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