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