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