ATLAS_2010_S8919674.cc

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