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