rivet is hosted by Hepforge, IPPP Durham
CMS_2013_I1224539_DIJET.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_DIJET : public Analysis {
00014   public:
00015 
00016     /// @name Constructors etc.
00017     //@{
00018 
00019     /// Constructor
00020     CMS_2013_I1224539_DIJET()
00021       : Analysis("CMS_2013_I1224539_DIJET"),
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       // Jet collections
00041       addProjection(FastJets(fs, FastJets::ANTIKT, 0.7), "JetsAK7");
00042       addProjection(FastJets(fs, FastJets::CAM, 0.8), "JetsCA8");
00043       addProjection(FastJets(fs, FastJets::CAM, 1.2), "JetsCA12");
00044 
00045       // Histograms
00046       for (size_t i = 0; i < N_PT_BINS_dj; ++i ) {
00047         _h_ungroomedAvgJetMass_dj[i] = bookHisto1D(i+1+0*N_PT_BINS_dj, 1, 1);
00048         _h_filteredAvgJetMass_dj[i] = bookHisto1D(i+1+1*N_PT_BINS_dj, 1, 1);
00049         _h_trimmedAvgJetMass_dj[i] = bookHisto1D(i+1+2*N_PT_BINS_dj, 1, 1);
00050         _h_prunedAvgJetMass_dj[i] = bookHisto1D(i+1+3*N_PT_BINS_dj, 1, 1);
00051       }
00052     }
00053 
00054 
00055     // Find the pT histogram bin index for value pt (in GeV), to hack a 2D histogram equivalent
00056     /// @todo Use a YODA axis/finder alg when available
00057     size_t findPtBin(double ptJ) {
00058       const double ptBins_dj[N_PT_BINS_dj+1] = { 220.0, 300.0, 450.0, 500.0, 600.0, 800.0, 1000.0, 1500.0};
00059       for (size_t ibin = 0; ibin < N_PT_BINS_dj; ++ibin) {
00060         if (inRange(ptJ, ptBins_dj[ibin], ptBins_dj[ibin+1])) return ibin;
00061       }
00062       return N_PT_BINS_dj;
00063     }
00064 
00065 
00066     /// Perform the per-event analysis
00067     void analyze(const Event& event) {
00068       const double weight = event.weight();
00069 
00070       // Look at events with >= 2 jets
00071       const PseudoJets& psjetsAK7 = applyProjection<FastJets>(event, "JetsAK7").pseudoJetsByPt( 50.0*GeV );
00072       if (psjetsAK7.size() < 2) vetoEvent;
00073 
00074       // Get the leading two jets and find their average pT
00075       const fastjet::PseudoJet& j0 = psjetsAK7[0];
00076       const fastjet::PseudoJet& j1 = psjetsAK7[1];
00077       double ptAvg = 0.5 * (j0.pt() + j1.pt());
00078 
00079       // Find the appropriate mean pT bin and escape if needed
00080       const size_t njetBin = findPtBin(ptAvg/GeV);
00081       if (njetBin >= N_PT_BINS_dj) vetoEvent;
00082 
00083       // Now run the substructure algs...
00084       fastjet::PseudoJet filtered0 = _filter(j0);
00085       fastjet::PseudoJet filtered1 = _filter(j1);
00086       fastjet::PseudoJet trimmed0 = _trimmer(j0);
00087       fastjet::PseudoJet trimmed1 = _trimmer(j1);
00088       fastjet::PseudoJet pruned0 = _pruner(j0);
00089       fastjet::PseudoJet pruned1 = _pruner(j1);
00090 
00091       // ... and fill the histograms
00092       _h_ungroomedAvgJetMass_dj[njetBin]->fill(0.5*(j0.m() + j1.m())/GeV, weight);
00093       _h_filteredAvgJetMass_dj[njetBin]->fill(0.5*(filtered0.m() + filtered1.m())/GeV, weight);
00094       _h_trimmedAvgJetMass_dj[njetBin]->fill(0.5*(trimmed0.m() + trimmed1.m())/GeV, weight);
00095       _h_prunedAvgJetMass_dj[njetBin]->fill(0.5*(pruned0.m() + pruned1.m())/GeV, weight);
00096     }
00097 
00098 
00099     /// Normalise histograms etc., after the run
00100     void finalize() {
00101       const double normalizationVal = 1000;
00102       for (size_t i = 0; i < N_PT_BINS_dj; ++i) {
00103         normalize(_h_ungroomedAvgJetMass_dj[i], normalizationVal);
00104         normalize(_h_filteredAvgJetMass_dj[i], normalizationVal);
00105         normalize(_h_trimmedAvgJetMass_dj[i], normalizationVal);
00106         normalize(_h_prunedAvgJetMass_dj[i], normalizationVal);
00107       }
00108     }
00109 
00110     //@}
00111 
00112 
00113   private:
00114 
00115     /// @name FastJet grooming tools (configured in constructor init list)
00116     //@{
00117     const fastjet::Filter _filter;
00118     const fastjet::Filter _trimmer;
00119     const fastjet::Pruner _pruner;
00120     //@}
00121 
00122 
00123     /// @name Histograms
00124     //@{
00125     enum { PT_220_300_dj=0, PT_300_450_dj, PT_450_500_dj, PT_500_600_dj,
00126            PT_600_800_dj, PT_800_1000_dj, PT_1000_1500_dj, N_PT_BINS_dj } BINS_dj;
00127     Histo1DPtr _h_ungroomedJet0pt, _h_ungroomedJet1pt;
00128     Histo1DPtr _h_ungroomedAvgJetMass_dj[N_PT_BINS_dj];
00129     Histo1DPtr _h_filteredAvgJetMass_dj[N_PT_BINS_dj];
00130     Histo1DPtr _h_trimmedAvgJetMass_dj[N_PT_BINS_dj];
00131     Histo1DPtr _h_prunedAvgJetMass_dj[N_PT_BINS_dj];
00132     //@}
00133 
00134 
00135   };
00136 
00137 
00138 
00139   // The hook for the plugin system
00140   DECLARE_RIVET_PLUGIN(CMS_2013_I1224539_DIJET);
00141 
00142 
00143 }