CMS_2011_S8884919.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/RivetYODA.hh" 00004 #include "Rivet/Projections/ChargedFinalState.hh" 00005 #include "Rivet/Tools/ParticleIdUtils.hh" 00006 #include "Rivet/Projections/Beam.hh" 00007 using namespace std; 00008 00009 namespace Rivet { 00010 00011 class CMS_2011_S8884919 : public Analysis { 00012 public: 00013 00014 CMS_2011_S8884919() 00015 : Analysis("CMS_2011_S8884919") 00016 { } 00017 00018 00019 void init() { 00020 ChargedFinalState cfs(-2.4, 2.4, 0.0*GeV); 00021 addProjection(cfs, "CFS"); 00022 00023 // eta bins 00024 _etabins.push_back(0.5); 00025 _etabins.push_back(1.0); 00026 _etabins.push_back(1.5); 00027 _etabins.push_back(2.0); 00028 _etabins.push_back(2.4) ; 00029 00030 if (fuzzyEquals(sqrtS()/GeV, 900)) { 00031 for (size_t ietabin=0; ietabin < _etabins.size(); ietabin++) { 00032 _h_dNch_dn.push_back( bookHisto1D( 2 + ietabin, 1, 1) ); 00033 } 00034 _h_dNch_dn_pt500_eta24 = bookHisto1D(20, 1, 1); 00035 _h_dmpt_dNch_eta24 = bookProfile1D(23, 1, 1); 00036 } 00037 00038 if (fuzzyEquals(sqrtS()/GeV, 2360)) { 00039 for (size_t ietabin=0; ietabin < _etabins.size(); ietabin++) { 00040 _h_dNch_dn.push_back( bookHisto1D(7 + ietabin, 1, 1) ); 00041 } 00042 _h_dNch_dn_pt500_eta24 = bookHisto1D(21, 1, 1); 00043 _h_dmpt_dNch_eta24 = bookProfile1D(24, 1, 1); 00044 } 00045 00046 if (fuzzyEquals(sqrtS()/GeV, 7000)) { 00047 for (size_t ietabin=0; ietabin < _etabins.size(); ietabin++) { 00048 _h_dNch_dn.push_back( bookHisto1D(12 + ietabin, 1, 1) ); 00049 } 00050 _h_dNch_dn_pt500_eta24 = bookHisto1D(22, 1, 1); 00051 _h_dmpt_dNch_eta24 = bookProfile1D(25, 1, 1); 00052 } 00053 } 00054 00055 00056 void analyze(const Event& event) { 00057 const double weight = event.weight(); 00058 00059 // Get the charged particles 00060 const ChargedFinalState& charged = applyProjection<ChargedFinalState>(event, "CFS"); 00061 00062 // Resetting the multiplicity for the event to 0; 00063 vector<int> _nch_in_Evt; 00064 vector<int> _nch_in_Evt_pt500; 00065 _nch_in_Evt.assign(_etabins.size(), 0); 00066 _nch_in_Evt_pt500.assign(_etabins.size(), 0); 00067 double sumpt = 0; 00068 00069 // Loop over particles in event 00070 foreach (const Particle& p, charged.particles()) { 00071 // Selecting only charged hadrons 00072 if (! PID::isHadron(p.pdgId())) continue; 00073 00074 double pT = p.momentum().pT(); 00075 double eta = p.momentum().eta(); 00076 sumpt += pT; 00077 for (int ietabin = _etabins.size()-1; ietabin >= 0; --ietabin) { 00078 if (fabs(eta) <= _etabins[ietabin]){ 00079 ++(_nch_in_Evt[ietabin]); 00080 if (pT > 0.5/GeV) ++(_nch_in_Evt_pt500[ietabin]); 00081 } 00082 else 00083 break; 00084 } 00085 } 00086 00087 // Filling multiplicity-dependent histogramms 00088 for (size_t ietabin = 0; ietabin < _etabins.size(); ietabin++) { 00089 _h_dNch_dn[ietabin]->fill(_nch_in_Evt[ietabin], weight); 00090 } 00091 00092 // Do only if eta bins are the needed ones 00093 if (_etabins[4] == 2.4 && _etabins[0] == 0.5) { 00094 if (_nch_in_Evt[4] != 0) { 00095 _h_dmpt_dNch_eta24->fill(_nch_in_Evt[4], sumpt/GeV / _nch_in_Evt[4], weight); 00096 } 00097 _h_dNch_dn_pt500_eta24->fill(_nch_in_Evt_pt500[4], weight); 00098 } else { 00099 MSG_WARNING("You changed the number of eta bins, but forgot to propagate it everywhere !!"); 00100 } 00101 } 00102 00103 00104 void finalize() { 00105 for (size_t ietabin = 0; ietabin < _etabins.size(); ietabin++){ 00106 normalize(_h_dNch_dn[ietabin]); 00107 } 00108 normalize(_h_dNch_dn_pt500_eta24); 00109 } 00110 00111 00112 private: 00113 00114 vector<Histo1DPtr> _h_dNch_dn; 00115 Histo1DPtr _h_dNch_dn_pt500_eta24; 00116 Profile1DPtr _h_dmpt_dNch_eta24; 00117 00118 vector<double> _etabins; 00119 00120 }; 00121 00122 00123 // The hook for the plugin system 00124 DECLARE_RIVET_PLUGIN(CMS_2011_S8884919); 00125 00126 } Generated on Fri Dec 21 2012 15:03:40 for The Rivet MC analysis system by ![]() |