Public Member Functions | |
CMS_2011_S8884919 () | |
Private Member Functions | |
void | analyze (const Event &) |
void | finalize () |
void | init () |
Private Attributes | |
vector< double > | _etabins |
AIDA::IProfile1D * | _h_dmpt_dNch_eta24 |
vector< AIDA::IHistogram1D * > | _h_dNch_dn |
AIDA::IHistogram1D * | _h_dNch_dn_pt500_eta24 |
AIDA::IHistogram1D * | _h_KNO_eta05 |
AIDA::IHistogram1D * | _h_KNO_eta24 |
vector< int > | _nch_in_Evt |
vector< int > | _nch_in_Evt_pt500 |
double | _Nevt_after_cuts |
double | _sumpt |
int | nch_max |
Definition at line 13 of file CMS_2011_S8884919.cc.
Rivet::CMS_2011_S8884919::CMS_2011_S8884919 | ( | ) | [inline] |
Definition at line 16 of file CMS_2011_S8884919.cc.
: Analysis("CMS_2011_S8884919") { setBeams(PROTON, PROTON); setNeedsCrossSection(false); }
void Rivet::CMS_2011_S8884919::analyze | ( | const Event & | event | ) | [private] |
Definition at line 194 of file CMS_2011_S8884919.cc.
References eta(), event(), AlCaHLTBitMon_ParallelJobs::p, dqm::qstatus::WARNING, and CommonMethods::weight().
{ const double weight = event.weight(); //charge particles const ChargedFinalState& charged = applyProjection<ChargedFinalState>(event, "CFS"); //This cut is not needed /*if (charged.particles().size()<1) { vetoEvent; }*/ _Nevt_after_cuts += weight; //resetting the multiplicity for the event to 0; _nch_in_Evt.assign(_etabins.size() , 0); _nch_in_Evt_pt500.assign(_etabins.size() , 0); _sumpt = 0; //std::cout << charged.size() << std::endl; //Loop over particles in event foreach (const Particle& p, charged.particles()) { //selecting only charged hadrons if(! PID::isHadron(p.pdgId())) continue; double pT = p.momentum().pT(); double eta = p.momentum().eta(); _sumpt+=pT; //cout << "eta : " << eta << " pT : " << pT << endl; for (int ietabin=_etabins.size()-1; ietabin >= 0 ; --ietabin){ //cout << " etabin : " << _etabins[ietabin] << " res : " << (fabs(eta) <= _etabins[ietabin]) << endl; if (fabs(eta) <= _etabins[ietabin]){ ++(_nch_in_Evt[ietabin]); if(pT>0.5) ++(_nch_in_Evt_pt500[ietabin]); } //break loop to go faster else break; } } //filling mutliplicity dependent histogramms for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){ _h_dNch_dn[ietabin]->fill(_nch_in_Evt[ietabin] , weight); } //Do only if eta bins are the needed ones if(_etabins[4] == 2.4 && _etabins[0] == 0.5){ if(_nch_in_Evt[4] != 0) _h_dmpt_dNch_eta24->fill(_nch_in_Evt[4] , _sumpt / _nch_in_Evt[4] , weight); _h_dNch_dn_pt500_eta24->fill(_nch_in_Evt_pt500[4] , weight); /*if( _nch_in_Evt[4] < nch_max ){ (_v_dNch_dn_binning1_eta05[_nch_in_Evt[0]])+=weight; (_v_dNch_dn_binning1_eta24[_nch_in_Evt[0]])+=weight; }*/ } else getLog() << Log::WARNING << "You changed the number of eta bins, but forgot to propagate it everywhere !! " << endl; }
void Rivet::CMS_2011_S8884919::finalize | ( | void | ) | [private] |
Definition at line 267 of file CMS_2011_S8884919.cc.
References dtDQMClient_cfg::INFO.
{ getLog() << Log::INFO << "Number of events after event selection: " << _Nevt_after_cuts << endl; /*makeMoments(_v_dNch_dn_binning1_eta05 , _moments_eta05); makeMoments(_v_dNch_dn_binning1_eta24 , _moments_eta24); for(int m = 0 ; m < _moments_eta05.n ; ++m){ _h_cq_eta05->fill(m , _moments_eta05.cq[m]); _h_cq_eta24->fill(m , _moments_eta24.cq[m]); }*/ for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){ normalize(_h_dNch_dn[ietabin]); } normalize(_h_dNch_dn_pt500_eta24); }
void Rivet::CMS_2011_S8884919::init | ( | void | ) | [private] |
Definition at line 130 of file CMS_2011_S8884919.cc.
References lumiQTWidget::t.
{ ChargedFinalState cfs(-2.4, 2.4, 0.0*GeV); addProjection(cfs, "CFS"); addProjection(Beam(), "Beam"); nch_max = 400; _Nevt_after_cuts = 0; //eta bins _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) ; _h_dNch_dn.reserve( _etabins.size() ); ostringstream t(""); if(fuzzyEquals(sqrtS(), 900, 1E-3)){ for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){ t.str("") ; t << "$|\\eta| <$ " << _etabins[ietabin] << " , $\\sqrt(s)$ = 0.9 TeV" ; _h_dNch_dn.push_back( bookHistogram1D( 2 + ietabin, 1, 1 , t.str() , "n" , "$P_{n}$") ); } _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}$"); //_h_cq_eta05 = bookHistogram1D( 17 , 1, 1 , "$|\\eta| <$ 0.5 , $\\sqrt(s)$ = 0.9 TeV" , "q" , "$C_{q}$"); //_h_cq_eta24 = bookHistogram1D( 17 , 1, 2 , "$|\\eta| <$ 2.4 , $\\sqrt(s)$ = 0.9 TeV" , "q" , "$C_{q}$"); _h_dmpt_dNch_eta24 = bookProfile1D( 23 , 1, 1 , "$|\\eta| <$ 2.4 , $\\sqrt(s)$ = 0.9 TeV" , "n" , "$< p_{T}> [ GeV/c ]$"); } if(fuzzyEquals(sqrtS(), 2360, 1E-3)){ for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){ t.str("") ; t << "$|\\eta| <$ " << _etabins[ietabin] << " , $\\sqrt(s)$ = 2.36 TeV" ; _h_dNch_dn.push_back( bookHistogram1D( 7 + ietabin, 1, 1 , t.str() , "n" , "$P_{n}$") ); } _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}$"); //_h_cq_eta05 = bookHistogram1D( 18 , 1, 1 , "$|\\eta| <$ 0.5 , $\\sqrt(s)$ = 2.36 TeV" , "q" , "$C_{q}$"); //_h_cq_eta24 = bookHistogram1D( 18 , 1, 2 , "$|\\eta| <$ 2.4 , $\\sqrt(s)$ = 2.36 TeV" , "q" , "$C_{q}$"); _h_dmpt_dNch_eta24 = bookProfile1D( 24 , 1, 1 , "$|\\eta| <$ 2.4 , $\\sqrt(s)$ = 2.36 TeV" , "n" , "$< p_{T}> [ GeV/c ]$"); } if(fuzzyEquals(sqrtS(), 7000, 1E-3)){ for (unsigned ietabin=0; ietabin < _etabins.size(); ietabin++){ t.str("") ; t << "$|\\eta| <$ " << _etabins[ietabin] << " , $\\sqrt(s)$ = 7 TeV" ; _h_dNch_dn.push_back( bookHistogram1D( 12 + ietabin, 1, 1 , t.str() , "n" , "$P_{n}$") ); } _h_dNch_dn_pt500_eta24 = bookHistogram1D( 22 , 1, 1 , "$p_{T} >$ 500 GeV/c , $|\\eta| <$ 2.4 , $\\sqrt(s)$ = 7 TeV" , "n" , "$P_{n}$"); //_h_cq_eta05 = bookHistogram1D( 19 , 1, 1 , "$|\\eta| <$ 0.5 , $\\sqrt(s)$ = 7 TeV" , "q" , "$C_{q}$"); //_h_cq_eta24 = bookHistogram1D( 19 , 1, 2 , "$|\\eta| <$ 2.4 , $\\sqrt(s)$ = 7 TeV" , "q" , "$C_{q}$"); _h_dmpt_dNch_eta24 = bookProfile1D( 25 , 1, 1 , "$|\\eta| <$2.4 , $\\sqrt(s)$ = 7 TeV" , "n" , "$< p_{T}> [ GeV/c ]$"); } //_v_dNch_dn_binning1_eta05.assign(nch_max+1,0); //_v_dNch_dn_binning1_eta24.assign(nch_max+1,0); }
vector<double> Rivet::CMS_2011_S8884919::_etabins [private] |
Definition at line 55 of file CMS_2011_S8884919.cc.
AIDA::IProfile1D* Rivet::CMS_2011_S8884919::_h_dmpt_dNch_eta24 [private] |
Definition at line 49 of file CMS_2011_S8884919.cc.
vector<AIDA::IHistogram1D*> Rivet::CMS_2011_S8884919::_h_dNch_dn [private] |
Definition at line 39 of file CMS_2011_S8884919.cc.
AIDA::IHistogram1D* Rivet::CMS_2011_S8884919::_h_dNch_dn_pt500_eta24 [private] |
Definition at line 44 of file CMS_2011_S8884919.cc.
AIDA::IHistogram1D* Rivet::CMS_2011_S8884919::_h_KNO_eta05 [private] |
Definition at line 51 of file CMS_2011_S8884919.cc.
AIDA::IHistogram1D* Rivet::CMS_2011_S8884919::_h_KNO_eta24 [private] |
Definition at line 52 of file CMS_2011_S8884919.cc.
vector<int> Rivet::CMS_2011_S8884919::_nch_in_Evt [private] |
Definition at line 56 of file CMS_2011_S8884919.cc.
vector<int> Rivet::CMS_2011_S8884919::_nch_in_Evt_pt500 [private] |
Definition at line 57 of file CMS_2011_S8884919.cc.
double Rivet::CMS_2011_S8884919::_Nevt_after_cuts [private] |
Definition at line 54 of file CMS_2011_S8884919.cc.
double Rivet::CMS_2011_S8884919::_sumpt [private] |
Definition at line 58 of file CMS_2011_S8884919.cc.
int Rivet::CMS_2011_S8884919::nch_max [private] |
Definition at line 59 of file CMS_2011_S8884919.cc.