rivet is hosted by Hepforge, IPPP Durham
ATLAS_2012_I1082009.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 #include "Rivet/Projections/UnstableFinalState.hh"
00010 
00011 namespace Rivet {
00012 
00013 
00014   class ATLAS_2012_I1082009 : public Analysis {
00015   public:
00016 
00017     /// @name Constructors etc.
00018     //@{
00019 
00020     /// Constructor
00021     ATLAS_2012_I1082009()
00022       : Analysis("ATLAS_2012_I1082009"),
00023         _weight25_30(0.),_weight30_40(0.),_weight40_50(0.),
00024         _weight50_60(0.),_weight60_70(0.),_weight25_70(0.)
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       // Input for the jets: No neutrinos, no muons
00039       VetoedFinalState veto;
00040       veto.addVetoPairId(PID::MUON);
00041       veto.vetoNeutrinos();
00042       FastJets jets(veto, FastJets::ANTIKT, 0.6);
00043       addProjection(jets, "jets");
00044       // unstable final-state for D*
00045       addProjection(UnstableFinalState(), "UFS");
00046 
00047       _h_pt25_30 = bookHisto1D( 8,1,1);
00048       _h_pt30_40 = bookHisto1D( 9,1,1);
00049       _h_pt40_50 = bookHisto1D(10,1,1);
00050       _h_pt50_60 = bookHisto1D(11,1,1);
00051       _h_pt60_70 = bookHisto1D(12,1,1);
00052       _h_pt25_70 = bookHisto1D(13,1,1);
00053     }
00054 
00055 
00056     /// Perform the per-event analysis
00057     void analyze(const Event& event) {
00058       const double weight = event.weight();
00059 
00060       // get the jets
00061       Jets jets;
00062       foreach (const Jet& jet, applyProjection<FastJets>(event, "jets").jetsByPt(25.0*GeV)) {
00063         if ( fabs(jet.eta()) < 2.5 ) jets.push_back(jet);
00064       }
00065       // get the D* mesons
00066       const UnstableFinalState& ufs = applyProjection<UnstableFinalState>(event, "UFS");
00067       Particles Dstar;
00068       foreach (const Particle& p, ufs.particles()) {
00069         const int id = abs(p.pdgId());
00070         if(id==413) Dstar.push_back(p);
00071       }
00072 
00073       // loop over the jobs
00074       foreach (const Jet& jet, jets ) {
00075         double perp = jet.momentum().perp();
00076         bool found = false;
00077         double z(0.);
00078         if(perp<25.||perp>70.) continue;
00079         foreach(const Particle & p, Dstar) {
00080           if(p.momentum().perp()<7.5) continue;
00081           if(deltaR(p, jet.momentum())<0.6) {
00082             Vector3 axis = jet.momentum().vector3().unit();
00083             z = axis.dot(p.momentum().vector3())/jet.momentum().E();
00084             if(z<0.3) continue;
00085             found = true;
00086             break;
00087           }
00088         }
00089         _weight25_70 += weight;
00090         if(found) _h_pt25_70->fill(z,weight);
00091         if(perp>=25.&&perp<30.) {
00092           _weight25_30 += weight;
00093           if(found) _h_pt25_30->fill(z,weight);
00094         }
00095         else if(perp>=30.&&perp<40.) {
00096           _weight30_40 += weight;
00097           if(found) _h_pt30_40->fill(z,weight);
00098         }
00099         else if(perp>=40.&&perp<50.) {
00100           _weight40_50 += weight;
00101           if(found) _h_pt40_50->fill(z,weight);
00102         }
00103         else if(perp>=50.&&perp<60.) {
00104           _weight50_60 += weight;
00105           if(found) _h_pt50_60->fill(z,weight);
00106         }
00107         else if(perp>=60.&&perp<70.) {
00108           _weight60_70 += weight;
00109           if(found) _h_pt60_70->fill(z,weight);
00110         }
00111       }
00112     }
00113 
00114 
00115     /// Normalise histograms etc., after the run
00116     void finalize() {
00117       scale(_h_pt25_30,1./_weight25_30);
00118       scale(_h_pt30_40,1./_weight30_40);
00119       scale(_h_pt40_50,1./_weight40_50);
00120       scale(_h_pt50_60,1./_weight50_60);
00121       scale(_h_pt60_70,1./_weight60_70);
00122       scale(_h_pt25_70,1./_weight25_70);
00123     }
00124 
00125     //@}
00126 
00127 
00128   private:
00129 
00130     /// @name Histograms
00131     //@{
00132     double _weight25_30,_weight30_40,_weight40_50;
00133     double _weight50_60,_weight60_70,_weight25_70;
00134 
00135     Histo1DPtr _h_pt25_30;
00136     Histo1DPtr _h_pt30_40;
00137     Histo1DPtr _h_pt40_50;
00138     Histo1DPtr _h_pt50_60;
00139     Histo1DPtr _h_pt60_70;
00140     Histo1DPtr _h_pt25_70;
00141     //@}
00142 
00143   };
00144 
00145   // The hook for the plugin system
00146   DECLARE_RIVET_PLUGIN(ATLAS_2012_I1082009);
00147 
00148 }