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     /// @name Construction
00015     //@{
00016     /// Constructor
00017     D0_2009_S8202443() : Analysis("D0_2009_S8202443"),
00018         _sum_of_weights(0.0), _sum_of_weights_constrained(0.0)
00019     {
00020     }
00021 
00022     //@}
00023 
00024 
00025     /// @name Analysis methods
00026     //@{
00027 
00028     /// Book histograms
00029     void init() {
00030       FinalState fs;
00031       // Leptons in constrained tracking acceptance
00032       vector<pair<double, double> > etaRanges;
00033       etaRanges.push_back(make_pair(-2.5, -1.5));
00034       etaRanges.push_back(make_pair(-1.1, 1.1));
00035       etaRanges.push_back(make_pair(1.5, 2.5));
00036       ZFinder zfinder_constrained(fs, etaRanges, 25.0*GeV, PID::ELECTRON,
00037                                   65.0*GeV, 115.0*GeV, 0.2, true, true);
00038       addProjection(zfinder_constrained, "ZFinderConstrained");
00039       FastJets conefinder_constrained(zfinder_constrained.remainingFinalState(),
00040                                       FastJets::D0ILCONE, 0.5);
00041       addProjection(conefinder_constrained, "ConeFinderConstrained");
00042 
00043       // Unconstrained leptons
00044       ZFinder zfinder(fs, -MAXRAPIDITY, MAXRAPIDITY, 0.0*GeV, PID::ELECTRON,
00045                       65.0*GeV, 115.0*GeV, 0.2, true, true);
00046       addProjection(zfinder, "ZFinder");
00047       FastJets conefinder(zfinder.remainingFinalState(), FastJets::D0ILCONE, 0.5);
00048       addProjection(conefinder, "ConeFinder");
00049 
00050       _h_jet1_pT_constrained = bookHisto1D(1, 1, 1);
00051       _h_jet2_pT_constrained = bookHisto1D(3, 1, 1);
00052       _h_jet3_pT_constrained = bookHisto1D(5, 1, 1);
00053       _h_jet1_pT = bookHisto1D(2, 1, 1);
00054       _h_jet2_pT = bookHisto1D(4, 1, 1);
00055       _h_jet3_pT = bookHisto1D(6, 1, 1);
00056     }
00057 
00058 
00059 
00060     // Do the analysis
00061     void analyze(const Event& e) {
00062       double weight = e.weight();
00063 
00064       // unconstrained electrons first
00065       const ZFinder& zfinder = applyProjection<ZFinder>(e, "ZFinder");
00066       if (zfinder.bosons().size()==1) {
00067         _sum_of_weights += weight;
00068         const JetAlg& jetpro = applyProjection<JetAlg>(e, "ConeFinder");
00069         const Jets& jets = jetpro.jetsByPt(20.0*GeV);
00070         Jets jets_cut;
00071         foreach (const Jet& j, jets) {
00072           if (fabs(j.momentum().pseudorapidity()) < 2.5) {
00073             jets_cut.push_back(j);
00074           }
00075         }
00076 
00077         if (jets_cut.size()>0) {
00078           _h_jet1_pT->fill(jets_cut[0].pT()/GeV, weight);
00079         }
00080         if (jets_cut.size()>1) {
00081           _h_jet2_pT->fill(jets_cut[1].pT()/GeV, weight);
00082         }
00083         if (jets_cut.size()>2) {
00084           _h_jet3_pT->fill(jets_cut[2].pT()/GeV, weight);
00085         }
00086       }
00087       else {
00088         MSG_DEBUG("no unique lepton pair found.");
00089       }
00090 
00091 
00092       // constrained electrons
00093       const ZFinder& zfinder_constrained = applyProjection<ZFinder>(e, "ZFinderConstrained");
00094       if (zfinder_constrained.bosons().size()==1) {
00095         _sum_of_weights_constrained += weight;
00096         const JetAlg& jetpro = applyProjection<JetAlg>(e, "ConeFinderConstrained");
00097         const Jets& jets = jetpro.jetsByPt(20.0*GeV);
00098         Jets jets_cut;
00099         foreach (const Jet& j, jets) {
00100           if (fabs(j.momentum().pseudorapidity()) < 2.5) {
00101             jets_cut.push_back(j);
00102           }
00103         }
00104 
00105         if (jets_cut.size()>0) {
00106           _h_jet1_pT_constrained->fill(jets_cut[0].pT()/GeV, weight);
00107         }
00108         if (jets_cut.size()>1) {
00109           _h_jet2_pT_constrained->fill(jets_cut[1].pT()/GeV, weight);
00110         }
00111         if (jets_cut.size()>2) {
00112           _h_jet3_pT_constrained->fill(jets_cut[2].pT()/GeV, weight);
00113         }
00114       }
00115       else {
00116         MSG_DEBUG("no unique lepton pair found.");
00117         vetoEvent;
00118       }
00119     }
00120 
00121 
00122 
00123     // Finalize
00124     void finalize() {
00125       scale(_h_jet1_pT, 1.0/_sum_of_weights);
00126       scale(_h_jet2_pT, 1.0/_sum_of_weights);
00127       scale(_h_jet3_pT, 1.0/_sum_of_weights);
00128       scale(_h_jet1_pT_constrained, 1.0/_sum_of_weights_constrained);
00129       scale(_h_jet2_pT_constrained, 1.0/_sum_of_weights_constrained);
00130       scale(_h_jet3_pT_constrained, 1.0/_sum_of_weights_constrained);
00131     }
00132 
00133     //@}
00134 
00135 
00136   private:
00137 
00138     /// @name Histograms
00139     //@{
00140     Histo1DPtr _h_jet1_pT;
00141     Histo1DPtr _h_jet2_pT;
00142     Histo1DPtr _h_jet3_pT;
00143     Histo1DPtr _h_jet1_pT_constrained;
00144     Histo1DPtr _h_jet2_pT_constrained;
00145     Histo1DPtr _h_jet3_pT_constrained;
00146     //@}
00147 
00148     double _sum_of_weights, _sum_of_weights_constrained;
00149 
00150   };
00151 
00152 
00153 
00154   // The hook for the plugin system
00155   DECLARE_RIVET_PLUGIN(D0_2009_S8202443);
00156 
00157 }