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