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