CMS_2010_S8656010.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 00007 namespace Rivet { 00008 00009 00010 class CMS_2010_S8656010 : public Analysis { 00011 public: 00012 00013 CMS_2010_S8656010() : Analysis("CMS_2010_S8656010") {} 00014 00015 00016 void init() { 00017 ChargedFinalState cfs(-2.5, 2.5, 0.0*GeV); 00018 addProjection(cfs, "CFS"); 00019 00020 for (int d=1; d<=3; d++) { 00021 for (int y=1; y<=4; y++) { 00022 _h_dNch_dpT.push_back(bookHisto1D(d, 1, y)); 00023 } 00024 } 00025 00026 _h_dNch_dpT_all = bookHisto1D(4, 1, 1); 00027 _h_dNch_dEta = bookHisto1D(5, 1, 1); 00028 } 00029 00030 00031 void analyze(const Event& event) { 00032 const double weight = event.weight(); 00033 00034 //charged particles 00035 const ChargedFinalState& charged = applyProjection<ChargedFinalState>(event, "CFS"); 00036 00037 foreach (const Particle& p, charged.particles()) { 00038 //selecting only charged hadrons 00039 if (! PID::isHadron(p.pdgId())) continue; 00040 00041 const double pT = p.momentum().pT(); 00042 const double eta = p.momentum().eta(); 00043 00044 // The data is actually a duplicated folded distribution. This should mimic it. 00045 _h_dNch_dEta->fill(eta, 0.5*weight); 00046 _h_dNch_dEta->fill(-eta, 0.5*weight); 00047 if (fabs(eta) < 2.4 && pT > 0.1*GeV) { 00048 if (pT < 6.0*GeV) { 00049 _h_dNch_dpT_all->fill(pT/GeV, weight/(pT/GeV)); 00050 if (pT < 2.0*GeV) { 00051 int ietabin = int(fabs(eta)/0.2); 00052 _h_dNch_dpT[ietabin]->fill(pT/GeV, weight); 00053 } 00054 } 00055 } 00056 } 00057 } 00058 00059 00060 void finalize() { 00061 const double normfac = 1.0/sumOfWeights(); // Normalizing to unit eta is automatic 00062 // The pT distributions in bins of eta must be normalized to unit eta. This is a factor of 2 00063 // for the |eta| times 0.2 (eta range). 00064 // The pT distributions over all eta are normalized to unit eta (2.0*2.4) and by 1/2*pi*pT. 00065 // The 1/pT part is taken care of in the filling. The 1/2pi is taken care of here. 00066 const double normpT = normfac/(2.0*0.2); 00067 const double normpTall = normfac/(2.0*M_PI*2.0*2.4); 00068 00069 for (size_t ietabin=0; ietabin < _h_dNch_dpT.size(); ietabin++){ 00070 scale(_h_dNch_dpT[ietabin], normpT); 00071 } 00072 scale(_h_dNch_dpT_all, normpTall); 00073 scale(_h_dNch_dEta, normfac); 00074 } 00075 00076 00077 private: 00078 00079 std::vector<Histo1DPtr> _h_dNch_dpT; 00080 Histo1DPtr _h_dNch_dpT_all; 00081 Histo1DPtr _h_dNch_dEta; 00082 00083 }; 00084 00085 00086 00087 // The hook for the plugin system 00088 DECLARE_RIVET_PLUGIN(CMS_2010_S8656010); 00089 00090 } Generated on Fri Dec 21 2012 15:03:40 for The Rivet MC analysis system by ![]() |