D0_2008_S6879055.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/RivetYODA.hh" 00004 #include "Rivet/Tools/Logging.hh" 00005 #include "Rivet/Projections/ZFinder.hh" 00006 #include "Rivet/Projections/FastJets.hh" 00007 00008 namespace Rivet { 00009 00010 00011 /// @brief D0 measurement of the ratio \f$ \sigma(Z/\gamma^* + n \text{ jets})/\sigma(Z/\gamma^*) \f$ 00012 class D0_2008_S6879055 : public Analysis { 00013 public: 00014 00015 /// Default constructor. 00016 D0_2008_S6879055() : Analysis("D0_2008_S6879055") 00017 { 00018 } 00019 00020 00021 /// @name Analysis methods 00022 //@{ 00023 00024 // Book histograms 00025 void init() { 00026 FinalState fs; 00027 ZFinder zfinder(fs, -MAXRAPIDITY, MAXRAPIDITY, 0.0*GeV, ELECTRON, 00028 40.0*GeV, 200.0*GeV, 0.2, true, true); 00029 addProjection(zfinder, "ZFinder"); 00030 00031 FastJets conefinder(zfinder.remainingFinalState(), FastJets::D0ILCONE, 0.5); 00032 addProjection(conefinder, "ConeFinder"); 00033 00034 _crossSectionRatio = bookHisto1D(1, 1, 1); 00035 _pTjet1 = bookHisto1D(2, 1, 1); 00036 _pTjet2 = bookHisto1D(3, 1, 1); 00037 _pTjet3 = bookHisto1D(4, 1, 1); 00038 } 00039 00040 00041 00042 /// Do the analysis 00043 void analyze(const Event& event) { 00044 const double weight = event.weight(); 00045 00046 00047 00048 const ZFinder& zfinder = applyProjection<ZFinder>(event, "ZFinder"); 00049 if (zfinder.bosons().size()!=1) { 00050 vetoEvent; 00051 } 00052 00053 FourMomentum e0 = zfinder.constituents()[0].momentum(); 00054 FourMomentum e1 = zfinder.constituents()[1].momentum(); 00055 const double e0eta = e0.eta(); 00056 const double e0phi = e0.phi(); 00057 const double e1eta = e1.eta(); 00058 const double e1phi = e1.phi(); 00059 00060 vector<FourMomentum> finaljet_list; 00061 foreach (const Jet& j, applyProjection<JetAlg>(event, "ConeFinder").jetsByPt(20.0*GeV)) { 00062 const double jeta = j.momentum().eta(); 00063 const double jphi = j.momentum().phi(); 00064 if (fabs(jeta) < 2.5) { 00065 if (deltaR(e0eta, e0phi, jeta, jphi) > 0.4 && 00066 deltaR(e1eta, e1phi, jeta, jphi) > 0.4) { 00067 finaljet_list.push_back(j.momentum()); 00068 } 00069 } 00070 } 00071 00072 // For normalisation of crossSection data (includes events with no jets passing cuts) 00073 _crossSectionRatio->fill(0, weight); 00074 00075 // Fill jet pT and multiplicities 00076 if (finaljet_list.size() >= 1) { 00077 _crossSectionRatio->fill(1, weight); 00078 _pTjet1->fill(finaljet_list[0].pT(), weight); 00079 } 00080 if (finaljet_list.size() >= 2) { 00081 _crossSectionRatio->fill(2, weight); 00082 _pTjet2->fill(finaljet_list[1].pT(), weight); 00083 } 00084 if (finaljet_list.size() >= 3) { 00085 _crossSectionRatio->fill(3, weight); 00086 _pTjet3->fill(finaljet_list[2].pT(), weight); 00087 } 00088 if (finaljet_list.size() >= 4) { 00089 _crossSectionRatio->fill(4, weight); 00090 } 00091 } 00092 00093 00094 00095 /// Finalize 00096 void finalize() { 00097 // Now divide by the inclusive result 00098 scale(_crossSectionRatio,1.0/_crossSectionRatio->bin(0).area()); 00099 00100 // Normalise jet pTs to integrals of data 00101 // NB. There is no other way to do this, because these quantities are not 00102 // detector-corrected 00103 normalize(_pTjet1, 10439.0); // fixed norm OK 00104 normalize(_pTjet2, 1461.5); // fixed norm OK 00105 normalize(_pTjet3, 217.0); // fixed norm OK 00106 } 00107 00108 //@} 00109 00110 00111 private: 00112 00113 /// @name Histograms 00114 //@{ 00115 Histo1DPtr _crossSectionRatio; 00116 Histo1DPtr _pTjet1; 00117 Histo1DPtr _pTjet2; 00118 Histo1DPtr _pTjet3; 00119 //@} 00120 00121 }; 00122 00123 00124 00125 // The hook for the plugin system 00126 DECLARE_RIVET_PLUGIN(D0_2008_S6879055); 00127 00128 } Generated on Fri Dec 21 2012 15:03:40 for The Rivet MC analysis system by ![]() |