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/IdentifiedFinalState.hh"
00004 #include "Rivet/Projections/VetoedFinalState.hh"
00005 #include "Rivet/Projections/MissingMomentum.hh"
00006 #include "Rivet/Projections/FastJets.hh"
00007 #include "Rivet/Projections/ClusteredPhotons.hh"
00008 #include "Rivet/Projections/LeadingParticlesFinalState.hh"
00009 
00010 namespace Rivet {
00011 
00012 
00013   class ATLAS_2010_S8919674 : public Analysis {
00014   public:
00015 
00016     /// @name Constructors etc.
00017     //@{
00018 
00019     /// Constructor
00020     ATLAS_2010_S8919674()
00021       : Analysis("ATLAS_2010_S8919674")
00022     {    }
00023 
00024     //@}
00025 
00026 
00027   public:
00028 
00029     /// @name Analysis methods
00030     //@{
00031 
00032     /// Book histograms and initialise projections before the run
00033     void init() {
00034 
00035       /// Initialise and register projections (selections on the final state)
00036       // projection to find the electrons
00037       std::vector<std::pair<double, double> > eta_e;
00038       eta_e.push_back(make_pair(-2.47,-1.52));
00039       eta_e.push_back(make_pair(-1.37,1.37));
00040       eta_e.push_back(make_pair(1.52,2.47));
00041       IdentifiedFinalState elecs(eta_e, 20.0*GeV);
00042       elecs.acceptIdPair(PID::ELECTRON);
00043       addProjection(elecs, "elecs");
00044       // projection for finding the photons which have to be clustered into
00045       // the lepton later
00046       ClusteredPhotons cphotons_e(FinalState(), elecs, 0.1);
00047       addProjection(cphotons_e, "cphotons_e");
00048 
00049       // projection to find the muons
00050       std::vector<std::pair<double, double> > eta_m;
00051       eta_m.push_back(make_pair(-2.4,2.4));
00052       IdentifiedFinalState muons(eta_m, 20.0*GeV);
00053       muons.acceptIdPair(PID::MUON);
00054       addProjection(muons, "muons");
00055       // projection for finding the photons which have to be clustered into
00056       // the lepton later
00057       ClusteredPhotons cphotons_m(FinalState(), muons, 0.1);
00058       addProjection(cphotons_m, "cphotons_m");
00059 
00060       // Leading neutrinos for Etmiss
00061       FinalState fs;
00062       LeadingParticlesFinalState muon_neutrino(fs);
00063       muon_neutrino.addParticleIdPair(PID::NU_MU);
00064       muon_neutrino.setLeadingOnly(true);
00065       addProjection(muon_neutrino, "muon_neutrino");
00066       LeadingParticlesFinalState elec_neutrino(fs);
00067       elec_neutrino.addParticleIdPair(PID::NU_E);
00068       elec_neutrino.setLeadingOnly(true);
00069       addProjection(elec_neutrino, "elec_neutrino");
00070 
00071       // Input for the jets: No neutrinos, no muons, and no electron which
00072       // passed the electron cuts ("elecs" finalstate from above)
00073       VetoedFinalState veto;
00074       veto.addVetoOnThisFinalState(elecs);
00075       veto.addVetoPairId(PID::MUON);
00076       veto.vetoNeutrinos();
00077       FastJets jets(veto, FastJets::ANTIKT, 0.4);
00078       addProjection(jets, "jets");
00079 
00080 
00081 
00082       /// book histograms
00083       _h_el_njet_inclusive = bookHisto1D(1,1,1);
00084       _h_mu_njet_inclusive = bookHisto1D(2,1,1);
00085 
00086       _h_el_pT_jet1 = bookHisto1D(5,1,1);
00087       _h_mu_pT_jet1 = bookHisto1D(6,1,1);
00088 
00089       _h_el_pT_jet2 = bookHisto1D(7,1,1);
00090       _h_mu_pT_jet2 = bookHisto1D(8,1,1);
00091     }
00092 
00093 
00094     /// Perform the per-event analysis
00095     void analyze(const Event& event) {
00096       const double weight = event.weight();
00097 
00098       const FinalState& elecs = applyProjection<FinalState>(event, "elecs");
00099       Particles elec_neutrino=applyProjection<FinalState>(event, "elec_neutrino").particles();
00100       if (elecs.size()==1 && elec_neutrino.size()>0) {
00101         FourMomentum lepton=elecs.particles()[0].momentum();
00102         foreach (const Particle& photon,
00103                  applyProjection<FinalState>(event, "cphotons_e").particles()) {
00104           lepton+=photon.momentum();
00105         }
00106         FourMomentum p_miss = elec_neutrino[0].momentum();
00107         double mT=sqrt(2.0*lepton.pT()*p_miss.Et()*(1.0-cos(lepton.phi()-p_miss.phi())));
00108         if (p_miss.Et()>25.0*GeV && mT>40.0*GeV) {
00109           Jets jets;
00110           foreach (const Jet& jet, applyProjection<FastJets>(event, "jets").jetsByPt(20.0*GeV)) {
00111             if (fabs(jet.eta())<2.8 && deltaR(lepton, jet.momentum())>0.5) {
00112               jets.push_back(jet);
00113             }
00114           }
00115 
00116           _h_el_njet_inclusive->fill(0, weight);
00117           if (jets.size()>=1) {
00118             _h_el_njet_inclusive->fill(1, weight);
00119             _h_el_pT_jet1->fill(jets[0].pT(), weight);
00120           }
00121           if (jets.size()>=2) {
00122             _h_el_njet_inclusive->fill(2, weight);
00123             _h_el_pT_jet2->fill(jets[1].pT(), weight);
00124           }
00125           if (jets.size()>=3) {
00126             _h_el_njet_inclusive->fill(3, weight);
00127           }
00128         }
00129       }
00130 
00131       const FinalState& muons = applyProjection<FinalState>(event, "muons");
00132       Particles muon_neutrino=applyProjection<FinalState>(event, "muon_neutrino").particles();
00133       if (muons.size()==1 && muon_neutrino.size()>0) {
00134         FourMomentum lepton=muons.particles()[0].momentum();
00135         foreach (const Particle& photon,
00136                  applyProjection<FinalState>(event, "cphotons_m").particles()) {
00137           lepton+=photon.momentum();
00138         }
00139         FourMomentum p_miss = muon_neutrino[0].momentum();
00140         double mT=sqrt(2.0*lepton.pT()*p_miss.Et()*(1.0-cos(lepton.phi()-p_miss.phi())));
00141         if (p_miss.Et()>25.0*GeV && mT>40.0*GeV) {
00142           Jets jets;
00143           foreach (const Jet& jet, applyProjection<FastJets>(event, "jets").jetsByPt(20.0*GeV)) {
00144             if (fabs(jet.eta())<2.8 && deltaR(lepton, jet.momentum())>0.5) {
00145               jets.push_back(jet);
00146             }
00147           }
00148 
00149           _h_mu_njet_inclusive->fill(0, weight);
00150           if (jets.size()>=1) {
00151             _h_mu_njet_inclusive->fill(1, weight);
00152             _h_mu_pT_jet1->fill(jets[0].pT(), weight);
00153           }
00154           if (jets.size()>=2) {
00155             _h_mu_njet_inclusive->fill(2, weight);
00156             _h_mu_pT_jet2->fill(jets[1].pT(), weight);
00157           }
00158           if (jets.size()>=3) {
00159             _h_mu_njet_inclusive->fill(3, weight);
00160           }
00161           if (jets.size()>=4) {
00162             _h_mu_njet_inclusive->fill(4, weight);
00163           }
00164         }
00165       }
00166 
00167     }
00168 
00169 
00170     /// Normalise histograms etc., after the run
00171     void finalize() {
00172       double normfac=crossSection()/sumOfWeights();
00173       scale(_h_el_njet_inclusive, normfac);
00174       scale(_h_mu_njet_inclusive, normfac);
00175       scale(_h_el_pT_jet1, normfac);
00176       scale(_h_mu_pT_jet1, normfac);
00177       scale(_h_el_pT_jet2, normfac);
00178       scale(_h_mu_pT_jet2, normfac);
00179     }
00180 
00181     //@}
00182 
00183 
00184   private:
00185 
00186     /// @name Histograms
00187     //@{
00188 
00189     Histo1DPtr _h_el_njet_inclusive;
00190     Histo1DPtr _h_mu_njet_inclusive;
00191     Histo1DPtr _h_el_pT_jet1;
00192     Histo1DPtr _h_mu_pT_jet1;
00193     Histo1DPtr _h_el_pT_jet2;
00194     Histo1DPtr _h_mu_pT_jet2;
00195     //@}
00196 
00197   };
00198 
00199 
00200 
00201   // The hook for the plugin system
00202   DECLARE_RIVET_PLUGIN(ATLAS_2010_S8919674);
00203 
00204 }