CMS_2012_I1107658.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/FinalState.hh" 00004 #include "Rivet/Projections/ZFinder.hh" 00005 #include "Rivet/Projections/FinalState.hh" 00006 #include "Rivet/Projections/ChargedFinalState.hh" 00007 #include "Rivet/Projections/VetoedFinalState.hh" 00008 #include "Rivet/ParticleName.hh" 00009 00010 namespace Rivet { 00011 00012 00013 class CMS_2012_I1107658 : public Analysis { 00014 public: 00015 00016 /// Constructor 00017 CMS_2012_I1107658() : Analysis("CMS_2012_I1107658") {} 00018 00019 void init() { 00020 00021 FinalState fs; 00022 ZFinder zfinder(fs, -2.4, 2.4, 20.0*GeV, PID::MUON, 4.0*GeV, 140.0*GeV, 0.2, false, false); 00023 addProjection(zfinder, "ZFinder"); 00024 00025 ChargedFinalState cfs(-2.0, 2.0, 500*MeV); // For charged particles 00026 VetoedFinalState nonmuons(cfs); 00027 nonmuons.addVetoPairId(PID::MUON); 00028 addProjection(nonmuons, "nonmuons"); 00029 00030 _h_Nchg_towards_pTmumu = bookProfile1D(1, 1, 1); 00031 _h_Nchg_transverse_pTmumu = bookProfile1D(2, 1, 1); 00032 _h_Nchg_away_pTmumu = bookProfile1D(3, 1, 1); 00033 _h_pTsum_towards_pTmumu = bookProfile1D(4, 1, 1); 00034 _h_pTsum_transverse_pTmumu = bookProfile1D(5, 1, 1); 00035 _h_pTsum_away_pTmumu = bookProfile1D(6, 1, 1); 00036 _h_avgpT_towards_pTmumu = bookProfile1D(7, 1, 1); 00037 _h_avgpT_transverse_pTmumu = bookProfile1D(8, 1, 1); 00038 _h_avgpT_away_pTmumu = bookProfile1D(9, 1, 1); 00039 _h_Nchg_towards_plus_transverse_Mmumu = bookProfile1D(10, 1, 1); 00040 _h_pTsum_towards_plus_transverse_Mmumu = bookProfile1D(11, 1, 1); 00041 _h_avgpT_towards_plus_transverse_Mmumu = bookProfile1D(12, 1, 1); 00042 _h_Nchg_towards_zmass_81_101 = bookHisto1D(13, 1, 1); 00043 _h_Nchg_transverse_zmass_81_101 = bookHisto1D(14, 1, 1); 00044 _h_Nchg_away_zmass_81_101 = bookHisto1D(15, 1, 1); 00045 _h_pT_towards_zmass_81_101 = bookHisto1D(16, 1, 1); 00046 _h_pT_transverse_zmass_81_101 = bookHisto1D(17, 1, 1); 00047 _h_pT_away_zmass_81_101 = bookHisto1D(18, 1, 1); 00048 _h_Nchg_transverse_zpt_5 = bookHisto1D(19, 1, 1); 00049 _h_pT_transverse_zpt_5 = bookHisto1D(20, 1, 1); 00050 } 00051 00052 00053 /// Perform the per-event analysis 00054 void analyze(const Event& event) { 00055 const double weight = event.weight(); 00056 const ZFinder& zfinder = applyProjection<ZFinder>(event, "ZFinder"); 00057 00058 if (zfinder.bosons().size() != 1) vetoEvent; 00059 00060 double Zpt = zfinder.bosons()[0].pT()/GeV; 00061 double Zphi = zfinder.bosons()[0].momentum().phi(); 00062 double Zmass = zfinder.bosons()[0].momentum().mass()/GeV; 00063 00064 Particles particles = applyProjection<VetoedFinalState>(event, "nonmuons").particles(); 00065 00066 int nTowards = 0; 00067 int nTransverse = 0; 00068 int nAway = 0; 00069 double ptSumTowards = 0.0; 00070 double ptSumTransverse = 0.0; 00071 double ptSumAway = 0.0; 00072 00073 foreach (const Particle& p, particles) { 00074 double dphi = fabs(deltaPhi(Zphi, p.momentum().phi())); 00075 double pT = p.pT(); 00076 00077 if ( dphi < M_PI/3.0 ) { 00078 nTowards++; 00079 ptSumTowards += pT; 00080 if (Zmass > 81. && Zmass < 101.) _h_pT_towards_zmass_81_101->fill(pT, weight); 00081 } else if ( dphi < 2.*M_PI/3.0 ) { 00082 nTransverse++; 00083 ptSumTransverse += pT; 00084 if (Zmass > 81. && Zmass < 101.) _h_pT_transverse_zmass_81_101->fill(pT, weight); 00085 if (Zpt < 5.) _h_pT_transverse_zpt_5->fill(pT, weight); 00086 } else { 00087 nAway++; 00088 ptSumAway += pT; 00089 if (Zmass > 81. && Zmass < 101.) _h_pT_away_zmass_81_101->fill(pT, weight); 00090 } 00091 00092 } // Loop over particles 00093 00094 00095 const double area = 8./3.*M_PI; 00096 if (Zmass > 81. && Zmass < 101.) { 00097 _h_Nchg_towards_pTmumu-> fill(Zpt, 1./area * nTowards, weight); 00098 _h_Nchg_transverse_pTmumu-> fill(Zpt, 1./area * nTransverse, weight); 00099 _h_Nchg_away_pTmumu-> fill(Zpt, 1./area * nAway, weight); 00100 _h_pTsum_towards_pTmumu-> fill(Zpt, 1./area * ptSumTowards, weight); 00101 _h_pTsum_transverse_pTmumu-> fill(Zpt, 1./area * ptSumTransverse, weight); 00102 _h_pTsum_away_pTmumu-> fill(Zpt, 1./area * ptSumAway, weight); 00103 if (nTowards > 0) _h_avgpT_towards_pTmumu-> fill(Zpt, ptSumTowards/nTowards, weight); 00104 if (nTransverse > 0) _h_avgpT_transverse_pTmumu-> fill(Zpt, ptSumTransverse/nTransverse, weight); 00105 if (nAway > 0) _h_avgpT_away_pTmumu-> fill(Zpt, ptSumAway/nAway, weight); 00106 _h_Nchg_towards_zmass_81_101-> fill(nTowards, weight); 00107 _h_Nchg_transverse_zmass_81_101->fill(nTransverse, weight); 00108 _h_Nchg_away_zmass_81_101-> fill(nAway, weight); 00109 } 00110 00111 if (Zpt < 5.) { 00112 _h_Nchg_towards_plus_transverse_Mmumu->fill(Zmass, (nTowards + nTransverse)/(2.*area), weight); 00113 _h_pTsum_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(2.*area), weight); 00114 if ((nTowards + nTransverse) > 0) _h_avgpT_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(nTowards + nTransverse), weight); 00115 _h_Nchg_transverse_zpt_5->fill(nTransverse, weight); 00116 } 00117 00118 } 00119 00120 00121 /// Normalise histograms etc., after the run 00122 void finalize() { 00123 if (integral(_h_Nchg_towards_zmass_81_101) > 0) scale(_h_pT_towards_zmass_81_101, 1.0/integral(_h_Nchg_towards_zmass_81_101)); 00124 if (integral(_h_Nchg_transverse_zmass_81_101) > 0) scale(_h_pT_transverse_zmass_81_101, 1.0/integral(_h_Nchg_transverse_zmass_81_101)); 00125 if (integral(_h_Nchg_away_zmass_81_101) > 0) scale(_h_pT_away_zmass_81_101, 1.0/integral(_h_Nchg_away_zmass_81_101)); 00126 if (integral(_h_Nchg_transverse_zpt_5) > 0) scale(_h_pT_transverse_zpt_5, 1.0/integral(_h_Nchg_transverse_zpt_5)); 00127 normalize(_h_Nchg_towards_zmass_81_101); 00128 normalize(_h_Nchg_transverse_zmass_81_101); 00129 normalize(_h_Nchg_away_zmass_81_101); 00130 normalize(_h_Nchg_transverse_zpt_5); 00131 } 00132 00133 00134 private: 00135 00136 00137 /// @name Histogram objects 00138 //@{ 00139 Profile1DPtr _h_Nchg_towards_pTmumu; 00140 Profile1DPtr _h_Nchg_transverse_pTmumu; 00141 Profile1DPtr _h_Nchg_away_pTmumu; 00142 Profile1DPtr _h_pTsum_towards_pTmumu; 00143 Profile1DPtr _h_pTsum_transverse_pTmumu; 00144 Profile1DPtr _h_pTsum_away_pTmumu; 00145 Profile1DPtr _h_avgpT_towards_pTmumu; 00146 Profile1DPtr _h_avgpT_transverse_pTmumu; 00147 Profile1DPtr _h_avgpT_away_pTmumu; 00148 Profile1DPtr _h_Nchg_towards_plus_transverse_Mmumu; 00149 Profile1DPtr _h_pTsum_towards_plus_transverse_Mmumu; 00150 Profile1DPtr _h_avgpT_towards_plus_transverse_Mmumu; 00151 Histo1DPtr _h_Nchg_towards_zmass_81_101; 00152 Histo1DPtr _h_Nchg_transverse_zmass_81_101; 00153 Histo1DPtr _h_Nchg_away_zmass_81_101; 00154 Histo1DPtr _h_pT_towards_zmass_81_101; 00155 Histo1DPtr _h_pT_transverse_zmass_81_101; 00156 Histo1DPtr _h_pT_away_zmass_81_101; 00157 Histo1DPtr _h_Nchg_transverse_zpt_5; 00158 Histo1DPtr _h_pT_transverse_zpt_5; 00159 //@} 00160 00161 }; 00162 00163 // This global object acts as a hook for the plugin system 00164 DECLARE_RIVET_PLUGIN(CMS_2012_I1107658); 00165 00166 } Generated on Fri Oct 25 2013 12:41:45 for The Rivet MC analysis system by ![]() |