rivet is hosted by Hepforge, IPPP Durham
D0_2009_S8202443.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/ZFinder.hh"
00005 #include "Rivet/Projections/FastJets.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   /// @brief D0 Z + jet + \f$ X \f$ cross-section / \f$ p_\perp \f$ distributions
00011   class D0_2009_S8202443 : public Analysis {
00012   public:
00013 
00014     /// Constructor
00015     D0_2009_S8202443()
00016       : Analysis("D0_2009_S8202443"),
00017         _sum_of_weights(0), _sum_of_weights_constrained(0)
00018     {    }
00019 
00020 
00021     /// @name Analysis methods
00022     //@{
00023 
00024     /// Book histograms
00025     void init() {
00026       FinalState fs;
00027       // Leptons in constrained tracking acceptance
00028       Cut cuts = (Cuts::abseta < 1.1 || Cuts::absetaIn(1.5, 2.5)) && Cuts::pT > 25*GeV;
00029       ZFinder zfinder_constrained(fs, cuts, PID::ELECTRON, 65*GeV, 115*GeV, 0.2, ZFinder::CLUSTERNODECAY, ZFinder::TRACK);
00030       addProjection(zfinder_constrained, "ZFinderConstrained");
00031       FastJets conefinder_constrained(zfinder_constrained.remainingFinalState(), FastJets::D0ILCONE, 0.5);
00032       addProjection(conefinder_constrained, "ConeFinderConstrained");
00033 
00034       // Unconstrained leptons
00035       ZFinder zfinder(fs, Cuts::open(), PID::ELECTRON, 65*GeV, 115*GeV, 0.2, ZFinder::CLUSTERNODECAY, ZFinder::TRACK);
00036       addProjection(zfinder, "ZFinder");
00037       FastJets conefinder(zfinder.remainingFinalState(), FastJets::D0ILCONE, 0.5);
00038       addProjection(conefinder, "ConeFinder");
00039 
00040       _h_jet1_pT_constrained = bookHisto1D(1, 1, 1);
00041       _h_jet2_pT_constrained = bookHisto1D(3, 1, 1);
00042       _h_jet3_pT_constrained = bookHisto1D(5, 1, 1);
00043       _h_jet1_pT = bookHisto1D(2, 1, 1);
00044       _h_jet2_pT = bookHisto1D(4, 1, 1);
00045       _h_jet3_pT = bookHisto1D(6, 1, 1);
00046     }
00047 
00048 
00049 
00050     // Do the analysis
00051     void analyze(const Event& e) {
00052       double weight = e.weight();
00053 
00054       // Unconstrained electrons
00055       const ZFinder& zfinder = applyProjection<ZFinder>(e, "ZFinder");
00056       if (zfinder.bosons().size() == 0) {
00057         MSG_DEBUG("No unique lepton pair found.");
00058         vetoEvent;
00059       }
00060       _sum_of_weights += weight;
00061       const Jets jets_cut = applyProjection<JetAlg>(e, "ConeFinder").jetsByPt(Cuts::pT > 20*GeV && Cuts::abseta < 2.5);
00062       if (jets_cut.size() > 0)
00063         _h_jet1_pT->fill(jets_cut[0].pT()/GeV, weight);
00064       if (jets_cut.size() > 1)
00065         _h_jet2_pT->fill(jets_cut[1].pT()/GeV, weight);
00066       if (jets_cut.size() > 2)
00067         _h_jet3_pT->fill(jets_cut[2].pT()/GeV, weight);
00068 
00069 
00070       // Constrained electrons
00071       const ZFinder& zfinder_constrained = applyProjection<ZFinder>(e, "ZFinderConstrained");
00072       if (zfinder_constrained.bosons().size() == 0) {
00073         MSG_DEBUG("No unique constrained lepton pair found.");
00074         return; // Not really a "veto", since if we got this far there is an unconstrained Z
00075       }
00076       _sum_of_weights_constrained += weight;
00077       const Jets& jets_constrained = applyProjection<JetAlg>(e, "ConeFinderConstrained").jetsByPt(20*GeV);
00078       /// @todo Replace this explicit selection with a Cut
00079       Jets jets_cut_constrained;
00080       foreach (const Jet& j, jets_constrained) {
00081         if (j.abseta() < 2.5) jets_cut_constrained.push_back(j);
00082       }
00083       if (jets_cut_constrained.size() > 0)
00084         _h_jet1_pT_constrained->fill(jets_cut_constrained[0].pT()/GeV, weight);
00085       if (jets_cut_constrained.size() > 1)
00086         _h_jet2_pT_constrained->fill(jets_cut_constrained[1].pT()/GeV, weight);
00087       if (jets_cut_constrained.size() > 2)
00088         _h_jet3_pT_constrained->fill(jets_cut_constrained[2].pT()/GeV, weight);
00089     }
00090 
00091 
00092     // Finalize
00093     void finalize() {
00094       scale(_h_jet1_pT, 1/_sum_of_weights);
00095       scale(_h_jet2_pT, 1/_sum_of_weights);
00096       scale(_h_jet3_pT, 1/_sum_of_weights);
00097       scale(_h_jet1_pT_constrained, 1/_sum_of_weights_constrained);
00098       scale(_h_jet2_pT_constrained, 1/_sum_of_weights_constrained);
00099       scale(_h_jet3_pT_constrained, 1/_sum_of_weights_constrained);
00100     }
00101 
00102     //@}
00103 
00104 
00105   private:
00106 
00107     /// @name Histograms
00108     //@{
00109     Histo1DPtr _h_jet1_pT;
00110     Histo1DPtr _h_jet2_pT;
00111     Histo1DPtr _h_jet3_pT;
00112     Histo1DPtr _h_jet1_pT_constrained;
00113     Histo1DPtr _h_jet2_pT_constrained;
00114     Histo1DPtr _h_jet3_pT_constrained;
00115     //@}
00116 
00117     double _sum_of_weights, _sum_of_weights_constrained;
00118 
00119   };
00120 
00121 
00122 
00123   // The hook for the plugin system
00124   DECLARE_RIVET_PLUGIN(D0_2009_S8202443);
00125 
00126 }