CMS_2013_I1224539_ZJET.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/FastJets.hh" 00005 #include "Rivet/Projections/WFinder.hh" 00006 #include "Rivet/Projections/ZFinder.hh" 00007 #include "fastjet/tools/Filter.hh" 00008 #include "fastjet/tools/Pruner.hh" 00009 00010 namespace Rivet { 00011 00012 00013 00014 00015 class CMS_2013_I1224539_ZJET : public Analysis { 00016 public: 00017 00018 /// @name Constructors etc. 00019 //@{ 00020 00021 /// Constructor 00022 CMS_2013_I1224539_ZJET() 00023 : Analysis("CMS_2013_I1224539_ZJET"), 00024 _filter(fastjet::Filter(fastjet::JetDefinition(fastjet::cambridge_algorithm, 0.3), fastjet::SelectorNHardest(3))), 00025 _trimmer(fastjet::Filter(fastjet::JetDefinition(fastjet::kt_algorithm, 0.2), fastjet::SelectorPtFractionMin(0.03))), 00026 _pruner(fastjet::Pruner(fastjet::cambridge_algorithm, 0.1, 0.5)) 00027 { } 00028 00029 //@} 00030 00031 00032 public: 00033 00034 /// @name Analysis methods 00035 //@{ 00036 00037 /// Book histograms and initialise projections before the run 00038 void init() { 00039 FinalState fs(Cuts::abseta < 2.4); 00040 addProjection(fs, "FS"); 00041 00042 // Find Zs with pT > 120 GeV 00043 ZFinder zfinder(fs, Cuts::abseta < 2.4 && Cuts::pT > 30*GeV, PID::ELECTRON, 80*GeV, 100*GeV, 00044 0.2, ZFinder::CLUSTERNODECAY, ZFinder::TRACK); 00045 00046 addProjection(zfinder, "ZFinder"); 00047 00048 // Z+jet jet collections 00049 addProjection(FastJets(zfinder.remainingFinalState(), FastJets::ANTIKT, 0.7), "JetsAK7_zj"); 00050 addProjection(FastJets(zfinder.remainingFinalState(), FastJets::CAM, 0.8), "JetsCA8_zj"); 00051 addProjection(FastJets(zfinder.remainingFinalState(), FastJets::CAM, 1.2), "JetsCA12_zj"); 00052 00053 // Histograms 00054 /// @note These are 2D histos rendered into slices 00055 const int zjetsOffset = 28; 00056 for (size_t i = 0; i < N_PT_BINS_vj; ++i ) { 00057 _h_ungroomedJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+0*N_PT_BINS_vj, 1, 1); 00058 _h_filteredJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+1*N_PT_BINS_vj,1,1); 00059 _h_trimmedJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+2*N_PT_BINS_vj,1,1); 00060 _h_prunedJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+3*N_PT_BINS_vj,1,1); 00061 _h_prunedJetMass_CA8_zj[i] = bookHisto1D(zjetsOffset+i+1+4*N_PT_BINS_vj,1,1); 00062 if (i > 0) _h_filteredJetMass_CA12_zj[i] = bookHisto1D(zjetsOffset+i+5*N_PT_BINS_vj,1,1); 00063 } 00064 } 00065 00066 00067 bool isBackToBack_zj(const ZFinder& zf, const fastjet::PseudoJet& psjet) { 00068 const FourMomentum& z = zf.bosons()[0].momentum(); 00069 const FourMomentum& l1 = zf.constituents()[0].momentum(); 00070 const FourMomentum& l2 = zf.constituents()[1].momentum(); 00071 /// @todo We should make FourMomentum know how to construct itself from a PseudoJet 00072 const FourMomentum jmom(psjet.e(), psjet.px(), psjet.py(), psjet.pz()); 00073 return (deltaPhi(z, jmom) > 2.0 && deltaR(l1, jmom) > 1.0 && deltaR(l2, jmom) > 1.0); 00074 } 00075 00076 00077 // Find the pT histogram bin index for value pt (in GeV), to hack a 2D histogram equivalent 00078 /// @todo Use a YODA axis/finder alg when available 00079 size_t findPtBin(double ptJ) { 00080 const double ptBins_vj[N_PT_BINS_vj+1] = { 125.0, 150.0, 220.0, 300.0, 450.0 }; 00081 for (size_t ibin = 0; ibin < N_PT_BINS_vj; ++ibin) { 00082 if (inRange(ptJ, ptBins_vj[ibin], ptBins_vj[ibin+1])) return ibin; 00083 } 00084 return N_PT_BINS_vj; 00085 } 00086 00087 00088 /// Perform the per-event analysis 00089 void analyze(const Event& event) { 00090 const double weight = event.weight(); 00091 00092 // Get the Z 00093 const ZFinder& zfinder = applyProjection<ZFinder>(event, "ZFinder"); 00094 if (zfinder.bosons().size() != 1) vetoEvent; 00095 const Particle& z = zfinder.bosons()[0]; 00096 const Particle& l1 = zfinder.constituents()[0]; 00097 const Particle& l2 = zfinder.constituents()[1]; 00098 00099 // Require a high-pT Z (and constituents) 00100 if (l1.pT() < 30*GeV || l2.pT() < 30*GeV || z.pT() < 120*GeV) vetoEvent; 00101 00102 // AK7 jets 00103 const PseudoJets& psjetsAK7_zj = applyProjection<FastJets>(event, "JetsAK7_zj").pseudoJetsByPt(50.0*GeV); 00104 if (!psjetsAK7_zj.empty()) { 00105 // Get the leading jet and make sure it's back-to-back with the Z 00106 const fastjet::PseudoJet& j0 = psjetsAK7_zj[0]; 00107 if (isBackToBack_zj(zfinder, j0)) { 00108 const size_t njetBin = findPtBin(j0.pt()/GeV); 00109 if (njetBin < N_PT_BINS_vj) { 00110 fastjet::PseudoJet filtered0 = _filter(j0); 00111 fastjet::PseudoJet trimmed0 = _trimmer(j0); 00112 fastjet::PseudoJet pruned0 = _pruner(j0); 00113 _h_ungroomedJetMass_AK7_zj[njetBin]->fill(j0.m()/GeV, weight); 00114 _h_filteredJetMass_AK7_zj[njetBin]->fill(filtered0.m()/GeV, weight); 00115 _h_trimmedJetMass_AK7_zj[njetBin]->fill(trimmed0.m()/GeV, weight); 00116 _h_prunedJetMass_AK7_zj[njetBin]->fill(pruned0.m()/GeV, weight); 00117 } 00118 } 00119 } 00120 00121 // CA8 jets 00122 const PseudoJets& psjetsCA8_zj = applyProjection<FastJets>(event, "JetsCA8_zj").pseudoJetsByPt(50.0*GeV); 00123 if (!psjetsCA8_zj.empty()) { 00124 // Get the leading jet and make sure it's back-to-back with the Z 00125 const fastjet::PseudoJet& j0 = psjetsCA8_zj[0]; 00126 if (isBackToBack_zj(zfinder, j0)) { 00127 const size_t njetBin = findPtBin(j0.pt()/GeV); 00128 if (njetBin < N_PT_BINS_vj) { 00129 fastjet::PseudoJet pruned0 = _pruner(j0); 00130 _h_prunedJetMass_CA8_zj[njetBin]->fill(pruned0.m()/GeV, weight); 00131 } 00132 } 00133 } 00134 00135 // CA12 jets 00136 const PseudoJets& psjetsCA12_zj = applyProjection<FastJets>(event, "JetsCA12_zj").pseudoJetsByPt(50.0*GeV); 00137 if (!psjetsCA12_zj.empty()) { 00138 // Get the leading jet and make sure it's back-to-back with the Z 00139 const fastjet::PseudoJet& j0 = psjetsCA12_zj[0]; 00140 if (isBackToBack_zj(zfinder, j0)) { 00141 const size_t njetBin = findPtBin(j0.pt()/GeV); 00142 if (njetBin>0 && njetBin < N_PT_BINS_vj) { 00143 fastjet::PseudoJet filtered0 = _filter(j0); 00144 _h_filteredJetMass_CA12_zj[njetBin]->fill( filtered0.m() / GeV, weight); 00145 } 00146 } 00147 } 00148 00149 } 00150 00151 00152 /// Normalise histograms etc., after the run 00153 void finalize() { 00154 const double normalizationVal = 1000; 00155 for (size_t i = 0; i < N_PT_BINS_vj; ++i ) { 00156 normalize( _h_ungroomedJetMass_AK7_zj[i], normalizationVal); 00157 normalize( _h_filteredJetMass_AK7_zj[i], normalizationVal); 00158 normalize( _h_trimmedJetMass_AK7_zj[i], normalizationVal); 00159 normalize( _h_prunedJetMass_AK7_zj[i], normalizationVal); 00160 normalize( _h_prunedJetMass_CA8_zj[i], normalizationVal); 00161 if (i > 0) normalize( _h_filteredJetMass_CA12_zj[i], normalizationVal); 00162 } 00163 } 00164 00165 //@} 00166 00167 00168 private: 00169 00170 /// @name FastJet grooming tools (configured in constructor init list) 00171 //@{ 00172 const fastjet::Filter _filter; 00173 const fastjet::Filter _trimmer; 00174 const fastjet::Pruner _pruner; 00175 //@} 00176 00177 00178 /// @name Histograms 00179 //@{ 00180 enum BINS_vj { PT_125_150_vj=0, PT_150_220_vj, PT_220_300_vj, PT_300_450_vj, N_PT_BINS_vj }; 00181 Histo1DPtr _h_ungroomedJetMass_AK7_zj[N_PT_BINS_vj]; 00182 Histo1DPtr _h_filteredJetMass_AK7_zj[N_PT_BINS_vj]; 00183 Histo1DPtr _h_trimmedJetMass_AK7_zj[N_PT_BINS_vj]; 00184 Histo1DPtr _h_prunedJetMass_AK7_zj[N_PT_BINS_vj]; 00185 Histo1DPtr _h_prunedJetMass_CA8_zj[N_PT_BINS_vj]; 00186 Histo1DPtr _h_filteredJetMass_CA12_zj[N_PT_BINS_vj]; 00187 //@} 00188 00189 }; 00190 00191 00192 // The hook for the plugin system 00193 DECLARE_RIVET_PLUGIN(CMS_2013_I1224539_ZJET); 00194 00195 } Generated on Thu Mar 10 2016 08:29:49 for The Rivet MC analysis system by ![]() |