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