rivet is hosted by Hepforge, IPPP Durham
ATLAS_2010_S8919674.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/VetoedFinalState.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 #include "Rivet/Projections/WFinder.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   /// W + jets jet multiplicities and pT
00011   class ATLAS_2010_S8919674 : public Analysis {
00012   public:
00013 
00014     /// @name Constructors etc.
00015     //@{
00016 
00017     /// Constructor
00018     ATLAS_2010_S8919674()
00019       : Analysis("ATLAS_2010_S8919674")
00020     {    }
00021 
00022     //@}
00023 
00024 
00025   public:
00026 
00027     /// @name Analysis methods
00028     //@{
00029 
00030     /// Book histograms and initialise projections before the run
00031     void init() {
00032 
00033       // Set up projections to find the electron and muon Ws
00034       FinalState fs;
00035       vector<pair<double, double> > eta_e;
00036       /// @todo Use C++11 uniform init list
00037       eta_e.push_back(make_pair(-2.47, -1.52));
00038       eta_e.push_back(make_pair(-1.37, 1.37));
00039       eta_e.push_back(make_pair(1.52, 2.47));
00040       WFinder wfinder_e(fs, eta_e, 20*GeV, PID::ELECTRON, 0*GeV, 1000*GeV, 25*GeV);
00041       addProjection(wfinder_e, "W_e");
00042       WFinder wfinder_mu(fs, -2.4, 2.4, 20*GeV, PID::MUON, 0*GeV, 1000*GeV, 25*GeV);
00043       addProjection(wfinder_mu, "W_mu");
00044 
00045       // Input for the jets: no neutrinos, no muons, and no electron which passed the electron cuts
00046       VetoedFinalState veto;
00047       veto.addVetoOnThisFinalState(wfinder_e);
00048       veto.addVetoOnThisFinalState(wfinder_mu);
00049       veto.addVetoPairId(PID::MUON);
00050       veto.vetoNeutrinos();
00051       FastJets jets(veto, FastJets::ANTIKT, 0.4);
00052       addProjection(jets, "jets");
00053 
00054       /// Book histograms
00055       _h_el_njet_inclusive = bookHisto1D(1,1,1);
00056       _h_mu_njet_inclusive = bookHisto1D(2,1,1);
00057       _h_el_pT_jet1 = bookHisto1D(5,1,1);
00058       _h_mu_pT_jet1 = bookHisto1D(6,1,1);
00059       _h_el_pT_jet2 = bookHisto1D(7,1,1);
00060       _h_mu_pT_jet2 = bookHisto1D(8,1,1);
00061     }
00062 
00063 
00064     /// Perform the per-event analysis
00065     void analyze(const Event& event) {
00066       const double weight = event.weight();
00067 
00068       const Jets& jets = applyProjection<FastJets>(event, "jets").jetsByPt(20.0*GeV);
00069 
00070       const WFinder& We = applyProjection<WFinder>(event, "W_e");
00071       if (We.bosons().size() == 1) {
00072         const FourMomentum& p_miss = We.constituentNeutrinos()[0].momentum();
00073         const FourMomentum& p_lept = We.constituentLeptons()[0].momentum();
00074         if (p_miss.Et() > 25*GeV && We.mT() > 40*GeV) {
00075           Jets js;
00076           foreach (const Jet& j, jets) {
00077             if (fabs(j.eta()) < 2.8 && deltaR(p_lept, j.momentum()) > 0.5) 
00078               js.push_back(j);
00079           }
00080           _h_el_njet_inclusive->fill(0, weight);
00081           if (js.size() >= 1) {
00082             _h_el_njet_inclusive->fill(1, weight);
00083             _h_el_pT_jet1->fill(js[0].pT(), weight);
00084           }
00085           if (js.size() >= 2) {
00086             _h_el_njet_inclusive->fill(2, weight);
00087             _h_el_pT_jet2->fill(js[1].pT(), weight);
00088           }
00089           if (js.size() >= 3) {
00090             _h_el_njet_inclusive->fill(3, weight);
00091           }
00092         }
00093       }
00094 
00095       const WFinder& Wm = applyProjection<WFinder>(event, "W_mu");
00096       if (Wm.bosons().size() == 1) {
00097         const FourMomentum& p_miss = Wm.constituentNeutrinos()[0].momentum();
00098         const FourMomentum& p_lept = Wm.constituentLeptons()[0].momentum();
00099         if (p_miss.Et() > 25*GeV && Wm.mT() > 40*GeV) {
00100           Jets js;
00101           foreach (const Jet& j, jets) {
00102             if (fabs(j.eta()) < 2.8 && deltaR(p_lept, j.momentum()) > 0.5) 
00103               js.push_back(j);
00104           }
00105           _h_mu_njet_inclusive->fill(0, weight);
00106           if (js.size() >= 1) {
00107             _h_mu_njet_inclusive->fill(1, weight);
00108             _h_mu_pT_jet1->fill(js[0].pT(), weight);
00109           }
00110           if (js.size() >= 2) {
00111             _h_mu_njet_inclusive->fill(2, weight);
00112             _h_mu_pT_jet2->fill(js[1].pT(), weight);
00113           }
00114           if (js.size() >= 3) {
00115             _h_mu_njet_inclusive->fill(3, weight);
00116           }
00117           if (js.size() >= 4) {
00118             _h_mu_njet_inclusive->fill(4, weight);
00119           }
00120         }
00121       }
00122 
00123     }
00124 
00125 
00126     /// Normalise histograms etc., after the run
00127     void finalize() {
00128       double normfac = crossSection()/sumOfWeights();
00129       scale(_h_el_njet_inclusive, normfac);
00130       scale(_h_mu_njet_inclusive, normfac);
00131       scale(_h_el_pT_jet1, normfac);
00132       scale(_h_mu_pT_jet1, normfac);
00133       scale(_h_el_pT_jet2, normfac);
00134       scale(_h_mu_pT_jet2, normfac);
00135     }
00136 
00137     //@}
00138 
00139 
00140   private:
00141 
00142     /// @name Histograms
00143     //@{
00144 
00145     Histo1DPtr _h_el_njet_inclusive;
00146     Histo1DPtr _h_mu_njet_inclusive;
00147     Histo1DPtr _h_el_pT_jet1;
00148     Histo1DPtr _h_mu_pT_jet1;
00149     Histo1DPtr _h_el_pT_jet2;
00150     Histo1DPtr _h_mu_pT_jet2;
00151     //@}
00152 
00153   };
00154 
00155 
00156 
00157   // The hook for the plugin system
00158   DECLARE_RIVET_PLUGIN(ATLAS_2010_S8919674);
00159 
00160 }