rivet is hosted by Hepforge, IPPP Durham
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   class CMS_2013_I1224539_ZJET : public Analysis {
00014   public:
00015 
00016     /// @name Constructors etc.
00017     //@{
00018 
00019     /// Constructor
00020     CMS_2013_I1224539_ZJET()
00021       : Analysis("CMS_2013_I1224539_ZJET"),
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 Zs with pT > 120 GeV
00041       ZFinder zfinder(fs, EtaIn(-2.4,2.4) & (Cuts::pT >= 30.0*GeV), PID::ELECTRON, 80*GeV, 100*GeV,
00042                       0.2, ZFinder::CLUSTERNODECAY, ZFinder::TRACK);
00043 
00044       addProjection(zfinder, "ZFinder");
00045 
00046       // Z+jet jet collections
00047       addProjection(FastJets(zfinder.remainingFinalState(), FastJets::ANTIKT, 0.7), "JetsAK7_zj");
00048       addProjection(FastJets(zfinder.remainingFinalState(), FastJets::CAM, 0.8), "JetsCA8_zj");
00049       addProjection(FastJets(zfinder.remainingFinalState(), FastJets::CAM, 1.2), "JetsCA12_zj");
00050 
00051       // Histograms
00052       /// @note These are 2D histos rendered into slices
00053       const int zjetsOffset = 28;
00054       for (size_t i = 0; i < N_PT_BINS_vj; ++i ) {
00055         _h_ungroomedJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+0*N_PT_BINS_vj, 1, 1);
00056         _h_filteredJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+1*N_PT_BINS_vj,1,1);
00057         _h_trimmedJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+2*N_PT_BINS_vj,1,1);
00058         _h_prunedJetMass_AK7_zj[i] = bookHisto1D(zjetsOffset+i+1+3*N_PT_BINS_vj,1,1);
00059         _h_prunedJetMass_CA8_zj[i] = bookHisto1D(zjetsOffset+i+1+4*N_PT_BINS_vj,1,1);
00060         if (i > 0) _h_filteredJetMass_CA12_zj[i] = bookHisto1D(zjetsOffset+i+5*N_PT_BINS_vj,1,1);
00061       }
00062     }
00063 
00064 
00065     bool isBackToBack_zj(const ZFinder& zf, const fastjet::PseudoJet& psjet) {
00066       const FourMomentum& z = zf.bosons()[0].momentum();
00067       const FourMomentum& l1 = zf.constituents()[0].momentum();
00068       const FourMomentum& l2 = zf.constituents()[1].momentum();
00069       /// @todo We should make FourMomentum know how to construct itself from a PseudoJet
00070       const FourMomentum jmom(psjet.e(), psjet.px(), psjet.py(), psjet.pz());
00071       return (deltaPhi(z, jmom) > 2.0 && deltaR(l1, jmom) > 1.0 && deltaR(l2, jmom) > 1.0);
00072     }
00073 
00074 
00075     // Find the pT histogram bin index for value pt (in GeV), to hack a 2D histogram equivalent
00076     /// @todo Use a YODA axis/finder alg when available
00077     size_t findPtBin(double ptJ) {
00078       const double ptBins_vj[N_PT_BINS_vj+1] = { 125.0, 150.0, 220.0, 300.0, 450.0 };
00079       for (size_t ibin = 0; ibin < N_PT_BINS_vj; ++ibin) {
00080         if (inRange(ptJ, ptBins_vj[ibin], ptBins_vj[ibin+1])) return ibin;
00081       }
00082       return N_PT_BINS_vj;
00083     }
00084 
00085 
00086     /// Perform the per-event analysis
00087     void analyze(const Event& event) {
00088       const double weight = event.weight();
00089 
00090       // Get the Z
00091       const ZFinder& zfinder = applyProjection<ZFinder>(event, "ZFinder");
00092       if (zfinder.bosons().size() != 1) vetoEvent;
00093       const Particle& z = zfinder.bosons()[0];
00094       const Particle& l1 = zfinder.constituents()[0];
00095       const Particle& l2 = zfinder.constituents()[1];
00096 
00097       // Require a high-pT Z (and constituents)
00098       if (l1.momentum().pT() < 30*GeV || l2.momentum().pT() < 30*GeV || z.momentum().pT() > 120*GeV) vetoEvent;
00099 
00100       // AK7 jets
00101       const PseudoJets& psjetsAK7_zj = applyProjection<FastJets>(event, "JetsAK7_zj").pseudoJetsByPt(50.0*GeV);
00102       if (!psjetsAK7_zj.empty()) {
00103         // Get the leading jet and make sure it's back-to-back with the Z
00104         const fastjet::PseudoJet& j0 = psjetsAK7_zj[0];
00105         if (isBackToBack_zj(zfinder, j0)) {
00106           const size_t njetBin = findPtBin(j0.pt()/GeV);
00107           if (njetBin < N_PT_BINS_vj) {
00108             fastjet::PseudoJet filtered0 = _filter(j0);
00109             fastjet::PseudoJet trimmed0 = _trimmer(j0);
00110             fastjet::PseudoJet pruned0 = _pruner(j0);
00111             _h_ungroomedJetMass_AK7_zj[njetBin]->fill(j0.m()/GeV, weight);
00112             _h_filteredJetMass_AK7_zj[njetBin]->fill(filtered0.m()/GeV, weight);
00113             _h_trimmedJetMass_AK7_zj[njetBin]->fill(trimmed0.m()/GeV, weight);
00114             _h_prunedJetMass_AK7_zj[njetBin]->fill(pruned0.m()/GeV, weight);
00115           }
00116         }
00117       }
00118 
00119       // CA8 jets
00120       const PseudoJets& psjetsCA8_zj = applyProjection<FastJets>(event, "JetsCA8_zj").pseudoJetsByPt(50.0*GeV);
00121       if (!psjetsCA8_zj.empty()) {
00122         // Get the leading jet and make sure it's back-to-back with the Z
00123         const fastjet::PseudoJet& j0 = psjetsCA8_zj[0];
00124         if (isBackToBack_zj(zfinder, j0)) {
00125           const size_t njetBin = findPtBin(j0.pt()/GeV);
00126           if (njetBin < N_PT_BINS_vj) {
00127             fastjet::PseudoJet pruned0 = _pruner(j0);
00128             _h_prunedJetMass_CA8_zj[njetBin]->fill(pruned0.m()/GeV, weight);
00129           }
00130         }
00131       }
00132 
00133       // CA12 jets
00134       const PseudoJets& psjetsCA12_zj = applyProjection<FastJets>(event, "JetsCA12_zj").pseudoJetsByPt(50.0*GeV);
00135       if (!psjetsCA12_zj.empty()) {
00136         // Get the leading jet and make sure it's back-to-back with the Z
00137         const fastjet::PseudoJet& j0 = psjetsCA12_zj[0];
00138         if (isBackToBack_zj(zfinder, j0)) {
00139           const size_t njetBin = findPtBin(j0.pt()/GeV);
00140           if (njetBin < N_PT_BINS_vj) {
00141             fastjet::PseudoJet filtered0 = _filter(j0);
00142             _h_filteredJetMass_CA12_zj[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_zj[i], normalizationVal);
00155         normalize( _h_filteredJetMass_AK7_zj[i], normalizationVal);
00156         normalize( _h_trimmedJetMass_AK7_zj[i], normalizationVal);
00157         normalize( _h_prunedJetMass_AK7_zj[i], normalizationVal);
00158         normalize( _h_prunedJetMass_CA8_zj[i], normalizationVal);
00159         if (i > 0) normalize( _h_filteredJetMass_CA12_zj[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 { PT_125_150_vj=0, PT_150_220_vj, PT_220_300_vj, PT_300_450_vj, N_PT_BINS_vj } BINS_vj;
00179     Histo1DPtr _h_ungroomedJetMass_AK7_zj[N_PT_BINS_vj];
00180     Histo1DPtr _h_filteredJetMass_AK7_zj[N_PT_BINS_vj];
00181     Histo1DPtr _h_trimmedJetMass_AK7_zj[N_PT_BINS_vj];
00182     Histo1DPtr _h_prunedJetMass_AK7_zj[N_PT_BINS_vj];
00183     Histo1DPtr _h_prunedJetMass_CA8_zj[N_PT_BINS_vj];
00184     Histo1DPtr _h_filteredJetMass_CA12_zj[N_PT_BINS_vj-1];
00185     //@}
00186 
00187   };
00188 
00189 
00190   // The hook for the plugin system
00191   DECLARE_RIVET_PLUGIN(CMS_2013_I1224539_ZJET);
00192 
00193 }