CDF_2006_S6653332.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/RivetYODA.hh" 00004 #include "Rivet/Tools/Logging.hh" 00005 #include "Rivet/Tools/ParticleIdUtils.hh" 00006 #include "Rivet/Projections/FinalState.hh" 00007 #include "Rivet/Projections/VetoedFinalState.hh" 00008 #include "Rivet/Projections/ChargedFinalState.hh" 00009 #include "Rivet/Projections/InvMassFinalState.hh" 00010 #include "Rivet/Projections/FastJets.hh" 00011 #include "Rivet/Projections/ChargedLeptons.hh" 00012 00013 namespace Rivet { 00014 00015 00016 /// @brief CDF Run II analysis: jet \f$ p_T \f$ and \f$ \eta \f$ 00017 /// distributions in Z + (b) jet production 00018 /// @author Lars Sonnenschein 00019 /// 00020 /// This CDF analysis provides \f$ p_T \f$ and \f$ \eta \f$ distributions of 00021 /// jets in Z + (b) jet production, before and after tagging. 00022 class CDF_2006_S6653332 : public Analysis { 00023 public: 00024 00025 /// Constructor 00026 CDF_2006_S6653332() 00027 : Analysis("CDF_2006_S6653332"), 00028 _Rjet(0.7), _JetPtCut(20.), _JetEtaCut(1.5), _Lep1PtCut(18.), _Lep2PtCut(10.), _LepEtaCut(1.1), 00029 _sumWeightsWithZ(0.0), _sumWeightsWithZJet(0.0) 00030 { } 00031 00032 00033 /// @name Analysis methods 00034 //@{ 00035 00036 void init() { 00037 const FinalState fs(-3.6, 3.6); 00038 addProjection(fs, "FS"); 00039 00040 // Create a final state with any e+e- or mu+mu- pair with 00041 // invariant mass 76 -> 106 GeV and ET > 20 (Z decay products) 00042 vector<pair<PdgId,PdgId> > vids; 00043 vids.push_back(make_pair(ELECTRON, POSITRON)); 00044 vids.push_back(make_pair(MUON, ANTIMUON)); 00045 FinalState fs2(-3.6, 3.6); 00046 InvMassFinalState invfs(fs2, vids, 66*GeV, 116*GeV); 00047 addProjection(invfs, "INVFS"); 00048 00049 // Make a final state without the Z decay products for jet clustering 00050 VetoedFinalState vfs(fs); 00051 vfs.addVetoOnThisFinalState(invfs); 00052 addProjection(vfs, "VFS"); 00053 addProjection(FastJets(vfs, FastJets::CDFMIDPOINT, 0.7), "Jets"); 00054 00055 // Book histograms 00056 _sigmaBJet = bookHisto1D(1, 1, 1); 00057 _ratioBJetToZ = bookHisto1D(2, 1, 1); 00058 _ratioBJetToJet = bookHisto1D(3, 1, 1); 00059 } 00060 00061 00062 /// Do the analysis 00063 void analyze(const Event& event) { 00064 // Check we have an l+l- pair that passes the kinematic cuts 00065 // Get the Z decay products (mu+mu- or e+e- pair) 00066 const InvMassFinalState& invMassFinalState = applyProjection<InvMassFinalState>(event, "INVFS"); 00067 const ParticleVector& ZDecayProducts = invMassFinalState.particles(); 00068 00069 // Make sure we have at least 2 Z decay products (mumu or ee) 00070 if (ZDecayProducts.size() < 2) vetoEvent; 00071 // 00072 double Lep1Pt = ZDecayProducts[0].momentum().perp(); 00073 double Lep2Pt = ZDecayProducts[1].momentum().perp(); 00074 double Lep1Eta = fabs(ZDecayProducts[0].momentum().rapidity()); 00075 double Lep2Eta = fabs(ZDecayProducts[1].momentum().rapidity()); 00076 00077 if (Lep1Eta > _LepEtaCut && Lep2Eta > _LepEtaCut) vetoEvent; 00078 00079 if (abs(ZDecayProducts[0].pdgId())==13 && (Lep1Eta > 1.0 && Lep2Eta > 1.)) { 00080 vetoEvent; 00081 } 00082 if (Lep1Pt < _Lep1PtCut && Lep2Pt < _Lep1PtCut) vetoEvent; 00083 00084 // 00085 00086 _sumWeightsWithZ += event.weight(); 00087 // @todo: write out a warning if there are more than two decay products 00088 FourMomentum Zmom = ZDecayProducts[0].momentum() + ZDecayProducts[1].momentum(); 00089 00090 // Put all b-quarks in a vector 00091 /// @todo Use jet contents rather than accessing quarks directly 00092 ParticleVector bquarks; 00093 /// @todo Use nicer looping 00094 for (GenEvent::particle_const_iterator p = event.genEvent().particles_begin(); 00095 p != event.genEvent().particles_end(); ++p) { 00096 if ( fabs((*p)->pdg_id()) == BQUARK ) { 00097 bquarks.push_back(Particle(**p)); 00098 } 00099 } 00100 00101 // Get jets 00102 const FastJets& jetpro = applyProjection<FastJets>(event, "Jets"); 00103 MSG_DEBUG("Jet multiplicity before any pT cut = " << jetpro.size()); 00104 00105 const PseudoJets& jets = jetpro.pseudoJetsByPt(); 00106 MSG_DEBUG("jetlist size = " << jets.size()); 00107 00108 int numBJet = 0; 00109 int numJet = 0; 00110 // for each b-jet plot the ET and the eta of the jet, normalise to the total cross section at the end 00111 // for each event plot N jet and pT(Z), normalise to the total cross section at the end 00112 for (PseudoJets::const_iterator jt = jets.begin(); jt != jets.end(); ++jt) { 00113 // select jets that pass the kinematic cuts 00114 if (jt->perp() > _JetPtCut && fabs(jt->rapidity()) <= _JetEtaCut) { 00115 ++numJet; 00116 // Does the jet contain a b-quark? 00117 /// @todo Use jet contents rather than accessing quarks directly 00118 00119 bool bjet = false; 00120 foreach (const Particle& bquark, bquarks) { 00121 if (deltaR(jt->rapidity(), jt->phi(), bquark.momentum().rapidity(),bquark.momentum().azimuthalAngle()) <= _Rjet) { 00122 bjet = true; 00123 break; 00124 } 00125 } // end loop around b-jets 00126 if (bjet) { 00127 numBJet++; 00128 } 00129 } 00130 } // end loop around jets 00131 00132 if (numJet > 0) _sumWeightsWithZJet += event.weight(); 00133 if (numBJet > 0) { 00134 _sigmaBJet->fill(1960.0,event.weight()); 00135 _ratioBJetToZ->fill(1960.0,event.weight()); 00136 _ratioBJetToJet->fill(1960.0,event.weight()); 00137 } 00138 00139 } 00140 00141 00142 /// Finalize 00143 void finalize() { 00144 MSG_DEBUG("Total sum of weights = " << sumOfWeights()); 00145 MSG_DEBUG("Sum of weights for Z production in mass range = " << _sumWeightsWithZ); 00146 MSG_DEBUG("Sum of weights for Z+jet production in mass range = " << _sumWeightsWithZJet); 00147 00148 scale(_sigmaBJet,crossSection()/sumOfWeights()); 00149 scale(_ratioBJetToZ,1.0/_sumWeightsWithZ); 00150 scale(_ratioBJetToJet,1.0/_sumWeightsWithZJet); 00151 } 00152 00153 //@} 00154 00155 00156 private: 00157 00158 /// @name Cuts and counters 00159 //@{ 00160 00161 double _Rjet; 00162 double _JetPtCut; 00163 double _JetEtaCut; 00164 double _Lep1PtCut; 00165 double _Lep2PtCut; 00166 double _LepEtaCut; 00167 00168 double _sumWeightsWithZ; 00169 double _sumWeightsWithZJet; 00170 00171 //@} 00172 00173 /// @name Histograms 00174 //@{ 00175 Histo1DPtr _sigmaBJet; 00176 Histo1DPtr _ratioBJetToZ; 00177 Histo1DPtr _ratioBJetToJet; 00178 //@} 00179 00180 }; 00181 00182 00183 00184 // The hook for the plugin system 00185 DECLARE_RIVET_PLUGIN(CDF_2006_S6653332); 00186 00187 } Generated on Fri Dec 21 2012 15:03:39 for The Rivet MC analysis system by ![]() |