CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
HiBasicGenTest Class Reference

#include <HiBasicGenTest.h>

Inheritance diagram for HiBasicGenTest:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &r, const edm::EventSetup &c) override
 
 HiBasicGenTest (const edm::ParameterSet &)
 
 ~HiBasicGenTest () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Private Attributes

MonitorElementb [3]
 
MonitorElementdnchdeta [3]
 
MonitorElementdnchdphi [3]
 
MonitorElementdnchdpt [3]
 
edm::EDGetTokenT< edm::HepMCProductgeneratorToken_
 
edm::ESHandle< ParticleDataTablepdt
 
MonitorElementrp
 

Detailed Description

Definition at line 21 of file HiBasicGenTest.h.

Constructor & Destructor Documentation

HiBasicGenTest::HiBasicGenTest ( const edm::ParameterSet iPSet)
explicit

Definition at line 17 of file HiBasicGenTest.cc.

References edm::ParameterSet::getParameter().

18 {
19  generatorToken_ = consumes<edm::HepMCProduct> (
20  iPSet.getParameter<edm::InputTag>("generatorLabel"));
21 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::HepMCProduct > generatorToken_
HiBasicGenTest::~HiBasicGenTest ( )
override

Definition at line 23 of file HiBasicGenTest.cc.

23 {}

Member Function Documentation

void HiBasicGenTest::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 59 of file HiBasicGenTest.cc.

References b, begin, ALCARECOTkAlJpsiMuMu_cff::charge, DEFINE_FWK_MODULE, end, PVValHelper::eta, edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), CaloTowersParam_cfi::mc, AlCaHLTBitMon_ParallelJobs::p, HiggsValidation_cfi::pdg_id, Pi, pi, and EnergyCorrector::pt.

60 {
61 
63  iEvent.getByToken(generatorToken_, mc);
64  const HepMC::GenEvent *evt = mc->GetEvent();
65  const HepMC::HeavyIon *hi = evt->heavy_ion();
66 
67  int cbin = 0;
68  double phi0 =0.;
69 
70  if(hi){
71 
72  double ip = hi->impact_parameter();
73  phi0 = hi->event_plane_angle();
74 
75  // fill reaction plane distribution
76  rp->Fill(phi0);
77 
78  // if the event is in one of the centrality bins of interest fill hists
79  int cbin=-1;
80  if(ip < 5.045) cbin=0;
81  else if (ip < 7.145 && ip > 5.045) cbin=1;
82  else if (ip < 15.202 && ip > 14.283) cbin=2;
83  if(cbin<0) return;
84 
85  // fill impact parameter distributions
86  b[cbin]->Fill(ip);
87  }
88 
89  // loop over particles
90  HepMC::GenEvent::particle_const_iterator begin = evt->particles_begin();
91  HepMC::GenEvent::particle_const_iterator end = evt->particles_end();
92  for(HepMC::GenEvent::particle_const_iterator it = begin; it != end; ++it){
93 
94  // only fill hists for status=1 particles
95  if((*it)->status() != 1) continue;
96 
97  // only fill hists for charged particles
98  int pdg_id = (*it)->pdg_id();
99  const ParticleData * part = pdt->particle(pdg_id);
100  int charge = static_cast<int>(part->charge());
101  if(charge==0) continue;
102 
103  float eta = (*it)->momentum().eta();
104  float phi = (*it)->momentum().phi();
105  float pt = (*it)->momentum().perp();
106 
107  dnchdeta[cbin]->Fill(eta);
108  dnchdpt[cbin]->Fill(pt);
109 
110  double pi = TMath::Pi();
111  double p = phi-phi0;
112  if(p > pi) p = p - 2*pi;
113  if(p < -1*pi) p = p + 2*pi;
114  dnchdphi[cbin]->Fill(p);
115 
116  }
117 
118  return;
119 
120 }
const double Pi
MonitorElement * rp
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
edm::EDGetTokenT< edm::HepMCProduct > generatorToken_
MonitorElement * dnchdpt[3]
edm::ESHandle< ParticleDataTable > pdt
void Fill(long long x)
const Double_t pi
MonitorElement * dnchdeta[3]
MonitorElement * dnchdphi[3]
MonitorElement * b[3]
#define end
Definition: vmac.h:39
HepPDT::ParticleData ParticleData
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
part
Definition: HCALResponse.h:20
#define begin
Definition: vmac.h:32
void HiBasicGenTest::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
override

Setting the DQM top directories

Booking the ME's

Definition at line 25 of file HiBasicGenTest.cc.

References b, DQMStore::IBooker::book1D(), DQMStore::IBooker::setCurrentFolder(), and DQMStore::IBooker::tag().

26  {
27 
29  ibooker.setCurrentFolder("Generator/Particles");
30 
32  for (int ibin = 0; ibin < 3; ibin++) {
33  dnchdeta[ibin] = ibooker.book1D(Form("dnchdeta%d", ibin), ";#eta;dN^{ch}/d#eta",
34  100, -6.0, 6.0);
35 
36  dnchdpt[ibin] = ibooker.book1D(Form("dnchdpt%d", ibin), ";p_{T};dN^{ch}/dp_{T}",
37  200, 0.0, 100.0);
38 
39  b[ibin] = ibooker.book1D(Form("b%d",ibin),";b[fm];events", 100, 0.0, 20.0);
40  dnchdphi[ibin] = ibooker.book1D(Form("dnchdphi%d", ibin), ";#phi;dN^{ch}/d#phi",
41  100, -3.2, 3.2);
42 
43  ibooker.tag(dnchdeta[ibin], 1+ibin*4);
44  ibooker.tag(dnchdpt[ibin], 2+ibin*4);
45  ibooker.tag(b[ibin], 3+ibin*4);
46  ibooker.tag(dnchdphi[ibin], 4+ibin*4);
47  }
48 
49  rp = ibooker.book1D("phi0", ";#phi_{RP};events", 100, -3.2, 3.2);
50  ibooker.tag(rp, 13);
51 }
MonitorElement * rp
MonitorElement * dnchdpt[3]
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * dnchdeta[3]
MonitorElement * dnchdphi[3]
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
void tag(MonitorElement *, unsigned int)
Definition: DQMStore.cc:283
MonitorElement * b[3]
void HiBasicGenTest::dqmBeginRun ( const edm::Run r,
const edm::EventSetup c 
)
override

Definition at line 53 of file HiBasicGenTest.cc.

References edm::EventSetup::getData().

54 {
55  iSetup.getData(pdt);
56 }
edm::ESHandle< ParticleDataTable > pdt

Member Data Documentation

MonitorElement* HiBasicGenTest::b[3]
private

Definition at line 37 of file HiBasicGenTest.h.

MonitorElement* HiBasicGenTest::dnchdeta[3]
private

Definition at line 35 of file HiBasicGenTest.h.

MonitorElement* HiBasicGenTest::dnchdphi[3]
private

Definition at line 38 of file HiBasicGenTest.h.

MonitorElement* HiBasicGenTest::dnchdpt[3]
private

Definition at line 36 of file HiBasicGenTest.h.

edm::EDGetTokenT<edm::HepMCProduct> HiBasicGenTest::generatorToken_
private

Definition at line 34 of file HiBasicGenTest.h.

edm::ESHandle< ParticleDataTable > HiBasicGenTest::pdt
private

Definition at line 41 of file HiBasicGenTest.h.

MonitorElement* HiBasicGenTest::rp
private

Definition at line 39 of file HiBasicGenTest.h.