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