rivet is hosted by Hepforge, IPPP Durham
D0_2010_S8566488.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Tools/BinnedHistogram.hh"
00004 #include "Rivet/Projections/FinalState.hh"
00005 #include "Rivet/Projections/FastJets.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   /// @brief D0 dijet invariant mass measurement
00011   class D0_2010_S8566488 : public Analysis {
00012   public:
00013 
00014     /// @name Constructors etc.
00015     //@{
00016 
00017     /// Constructor
00018     D0_2010_S8566488()
00019       : Analysis("D0_2010_S8566488")
00020     {    }
00021 
00022     //@}
00023 
00024 
00025   public:
00026 
00027     /// @name Analysis methods
00028     //@{
00029 
00030     /// Book histograms and initialise projections before the run
00031     void init() {
00032 
00033       FinalState fs;
00034       FastJets conefinder(fs, FastJets::D0ILCONE, 0.7);
00035       addProjection(conefinder, "ConeFinder");
00036 
00037       _h_m_dijet.addHistogram(0.0, 0.4, bookHisto1D(1, 1, 1));
00038       _h_m_dijet.addHistogram(0.4, 0.8, bookHisto1D(2, 1, 1));
00039       _h_m_dijet.addHistogram(0.8, 1.2, bookHisto1D(3, 1, 1));
00040       _h_m_dijet.addHistogram(1.2, 1.6, bookHisto1D(4, 1, 1));
00041       _h_m_dijet.addHistogram(1.6, 2.0, bookHisto1D(5, 1, 1));
00042       _h_m_dijet.addHistogram(2.0, 2.4, bookHisto1D(6, 1, 1));
00043     }
00044 
00045 
00046     /// Perform the per-event analysis
00047     void analyze(const Event& e) {
00048       const double weight = e.weight();
00049 
00050       const Jets& jets = applyProjection<JetAlg>(e, "ConeFinder").jetsByPt(40.0*GeV);
00051       if (jets.size() < 2) vetoEvent;
00052 
00053       FourMomentum j0(jets[0].momentum());
00054       FourMomentum j1(jets[1].momentum());
00055       double ymax = std::max(fabs(j0.rapidity()), fabs(j1.rapidity()));
00056       double mjj = FourMomentum(j0+j1).mass();
00057 
00058       _h_m_dijet.fill(ymax, mjj/TeV, weight);
00059     }
00060 
00061 
00062     /// Normalise histograms etc., after the run
00063     void finalize() {
00064       _h_m_dijet.scale(crossSection()/sumOfWeights(), this);
00065     }
00066 
00067     //@}
00068 
00069 
00070   private:
00071 
00072     // Data members like post-cuts event weight counters go here
00073 
00074 
00075   private:
00076 
00077     /// @name Histograms
00078     //@{
00079     BinnedHistogram<double> _h_m_dijet;
00080     //@}
00081 
00082   };
00083 
00084 
00085 
00086   // The hook for the plugin system
00087   DECLARE_RIVET_PLUGIN(D0_2010_S8566488);
00088 
00089 }