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].pT(); 00070 double Lep2Pt = ZDecayProducts[1].pT(); 00071 double Lep1Eta = ZDecayProducts[0].absrap(); ///< @todo This is y... should be abseta()? 00072 double Lep2Eta = ZDecayProducts[1].absrap(); ///< @todo This is y... should be abseta()? 00073 00074 if (Lep1Eta > _LepEtaCut && Lep2Eta > _LepEtaCut) vetoEvent; 00075 if (ZDecayProducts[0].abspid()==13 && Lep1Eta > 1. && Lep2Eta > 1.) vetoEvent; 00076 if (Lep1Pt < _Lep1PtCut && Lep2Pt < _Lep2PtCut) vetoEvent; 00077 00078 _sumWeightsWithZ += event.weight(); 00079 00080 /// @todo Write out a warning if there are more than two decay products 00081 FourMomentum Zmom = ZDecayProducts[0].momentum() + ZDecayProducts[1].momentum(); 00082 00083 // Put all b-quarks in a vector 00084 /// @todo Use jet contents rather than accessing quarks directly 00085 Particles bquarks; 00086 /// @todo Use nicer looping 00087 for (GenEvent::particle_const_iterator p = event.genEvent()->particles_begin(); p != event.genEvent()->particles_end(); ++p) { 00088 if ( std::abs((*p)->pdg_id()) == PID::BQUARK ) { 00089 bquarks.push_back(Particle(**p)); 00090 } 00091 } 00092 00093 // Get jets 00094 const FastJets& jetpro = applyProjection<FastJets>(event, "Jets"); 00095 MSG_DEBUG("Jet multiplicity before any pT cut = " << jetpro.size()); 00096 00097 const PseudoJets& jets = jetpro.pseudoJetsByPt(); 00098 MSG_DEBUG("jetlist size = " << jets.size()); 00099 00100 int numBJet = 0; 00101 int numJet = 0; 00102 // for each b-jet plot the ET and the eta of the jet, normalise to the total cross section at the end 00103 // for each event plot N jet and pT(Z), normalise to the total cross section at the end 00104 for (PseudoJets::const_iterator jt = jets.begin(); jt != jets.end(); ++jt) { 00105 // select jets that pass the kinematic cuts 00106 if (jt->perp() > _JetPtCut && fabs(jt->rapidity()) <= _JetEtaCut) { 00107 ++numJet; 00108 // Does the jet contain a b-quark? 00109 /// @todo Use jet contents rather than accessing quarks directly 00110 bool bjet = false; 00111 foreach (const Particle& bquark, bquarks) { 00112 if (deltaR(jt->rapidity(), jt->phi(), bquark.rapidity(), bquark.phi()) <= _Rjet) { 00113 bjet = true; 00114 break; 00115 } 00116 } // end loop around b-jets 00117 if (bjet) { 00118 numBJet++; 00119 } 00120 } 00121 } // end loop around jets 00122 00123 if (numJet > 0) _sumWeightsWithZJet += event.weight(); 00124 if (numBJet > 0) { 00125 _sigmaBJet->fill(1960.0,event.weight()); 00126 _ratioBJetToZ->fill(1960.0,event.weight()); 00127 _ratioBJetToJet->fill(1960.0,event.weight()); 00128 } 00129 00130 } 00131 00132 00133 /// Finalize 00134 void finalize() { 00135 MSG_DEBUG("Total sum of weights = " << sumOfWeights()); 00136 MSG_DEBUG("Sum of weights for Z production in mass range = " << _sumWeightsWithZ); 00137 MSG_DEBUG("Sum of weights for Z+jet production in mass range = " << _sumWeightsWithZJet); 00138 00139 scale(_sigmaBJet, crossSection()/sumOfWeights()); 00140 scale(_ratioBJetToZ, 1.0/_sumWeightsWithZ); 00141 scale(_ratioBJetToJet, 1.0/_sumWeightsWithZJet); 00142 } 00143 00144 //@} 00145 00146 00147 private: 00148 00149 /// @name Cuts and counters 00150 //@{ 00151 double _Rjet; 00152 double _JetPtCut; 00153 double _JetEtaCut; 00154 double _Lep1PtCut; 00155 double _Lep2PtCut; 00156 double _LepEtaCut; 00157 00158 double _sumWeightsWithZ; 00159 double _sumWeightsWithZJet; 00160 00161 //@} 00162 00163 /// @name Histograms 00164 //@{ 00165 Histo1DPtr _sigmaBJet; 00166 Histo1DPtr _ratioBJetToZ; 00167 Histo1DPtr _ratioBJetToJet; 00168 //@} 00169 00170 }; 00171 00172 00173 00174 // The hook for the plugin system 00175 DECLARE_RIVET_PLUGIN(CDF_2006_S6653332); 00176 00177 } Generated on Thu Mar 10 2016 08:29:48 for The Rivet MC analysis system by ![]() |