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