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