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