CMS_2011_S8884919.cc

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