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 00014 00015 /// Underlying event activity in the Drell-Yan process at 7 TeV 00016 class CMS_2012_I1107658 : public Analysis { 00017 public: 00018 00019 /// Constructor 00020 CMS_2012_I1107658() 00021 : Analysis("CMS_2012_I1107658") 00022 { } 00023 00024 00025 /// Initialization 00026 void init() { 00027 00028 /// @note Using a bare muon Z (but with a clustering radius!?) 00029 Cut cut = Cuts::abseta < 2.4 && Cuts::pT > 20*GeV; 00030 ZFinder zfinder(FinalState(), cut, PID::MUON, 4*GeV, 140*GeV, 0.2, ZFinder::NOCLUSTER); 00031 addProjection(zfinder, "ZFinder"); 00032 00033 ChargedFinalState cfs(-2, 2, 500*MeV); 00034 VetoedFinalState nonmuons(cfs); 00035 nonmuons.addVetoPairId(PID::MUON); 00036 addProjection(nonmuons, "nonmuons"); 00037 00038 _h_Nchg_towards_pTmumu = bookProfile1D(1, 1, 1); 00039 _h_Nchg_transverse_pTmumu = bookProfile1D(2, 1, 1); 00040 _h_Nchg_away_pTmumu = bookProfile1D(3, 1, 1); 00041 _h_pTsum_towards_pTmumu = bookProfile1D(4, 1, 1); 00042 _h_pTsum_transverse_pTmumu = bookProfile1D(5, 1, 1); 00043 _h_pTsum_away_pTmumu = bookProfile1D(6, 1, 1); 00044 _h_avgpT_towards_pTmumu = bookProfile1D(7, 1, 1); 00045 _h_avgpT_transverse_pTmumu = bookProfile1D(8, 1, 1); 00046 _h_avgpT_away_pTmumu = bookProfile1D(9, 1, 1); 00047 _h_Nchg_towards_plus_transverse_Mmumu = bookProfile1D(10, 1, 1); 00048 _h_pTsum_towards_plus_transverse_Mmumu = bookProfile1D(11, 1, 1); 00049 _h_avgpT_towards_plus_transverse_Mmumu = bookProfile1D(12, 1, 1); 00050 _h_Nchg_towards_zmass_81_101 = bookHisto1D(13, 1, 1); 00051 _h_Nchg_transverse_zmass_81_101 = bookHisto1D(14, 1, 1); 00052 _h_Nchg_away_zmass_81_101 = bookHisto1D(15, 1, 1); 00053 _h_pT_towards_zmass_81_101 = bookHisto1D(16, 1, 1); 00054 _h_pT_transverse_zmass_81_101 = bookHisto1D(17, 1, 1); 00055 _h_pT_away_zmass_81_101 = bookHisto1D(18, 1, 1); 00056 _h_Nchg_transverse_zpt_5 = bookHisto1D(19, 1, 1); 00057 _h_pT_transverse_zpt_5 = bookHisto1D(20, 1, 1); 00058 } 00059 00060 00061 /// Perform the per-event analysis 00062 void analyze(const Event& event) { 00063 const double weight = event.weight(); 00064 const ZFinder& zfinder = applyProjection<ZFinder>(event, "ZFinder"); 00065 00066 if (zfinder.bosons().size() != 1) vetoEvent; 00067 00068 double Zpt = zfinder.bosons()[0].pT()/GeV; 00069 double Zphi = zfinder.bosons()[0].phi(); 00070 double Zmass = zfinder.bosons()[0].mass()/GeV; 00071 00072 Particles particles = applyProjection<VetoedFinalState>(event, "nonmuons").particles(); 00073 00074 int nTowards = 0; 00075 int nTransverse = 0; 00076 int nAway = 0; 00077 double ptSumTowards = 0; 00078 double ptSumTransverse = 0; 00079 double ptSumAway = 0; 00080 00081 foreach (const Particle& p, particles) { 00082 double dphi = fabs(deltaPhi(Zphi, p.phi())); 00083 double pT = p.pT(); 00084 00085 if ( dphi < M_PI/3 ) { 00086 nTowards++; 00087 ptSumTowards += pT; 00088 if (Zmass > 81. && Zmass < 101.) _h_pT_towards_zmass_81_101->fill(pT, weight); 00089 } else if ( dphi < 2.*M_PI/3 ) { 00090 nTransverse++; 00091 ptSumTransverse += pT; 00092 if (Zmass > 81. && Zmass < 101.) _h_pT_transverse_zmass_81_101->fill(pT, weight); 00093 if (Zpt < 5.) _h_pT_transverse_zpt_5->fill(pT, weight); 00094 } else { 00095 nAway++; 00096 ptSumAway += pT; 00097 if (Zmass > 81. && Zmass < 101.) _h_pT_away_zmass_81_101->fill(pT, weight); 00098 } 00099 00100 } // Loop over particles 00101 00102 00103 const double area = 8./3.*M_PI; 00104 if (Zmass > 81. && Zmass < 101.) { 00105 _h_Nchg_towards_pTmumu-> fill(Zpt, 1./area * nTowards, weight); 00106 _h_Nchg_transverse_pTmumu-> fill(Zpt, 1./area * nTransverse, weight); 00107 _h_Nchg_away_pTmumu-> fill(Zpt, 1./area * nAway, weight); 00108 _h_pTsum_towards_pTmumu-> fill(Zpt, 1./area * ptSumTowards, weight); 00109 _h_pTsum_transverse_pTmumu-> fill(Zpt, 1./area * ptSumTransverse, weight); 00110 _h_pTsum_away_pTmumu-> fill(Zpt, 1./area * ptSumAway, weight); 00111 if (nTowards > 0) _h_avgpT_towards_pTmumu-> fill(Zpt, ptSumTowards/nTowards, weight); 00112 if (nTransverse > 0) _h_avgpT_transverse_pTmumu-> fill(Zpt, ptSumTransverse/nTransverse, weight); 00113 if (nAway > 0) _h_avgpT_away_pTmumu-> fill(Zpt, ptSumAway/nAway, weight); 00114 _h_Nchg_towards_zmass_81_101-> fill(nTowards, weight); 00115 _h_Nchg_transverse_zmass_81_101->fill(nTransverse, weight); 00116 _h_Nchg_away_zmass_81_101-> fill(nAway, weight); 00117 } 00118 00119 if (Zpt < 5.) { 00120 _h_Nchg_towards_plus_transverse_Mmumu->fill(Zmass, (nTowards + nTransverse)/(2.*area), weight); 00121 _h_pTsum_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(2.*area), weight); 00122 if ((nTowards + nTransverse) > 0) _h_avgpT_towards_plus_transverse_Mmumu->fill(Zmass, (ptSumTowards + ptSumTransverse)/(nTowards + nTransverse), weight); 00123 _h_Nchg_transverse_zpt_5->fill(nTransverse, weight); 00124 } 00125 00126 } 00127 00128 00129 /// Normalise histograms etc., after the run 00130 void finalize() { 00131 scale(_h_pT_towards_zmass_81_101, safediv(1, _h_Nchg_towards_zmass_81_101->integral(), 0)); 00132 scale(_h_pT_transverse_zmass_81_101, safediv(1, _h_Nchg_transverse_zmass_81_101->integral(), 0)); 00133 scale(_h_pT_away_zmass_81_101, safediv(1, _h_Nchg_away_zmass_81_101->integral(), 0)); 00134 scale(_h_pT_transverse_zpt_5, safediv(1, _h_Nchg_transverse_zpt_5->integral(), 0)); 00135 normalize(_h_Nchg_towards_zmass_81_101); 00136 normalize(_h_Nchg_transverse_zmass_81_101); 00137 normalize(_h_Nchg_away_zmass_81_101); 00138 normalize(_h_Nchg_transverse_zpt_5); 00139 } 00140 00141 00142 private: 00143 00144 00145 /// @name Histogram objects 00146 //@{ 00147 Profile1DPtr _h_Nchg_towards_pTmumu; 00148 Profile1DPtr _h_Nchg_transverse_pTmumu; 00149 Profile1DPtr _h_Nchg_away_pTmumu; 00150 Profile1DPtr _h_pTsum_towards_pTmumu; 00151 Profile1DPtr _h_pTsum_transverse_pTmumu; 00152 Profile1DPtr _h_pTsum_away_pTmumu; 00153 Profile1DPtr _h_avgpT_towards_pTmumu; 00154 Profile1DPtr _h_avgpT_transverse_pTmumu; 00155 Profile1DPtr _h_avgpT_away_pTmumu; 00156 Profile1DPtr _h_Nchg_towards_plus_transverse_Mmumu; 00157 Profile1DPtr _h_pTsum_towards_plus_transverse_Mmumu; 00158 Profile1DPtr _h_avgpT_towards_plus_transverse_Mmumu; 00159 Histo1DPtr _h_Nchg_towards_zmass_81_101; 00160 Histo1DPtr _h_Nchg_transverse_zmass_81_101; 00161 Histo1DPtr _h_Nchg_away_zmass_81_101; 00162 Histo1DPtr _h_pT_towards_zmass_81_101; 00163 Histo1DPtr _h_pT_transverse_zmass_81_101; 00164 Histo1DPtr _h_pT_away_zmass_81_101; 00165 Histo1DPtr _h_Nchg_transverse_zpt_5; 00166 Histo1DPtr _h_pT_transverse_zpt_5; 00167 //@} 00168 00169 }; 00170 00171 00172 // Hook for the plugin system 00173 DECLARE_RIVET_PLUGIN(CMS_2012_I1107658); 00174 00175 } Generated on Wed Oct 7 2015 12:09:12 for The Rivet MC analysis system by ![]() |