00001
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/RivetAIDA.hh"
00004 #include "Rivet/Projections/ChargedFinalState.hh"
00005 #include "Rivet/Projections/HadronicFinalState.hh"
00006 #include "Rivet/Tools/ParticleIdUtils.hh"
00007 #include "Rivet/Projections/Beam.hh"
00008 using namespace std;
00009
00010 namespace Rivet {
00011
00012
00013 class CMS_2011_S8884919 : public Analysis {
00014 public:
00015
00016 CMS_2011_S8884919() : Analysis("CMS_2011_S8884919") {
00017 setBeams(PROTON, PROTON);
00018 setNeedsCrossSection(false);
00019 }
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 private:
00037
00038
00039 vector<AIDA::IHistogram1D*> _h_dNch_dn;
00040
00041
00042
00043
00044 AIDA::IHistogram1D* _h_dNch_dn_pt500_eta24;
00045
00046
00047
00048
00049 AIDA::IProfile1D* _h_dmpt_dNch_eta24;
00050
00051 AIDA::IHistogram1D* _h_KNO_eta05;
00052 AIDA::IHistogram1D* _h_KNO_eta24;
00053
00054 double _Nevt_after_cuts;
00055 vector<double> _etabins;
00056 vector<int> _nch_in_Evt;
00057 vector<int> _nch_in_Evt_pt500;
00058 double _sumpt;
00059 int nch_max;
00060
00061
00062 void init();
00063 void analyze(const Event&);
00064 void finalize();
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076 };
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130 void CMS_2011_S8884919::init() {
00131 ChargedFinalState cfs(-2.4, 2.4, 0.0*GeV);
00132 addProjection(cfs, "CFS");
00133 addProjection(Beam(), "Beam");
00134
00135 nch_max = 400;
00136 _Nevt_after_cuts = 0;
00137
00138
00139 _etabins.push_back(0.5) ; _etabins.push_back(1.0) ; _etabins.push_back(1.5) ; _etabins.push_back(2.0) ; _etabins.push_back(2.4) ;
00140 _h_dNch_dn.reserve( _etabins.size() );
00141 ostringstream t("");
00142
00143 if(fuzzyEquals(sqrtS(), 900, 1E-3)){
00144
00145 for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){
00146 t.str("") ; t << "$|\\eta| <$ " << _etabins[ietabin] << " , $\\sqrt(s)$ = 0.9 TeV" ;
00147 _h_dNch_dn.push_back( bookHistogram1D( 2 + ietabin, 1, 1 , t.str() , "n" , "$P_{n}$") );
00148 }
00149
00150 _h_dNch_dn_pt500_eta24 = bookHistogram1D( 20 , 1, 1 , "$p_{T} >$ 500 GeV/c , $|\\eta| <$ 2.4 , $\\sqrt(s)$ = 0.9 TeV" , "n" , "$P_{n}$");
00151
00152
00153
00154
00155 _h_dmpt_dNch_eta24 = bookProfile1D( 23 , 1, 1 , "$|\\eta| <$ 2.4 , $\\sqrt(s)$ = 0.9 TeV" , "n" , "$< p_{T}> [ GeV/c ]$");
00156
00157 }
00158
00159 if(fuzzyEquals(sqrtS(), 2360, 1E-3)){
00160 for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){
00161 t.str("") ; t << "$|\\eta| <$ " << _etabins[ietabin] << " , $\\sqrt(s)$ = 2.36 TeV" ;
00162 _h_dNch_dn.push_back( bookHistogram1D( 7 + ietabin, 1, 1 , t.str() , "n" , "$P_{n}$") );
00163 }
00164
00165 _h_dNch_dn_pt500_eta24 = bookHistogram1D( 21 , 1, 1 , "$p_{T} >$ 500 GeV/c , $|\\eta| <$ 2.4 , $\\sqrt(s)$ = 2.36 TeV" , "n" , "$P_{n}$");
00166
00167
00168
00169
00170 _h_dmpt_dNch_eta24 = bookProfile1D( 24 , 1, 1 , "$|\\eta| <$ 2.4 , $\\sqrt(s)$ = 2.36 TeV" , "n" , "$< p_{T}> [ GeV/c ]$");
00171 }
00172
00173 if(fuzzyEquals(sqrtS(), 7000, 1E-3)){
00174 for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){
00175 t.str("") ; t << "$|\\eta| <$ " << _etabins[ietabin] << " , $\\sqrt(s)$ = 7 TeV" ;
00176 _h_dNch_dn.push_back( bookHistogram1D( 12 + ietabin, 1, 1 , t.str() , "n" , "$P_{n}$") );
00177 }
00178
00179 _h_dNch_dn_pt500_eta24 = bookHistogram1D( 22 , 1, 1 , "$p_{T} >$ 500 GeV/c , $|\\eta| <$ 2.4 , $\\sqrt(s)$ = 7 TeV" , "n" , "$P_{n}$");
00180
00181
00182
00183
00184 _h_dmpt_dNch_eta24 = bookProfile1D( 25 , 1, 1 , "$|\\eta| <$2.4 , $\\sqrt(s)$ = 7 TeV" , "n" , "$< p_{T}> [ GeV/c ]$");
00185
00186 }
00187
00188
00189
00190
00191 }
00192
00193
00194 void CMS_2011_S8884919::analyze(const Event& event) {
00195 const double weight = event.weight();
00196
00197
00198 const ChargedFinalState& charged = applyProjection<ChargedFinalState>(event, "CFS");
00199
00200
00201
00202
00203
00204
00205
00206 _Nevt_after_cuts += weight;
00207
00208
00209 _nch_in_Evt.assign(_etabins.size() , 0);
00210 _nch_in_Evt_pt500.assign(_etabins.size() , 0);
00211 _sumpt = 0;
00212
00213
00214
00215
00216 foreach (const Particle& p, charged.particles()) {
00217
00218
00219 if(! PID::isHadron(p.pdgId())) continue;
00220
00221 double pT = p.momentum().pT();
00222 double eta = p.momentum().eta();
00223
00224 _sumpt+=pT;
00225
00226
00227
00228 for (int ietabin=_etabins.size()-1; ietabin >= 0 ; --ietabin){
00229
00230 if (fabs(eta) <= _etabins[ietabin]){
00231 ++(_nch_in_Evt[ietabin]);
00232
00233 if(pT>0.5)
00234 ++(_nch_in_Evt_pt500[ietabin]);
00235 }
00236
00237 else
00238 break;
00239
00240 }
00241
00242 }
00243
00244
00245 for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){
00246 _h_dNch_dn[ietabin]->fill(_nch_in_Evt[ietabin] , weight);
00247 }
00248
00249
00250 if(_etabins[4] == 2.4 && _etabins[0] == 0.5){
00251 if(_nch_in_Evt[4] != 0) _h_dmpt_dNch_eta24->fill(_nch_in_Evt[4] , _sumpt / _nch_in_Evt[4] , weight);
00252
00253 _h_dNch_dn_pt500_eta24->fill(_nch_in_Evt_pt500[4] , weight);
00254
00255
00256
00257
00258
00259
00260 }
00261 else
00262 getLog() << Log::WARNING << "You changed the number of eta bins, but forgot to propagate it everywhere !! " << endl;
00263
00264 }
00265
00266
00267 void CMS_2011_S8884919::finalize() {
00268 getLog() << Log::INFO << "Number of events after event selection: " << _Nevt_after_cuts << endl;
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278 for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){
00279 normalize(_h_dNch_dn[ietabin]);
00280 }
00281 normalize(_h_dNch_dn_pt500_eta24);
00282 }
00283
00284
00285
00286 AnalysisBuilder<CMS_2011_S8884919> plugin_CMS_2011_S8884919;
00287
00288 }
00289