rivet is hosted by Hepforge, IPPP Durham
MC_HJETS.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analyses/MC_JetAnalysis.hh"
00003 #include "Rivet/Projections/ZFinder.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 #include "Rivet/Tools/Logging.hh"
00006 #include "Rivet/RivetYODA.hh"
00007 
00008 namespace Rivet {
00009 
00010 
00011   /// @brief MC validation analysis for higgs [-> tau tau] + jets events
00012   class MC_HJETS : public MC_JetAnalysis {
00013   public:
00014 
00015     /// Default constructor
00016     MC_HJETS()
00017       : MC_JetAnalysis("MC_HJETS", 4, "Jets")
00018     {    }
00019 
00020 
00021     /// @name Analysis methods
00022     //@{
00023 
00024     /// Book histograms
00025     void init() {
00026       FinalState fs;
00027       ZFinder hfinder(fs, -3.5, 3.5, 25.0*GeV, TAU, 115.0*GeV, 125.0*GeV, 0.0, false, false);
00028       addProjection(hfinder, "Hfinder");
00029       FastJets jetpro(hfinder.remainingFinalState(), FastJets::KT, 0.7);
00030       addProjection(jetpro, "Jets");
00031 
00032       _h_H_mass = bookHisto1D("H_mass", 50, 119.7, 120.3);
00033       _h_H_pT = bookHisto1D("H_pT", logspace(100, 1.0, 0.5*sqrtS()));
00034       _h_H_pT_peak = bookHisto1D("H_pT_peak", 25, 0.0, 25.0);
00035       _h_H_y = bookHisto1D("H_y", 40, -4.0, 4.0);
00036       _h_H_phi = bookHisto1D("H_phi", 25, 0.0, TWOPI);
00037       _h_H_jet1_deta = bookHisto1D("H_jet1_deta", 50, -5.0, 5.0);
00038       _h_H_jet1_dR = bookHisto1D("H_jet1_dR", 25, 0.5, 7.0);
00039       _h_lepton_pT = bookHisto1D("lepton_pT", logspace(100, 10.0, 0.25*sqrtS()));
00040       _h_lepton_eta = bookHisto1D("lepton_eta", 40, -4.0, 4.0);
00041 
00042       MC_JetAnalysis::init();
00043     }
00044 
00045 
00046 
00047     /// Do the analysis
00048     void analyze(const Event & e) {
00049       const ZFinder& hfinder = applyProjection<ZFinder>(e, "Hfinder");
00050       if (hfinder.bosons().size()!=1) {
00051         vetoEvent;
00052       }
00053       const double weight = e.weight();
00054 
00055       FourMomentum hmom(hfinder.bosons()[0].momentum());
00056       _h_H_mass->fill(hmom.mass(),weight);
00057       _h_H_pT->fill(hmom.pT(),weight);
00058       _h_H_pT_peak->fill(hmom.pT(),weight);
00059       _h_H_y->fill(hmom.rapidity(),weight);
00060       _h_H_phi->fill(hmom.azimuthalAngle(),weight);
00061       foreach (const Particle& l, hfinder.constituents()) {
00062         _h_lepton_pT->fill(l.momentum().pT(), weight);
00063         _h_lepton_eta->fill(l.momentum().eta(), weight);
00064       }
00065 
00066       const FastJets& jetpro = applyProjection<FastJets>(e, "Jets");
00067       const Jets& jets = jetpro.jetsByPt(20.0*GeV);
00068       if (jets.size() > 0) {
00069         _h_H_jet1_deta->fill(hmom.eta()-jets[0].momentum().eta(), weight);
00070         _h_H_jet1_dR->fill(deltaR(hmom, jets[0].momentum()), weight);
00071       }
00072 
00073       MC_JetAnalysis::analyze(e);
00074     }
00075 
00076 
00077     /// Finalize
00078     void finalize() {
00079       scale(_h_H_mass, crossSection()/sumOfWeights());
00080       scale(_h_H_pT, crossSection()/sumOfWeights());
00081       scale(_h_H_pT_peak, crossSection()/sumOfWeights());
00082       scale(_h_H_y, crossSection()/sumOfWeights());
00083       scale(_h_H_phi, crossSection()/sumOfWeights());
00084       scale(_h_H_jet1_deta, crossSection()/sumOfWeights());
00085       scale(_h_H_jet1_dR, crossSection()/sumOfWeights());
00086       scale(_h_lepton_pT, crossSection()/sumOfWeights());
00087       scale(_h_lepton_eta, crossSection()/sumOfWeights());
00088 
00089       MC_JetAnalysis::finalize();
00090     }
00091 
00092     //@}
00093 
00094 
00095   private:
00096 
00097     /// @name Histograms
00098     //@{
00099     Histo1DPtr _h_H_mass;
00100     Histo1DPtr _h_H_pT;
00101     Histo1DPtr _h_H_pT_peak;
00102     Histo1DPtr _h_H_y;
00103     Histo1DPtr _h_H_phi;
00104     Histo1DPtr _h_H_jet1_deta;
00105     Histo1DPtr _h_H_jet1_dR;
00106     Histo1DPtr _h_lepton_pT;
00107     Histo1DPtr _h_lepton_eta;
00108     //@}
00109 
00110   };
00111 
00112 
00113 
00114   // The hook for the plugin system
00115   DECLARE_RIVET_PLUGIN(MC_HJETS);
00116 
00117 }