CDF_2008_S7540469.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/IdentifiedFinalState.hh" 00005 #include "Rivet/Projections/VetoedFinalState.hh" 00006 #include "Rivet/Projections/FastJets.hh" 00007 00008 namespace Rivet { 00009 00010 00011 /// @brief Measurement differential Z/\f$ \gamma^* \f$ + jet + \f$ X \f$ cross sections 00012 /// @author Frank Siegert 00013 class CDF_2008_S7540469 : public Analysis { 00014 00015 public: 00016 00017 /// Constructor 00018 CDF_2008_S7540469() 00019 : Analysis("CDF_2008_S7540469") 00020 { } 00021 00022 00023 /// @name Analysis methods 00024 //@{ 00025 00026 /// Book histograms 00027 void init() { 00028 // Full final state 00029 FinalState fs(-5.0, 5.0); 00030 addProjection(fs, "FS"); 00031 00032 // Leading electrons in tracking acceptance 00033 IdentifiedFinalState elfs(Cuts::abseta < 5 && Cuts::pT > 25*GeV); 00034 elfs.acceptIdPair(PID::ELECTRON); 00035 addProjection(elfs, "LeadingElectrons"); 00036 00037 _h_jet_multiplicity = bookHisto1D(1, 1, 1); 00038 _h_jet_pT_cross_section_incl_1jet = bookHisto1D(2, 1, 1); 00039 _h_jet_pT_cross_section_incl_2jet = bookHisto1D(3, 1, 1); 00040 } 00041 00042 00043 /// Do the analysis 00044 void analyze(const Event & event) { 00045 const double weight = event.weight(); 00046 00047 // Skip if the event is empty 00048 const FinalState& fs = applyProjection<FinalState>(event, "FS"); 00049 if (fs.empty()) { 00050 MSG_DEBUG("Skipping event " << numEvents() << " because no final state pair found"); 00051 vetoEvent; 00052 } 00053 00054 // Find the Z candidates 00055 const FinalState & electronfs = applyProjection<FinalState>(event, "LeadingElectrons"); 00056 std::vector<std::pair<Particle, Particle> > Z_candidates; 00057 Particles all_els=electronfs.particles(); 00058 for (size_t i=0; i<all_els.size(); ++i) { 00059 for (size_t j=i+1; j<all_els.size(); ++j) { 00060 bool candidate=true; 00061 double mZ = FourMomentum(all_els[i].momentum()+all_els[j].momentum()).mass()/GeV; 00062 if (mZ < 66.0 || mZ > 116.0) { 00063 candidate = false; 00064 } 00065 double abs_eta_0 = fabs(all_els[i].eta()); 00066 double abs_eta_1 = fabs(all_els[j].eta()); 00067 if (abs_eta_1 < abs_eta_0) { 00068 double tmp = abs_eta_0; 00069 abs_eta_0 = abs_eta_1; 00070 abs_eta_1 = tmp; 00071 } 00072 if (abs_eta_0 > 1.0) { 00073 candidate = false; 00074 } 00075 if (!(abs_eta_1 < 1.0 || (inRange(abs_eta_1, 1.2, 2.8)))) { 00076 candidate = false; 00077 } 00078 if (candidate) { 00079 Z_candidates.push_back(make_pair(all_els[i], all_els[j])); 00080 } 00081 } 00082 } 00083 if (Z_candidates.size() != 1) { 00084 MSG_DEBUG("Skipping event " << numEvents() << " because no unique electron pair found "); 00085 vetoEvent; 00086 } 00087 00088 // Now build the jets on a FS without the electrons from the Z 00089 // (including their QED radiation) 00090 Particles jetparts; 00091 foreach (const Particle& p, fs.particles()) { 00092 bool copy = true; 00093 if (p.pid() == PID::PHOTON) { 00094 FourMomentum p_e0 = Z_candidates[0].first.momentum(); 00095 FourMomentum p_e1 = Z_candidates[0].second.momentum(); 00096 FourMomentum p_P = p.momentum(); 00097 if (deltaR(p_e0.eta(), p_e0.phi(), 00098 p_P.eta(), p_P.phi()) < 0.2) { 00099 copy = false; 00100 } 00101 if (deltaR(p_e1.eta(), p_e1.phi(), 00102 p_P.eta(), p_P.phi()) < 0.2) { 00103 copy = false; 00104 } 00105 } else { 00106 if (p.genParticle()->barcode() == Z_candidates[0].first.genParticle()->barcode()) { 00107 copy = false; 00108 } 00109 if (p.genParticle()->barcode() == Z_candidates[0].second.genParticle()->barcode()) { 00110 copy = false; 00111 } 00112 } 00113 if (copy) jetparts.push_back(p); 00114 } 00115 FastJets jetpro(FastJets::CDFMIDPOINT, 0.7); 00116 jetpro.calc(jetparts); 00117 00118 // Take jets with pt > 30, |eta| < 2.1: 00119 /// @todo Make this neater, using the JetAlg interface and the built-in sorting 00120 const Jets& jets = jetpro.jets(); 00121 Jets jets_cut; 00122 foreach (const Jet& j, jets) { 00123 if (j.pT()/GeV > 30.0 && j.abseta() < 2.1) { 00124 jets_cut.push_back(j); 00125 } 00126 } 00127 MSG_DEBUG("Num jets above 30 GeV = " << jets_cut.size()); 00128 00129 // Return if there are no jets: 00130 if (jets_cut.empty()) { 00131 MSG_DEBUG("No jets pass cuts "); 00132 vetoEvent; 00133 } 00134 00135 // Sort by pT: 00136 sort(jets_cut.begin(), jets_cut.end(), cmpMomByPt); 00137 00138 // cut on Delta R between jet and electrons 00139 foreach (const Jet& j, jets_cut) { 00140 Particle el = Z_candidates[0].first; 00141 if (deltaR(el.eta(), el.phi(), 00142 j.eta(), j.phi()) < 0.7) { 00143 vetoEvent; 00144 } 00145 el = Z_candidates[0].second; 00146 if (deltaR(el.eta(), el.phi(), 00147 j.eta(), j.phi()) < 0.7) { 00148 vetoEvent; 00149 } 00150 } 00151 00152 for (size_t njet=1; njet<=jets_cut.size(); ++njet) { 00153 _h_jet_multiplicity->fill(njet, weight); 00154 } 00155 foreach (const Jet& j, jets_cut) { 00156 if (jets_cut.size()>0) { 00157 _h_jet_pT_cross_section_incl_1jet->fill(j.pT(), weight); 00158 } 00159 if (jets_cut.size()>1) { 00160 _h_jet_pT_cross_section_incl_2jet->fill(j.pT(), weight); 00161 } 00162 } 00163 } 00164 00165 00166 /// Rescale histos 00167 void finalize() { 00168 const double invlumi = crossSection()/femtobarn/sumOfWeights(); 00169 scale(_h_jet_multiplicity, invlumi); 00170 scale(_h_jet_pT_cross_section_incl_1jet, invlumi); 00171 scale(_h_jet_pT_cross_section_incl_2jet, invlumi); 00172 } 00173 00174 //@} 00175 00176 private: 00177 00178 /// @name Histograms 00179 //@{ 00180 Histo1DPtr _h_jet_multiplicity; 00181 Histo1DPtr _h_jet_pT_cross_section_incl_1jet; 00182 Histo1DPtr _h_jet_pT_cross_section_incl_2jet; 00183 //@} 00184 00185 }; 00186 00187 00188 00189 // The hook for the plugin system 00190 DECLARE_RIVET_PLUGIN(CDF_2008_S7540469); 00191 00192 } Generated on Wed Oct 7 2015 12:09:11 for The Rivet MC analysis system by ![]() |