CMS 3D CMS Logo

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

#include <BasicHepMCValidation.h>

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

Classes

class  ParticleMonitor
 

Public Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
 BasicHepMCValidation (const edm::ParameterSet &)
 
void bookHistograms (DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &r, const edm::EventSetup &c) override
 
 ~BasicHepMCValidation () 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

MonitorElementBjorken_x
 
MonitorElementDeltaEcms
 
MonitorElementDeltaPx
 
MonitorElementDeltaPy
 
MonitorElementDeltaPz
 
edm::ESHandle< HepPDT::ParticleDataTablefPDGTable
 PDT table. More...
 
MonitorElementgenPtclNumber
 
MonitorElementgenPtclStatus
 
MonitorElementgenVrtxNumber
 
edm::InputTag hepmcCollection_
 
edm::EDGetTokenT< edm::HepMCProducthepmcCollectionToken_
 
MonitorElementlog10DeltaEcms
 
MonitorElementnEvt
 
MonitorElementotherPtclMomentum
 
MonitorElementotherPtclNumber
 other ME's More...
 
MonitorElementoutVrtxPtclNumber
 
MonitorElementoutVrtxStablePtclNumber
 
std::vector< ParticleMonitorparticles
 
MonitorElementparton1Id
 
MonitorElementparton2Id
 
MonitorElementpartonNumber
 
MonitorElementpartonpT
 
MonitorElementpdf_bbbar
 
MonitorElementpdf_ccbar
 
MonitorElementpdf_d
 
MonitorElementpdf_dbar
 
MonitorElementpdf_g
 
MonitorElementpdf_ssbar
 
MonitorElementpdf_u
 
MonitorElementpdf_ubar
 
MonitorElementscalePDF
 
MonitorElementstableChaNumber
 
MonitorElementstablePtclCharge
 
MonitorElementstablePtclEta
 
MonitorElementstablePtclNumber
 
MonitorElementstablePtclp
 
MonitorElementstablePtclPhi
 
MonitorElementstablePtclpT
 
MonitorElementstatus1ShortLived
 
MonitorElementunknownPDTNumber
 
MonitorElementvrtxRadius
 
MonitorElementvrtxZ
 
WeightManager wmanager_
 

Detailed Description

Definition at line 35 of file BasicHepMCValidation.h.

Constructor & Destructor Documentation

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

Definition at line 14 of file BasicHepMCValidation.cc.

References hepmcCollection_, and hepmcCollectionToken_.

14  :
15  wmanager_(iPSet,consumesCollector()),
16  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection"))
17 {
18  hepmcCollectionToken_=consumes<HepMCProduct>(hepmcCollection_);
19 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
edm::InputTag hepmcCollection_
BasicHepMCValidation::~BasicHepMCValidation ( )
override

Definition at line 21 of file BasicHepMCValidation.cc.

21 {}

Member Function Documentation

void BasicHepMCValidation::analyze ( edm::Event const &  ,
edm::EventSetup const &   
)
override

counters to zero for every event

Gathering the HepMCProduct information

PDF informations

Vertices

loop on vertex particles

Looping through the PARTICLES in the event

Particles

Status statistics

Stable particles

filling multiplicity ME's

Definition at line 166 of file BasicHepMCValidation.cc.

References funct::abs(), Bjorken_x, ALCARECOTkAlJpsiMuMu_cff::charge, DeltaEcms, DeltaPx, DeltaPy, DeltaPz, MonitorElement::Fill(), fPDGTable, GenParticle::GenParticle, genPtclNumber, genPtclStatus, genVrtxNumber, edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), hepmcCollectionToken_, mps_fire::i, log10DeltaEcms, logPdfBinsize, logPdfMax, logPdfMin, logQScaleBinsize, nEvt, otherPtclMomentum, otherPtclNumber, outVrtxPtclNumber, outVrtxStablePtclNumber, source_particleGun_cfi::ParticleID, particles, parton1Id, parton2Id, partonNumber, partonpT, pdf_bbbar, pdf_ccbar, pdf_d, pdf_dbar, pdf_g, pdf_ssbar, pdf_u, pdf_ubar, scalePDF, stableChaNumber, stablePtclCharge, stablePtclEta, stablePtclNumber, stablePtclp, stablePtclPhi, stablePtclpT, mps_update::status, status1ShortLived, unknownPDTNumber, vrtxRadius, vrtxZ, WeightManager::weight(), mps_merge::weight, and wmanager_.

167 {
169  int partonNum = 0;
170  //
171  int outVrtxStablePtclNum = 0; int stablePtclNum = 0; int otherPtclNum = 0; int unknownPDTNum = 0; int stableChaNum = 0;
172  //
173  double bjorken = 0.; double logPdf1 = 0.; double logPdf2 = 0.; double logQScale = 0.;
174  //
175  double etotal = 0. ; double pxtotal = 0.; double pytotal = 0.; double pztotal = 0.;
176 
179  iEvent.getByToken(hepmcCollectionToken_, evt);
180 
181  //Get EVENT
182  HepMC::GenEvent const *myGenEvent = evt->GetEvent();
183 
184  double weight = wmanager_.weight(iEvent);
185 
186  nEvt->Fill(0.5,weight);
187 
188  genPtclNumber->Fill(log10(myGenEvent->particles_size()),weight);
189  genVrtxNumber->Fill(log10(myGenEvent->vertices_size()),weight);
190 
192  HepMC::PdfInfo const *pdf = myGenEvent->pdf_info();
193  if(pdf){
194  bjorken = ((pdf->x1())/((pdf->x1())+(pdf->x2())));
195  logQScale = log10(pdf->scalePDF());
196  if(logQScale > logQScaleMax) logQScale = logQScaleMax - 0.5*logQScaleBinsize; // visualize overflow & underflow in the histograms
197  if(logQScale < logQScaleMin) logQScale = logQScaleMin + 0.5*logQScaleBinsize;
198  logPdf1 = log10(pdf->pdf1());
199  if(logPdf1 > logPdfMax) logPdf1 = logPdfMax - 0.5*logPdfBinsize;
200  if(logPdf1 < logPdfMin) logPdf1 = logPdfMin + 0.5*logPdfBinsize;
201  logPdf2 = log10(pdf->pdf2());
202  if(logPdf2 > logPdfMax) logPdf2 = logPdfMax - 0.5*logPdfBinsize;
203  if(logPdf2 < logPdfMin) logPdf2 = logPdfMin + 0.5*logPdfBinsize;
204  Bjorken_x->Fill(bjorken,weight);
205  scalePDF->Fill(logQScale,weight);
206  parton1Id->Fill((double)pdf->id1(),weight);
207  parton2Id->Fill((double)pdf->id2(),weight);
208  if(pdf->id1() == 2) pdf_u->Fill(logPdf1,weight);
209  if(pdf->id2() == 2) pdf_u->Fill(logPdf2,weight);
210  if(pdf->id1() == -2) pdf_ubar->Fill(logPdf1,weight);
211  if(pdf->id2() == -2) pdf_ubar->Fill(logPdf2,weight);
212  if(pdf->id1() == 1) pdf_d->Fill(logPdf1,weight);
213  if(pdf->id2() == 1) pdf_d->Fill(logPdf2,weight);
214  if(pdf->id1() == -1) pdf_dbar->Fill(logPdf1,weight);
215  if(pdf->id2() == -1) pdf_dbar->Fill(logPdf2,weight);
216  if(std::abs(pdf->id1()) == 3) pdf_ssbar->Fill(logPdf1,weight);
217  if(std::abs(pdf->id2()) == 3) pdf_ssbar->Fill(logPdf2,weight);
218  if(std::abs(pdf->id1()) == 4) pdf_ccbar->Fill(logPdf1,weight);
219  if(std::abs(pdf->id2()) == 4) pdf_ccbar->Fill(logPdf2,weight);
220  if(std::abs(pdf->id1()) == 5) pdf_bbbar->Fill(logPdf1,weight);
221  if(std::abs(pdf->id2()) == 5) pdf_bbbar->Fill(logPdf2,weight);
222  if(std::abs(pdf->id1()) == 21) pdf_g->Fill(logPdf1,weight);
223  if(std::abs(pdf->id2()) == 21) pdf_g->Fill(logPdf2,weight);
224  }
225 
226  //Looping through the VERTICES in the event
227  HepMC::GenEvent::vertex_const_iterator vrtxBegin = myGenEvent->vertices_begin();
228  HepMC::GenEvent::vertex_const_iterator vrtxEnd = myGenEvent->vertices_end();
229  unsigned int nvtx(0);
230  for(HepMC::GenEvent::vertex_const_iterator vrtxIt = vrtxBegin; vrtxIt!=vrtxEnd; ++vrtxIt)
231  {
233  HepMC::GenVertex const *vrtx = *vrtxIt;
234  outVrtxPtclNumber->Fill(vrtx->particles_out_size(),weight); //std::cout << "all " << vrtx->particles_out_size() << '\n';
235 
236  if(nvtx==0){
237  vrtxZ->Fill(vrtx->point3d().z(),weight);
238  vrtxRadius->Fill(vrtx->point3d().perp(),weight);
239  }
241  HepMC::GenVertex::particles_out_const_iterator vrtxPtclBegin = vrtx->particles_out_const_begin();
242  HepMC::GenVertex::particles_out_const_iterator vrtxPtclEnd = vrtx->particles_out_const_end();
243  outVrtxStablePtclNum = 0;
244  for(HepMC::GenVertex::particles_out_const_iterator vrtxPtclIt = vrtxPtclBegin; vrtxPtclIt != vrtxPtclEnd; ++vrtxPtclIt)
245  {
246  HepMC::GenParticle const *vrtxPtcl = *vrtxPtclIt;
247  if (vrtxPtcl->status() == 1){
248  ++outVrtxStablePtclNum; //std::cout << "stable " << outVrtxStablePtclNum << '\n';
249  }
250  }
251  outVrtxStablePtclNumber->Fill(outVrtxStablePtclNum,weight);
252  nvtx++;
253  }//vertices
254 
256  HepMC::GenEvent::particle_const_iterator ptclBegin = myGenEvent->particles_begin();
257  HepMC::GenEvent::particle_const_iterator ptclEnd = myGenEvent->particles_end();
258  for(HepMC::GenEvent::particle_const_iterator ptclIt = ptclBegin; ptclIt!=ptclEnd; ++ptclIt)
259  {
260 
262  HepMC::GenParticle const *ptcl = *ptclIt;
263  int Id = ptcl->pdg_id(); // std::cout << Id << '\n';
264  float Log_p = log10( ptcl->momentum().rho() );
265  double charge = 999.; // for the charge it's needed a HepPDT method
266  int status = ptcl->status();
267  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID(Id));
268  if(PData==nullptr) {
269  // std::cout << "Unknown id = " << Id << '\n';
270  ++unknownPDTNum;
271  }
272  else
273  charge = PData->charge();
274 
276  genPtclStatus->Fill((float)status,weight);
277 
279  if(ptcl->status() == 1){
280  ++stablePtclNum;
281  stablePtclPhi->Fill(ptcl->momentum().phi()/CLHEP::degree,weight); //std::cout << ptcl->polarization().phi() << '\n';
282  stablePtclEta->Fill(ptcl->momentum().pseudoRapidity(),weight);
283  stablePtclCharge->Fill(charge,weight); // std::cout << ptclData.charge() << '\n';
284  stablePtclp->Fill(Log_p,weight);
285  stablePtclpT->Fill(log10(ptcl->momentum().perp()),weight);
286  if (charge != 0. && charge != 999.) ++stableChaNum;
287  if ( std::abs(Id) == 1 ) status1ShortLived->Fill(1,weight);
288  if ( std::abs(Id) == 2 ) status1ShortLived->Fill(2,weight);
289  if ( std::abs(Id) == 3 ) status1ShortLived->Fill(3,weight);
290  if ( std::abs(Id) == 4 ) status1ShortLived->Fill(4,weight);
291  if ( std::abs(Id) == 5 ) status1ShortLived->Fill(5,weight);
292  if ( std::abs(Id) == 6 ) status1ShortLived->Fill(6,weight);
293  if ( Id == 21 ) status1ShortLived->Fill(7,weight);
294  if ( std::abs(Id) == 15 ) status1ShortLived->Fill(8,weight);
295  if ( Id == 23 ) status1ShortLived->Fill(9,weight);
296  if ( std::abs(Id) == 24 ) status1ShortLived->Fill(10,weight);
297  if ( std::abs(Id) == 7 || std::abs(Id) == 8 || std::abs(Id) == 17 || (std::abs(Id) >= 25 && std::abs(Id) <= 99) ) status1ShortLived->Fill(11,weight);
298  etotal += ptcl->momentum().e();
299  pxtotal += ptcl->momentum().px();
300  pytotal += ptcl->momentum().py();
301  pztotal += ptcl->momentum().pz();
302  }
303 
304  if (abs(Id) < 6 || abs(Id) == 22){
305  ++partonNum; partonpT->Fill(Log_p,weight);
306  }
307 
308  bool indentified=false;
309  for(unsigned int i=0;i<particles.size();i++){if(particles.at(i).Fill(ptcl,weight)){indentified=true; break;}}
310  if(!indentified){ ++otherPtclNum; otherPtclMomentum->Fill(Log_p,weight);}
311  }//event particles
312 
313  // set a default sqrt(s) and then check in the event
314  double ecms = 13000.;
315  if ( myGenEvent->valid_beam_particles() ) {
316  ecms = myGenEvent->beam_particles().first->momentum().e()+myGenEvent->beam_particles().second->momentum().e();
317  }
318  log10DeltaEcms->Fill(log10(fabs(etotal-ecms)),weight);
319  DeltaEcms->Fill(etotal-ecms,weight);
320  DeltaPx->Fill(pxtotal,weight);
321  DeltaPy->Fill(pytotal,weight);
322  DeltaPz->Fill(pztotal,weight);
323 
324 
326  stablePtclNumber->Fill(log10(stablePtclNum+0.1),weight);
327  stableChaNumber->Fill(log10(stableChaNum+0.1),weight);
328  otherPtclNumber->Fill(log10(otherPtclNum+0.1),weight);
329  unknownPDTNumber->Fill(log10(unknownPDTNum+0.1),weight);
330  //
331  partonNumber->Fill(partonNum,weight);
332  for(unsigned int i=0;i<particles.size();i++){particles.at(i).FillCount(weight);};
333 
334 }//analyze
MonitorElement * genPtclStatus
MonitorElement * DeltaEcms
MonitorElement * genVrtxNumber
MonitorElement * unknownPDTNumber
MonitorElement * pdf_bbbar
MonitorElement * outVrtxPtclNumber
MonitorElement * pdf_ccbar
MonitorElement * otherPtclNumber
other ME&#39;s
MonitorElement * Bjorken_x
MonitorElement * parton2Id
MonitorElement * otherPtclMomentum
MonitorElement * stablePtclCharge
Definition: weight.py:1
MonitorElement * partonNumber
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
void Fill(long long x)
double logPdfMin
int iEvent
Definition: GenABIO.cc:224
MonitorElement * vrtxRadius
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * log10DeltaEcms
MonitorElement * stablePtclPhi
MonitorElement * stablePtclEta
double logPdfBinsize
HepPDT::ParticleData ParticleData
MonitorElement * stablePtclp
MonitorElement * status1ShortLived
MonitorElement * stablePtclpT
MonitorElement * stableChaNumber
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
double logPdfMax
MonitorElement * stablePtclNumber
double logQScaleBinsize
double logQScaleMin
double logQScaleMax
MonitorElement * genPtclNumber
MonitorElement * outVrtxStablePtclNumber
double weight(const edm::Event &)
std::vector< ParticleMonitor > particles
MonitorElement * pdf_ssbar
MonitorElement * parton1Id
void BasicHepMCValidation::bookHistograms ( DQMStore::IBooker i,
edm::Run const &  ,
edm::EventSetup const &   
)
override

Setting the DQM top directories

Booking the ME's multiplicity

other

Definition at line 31 of file BasicHepMCValidation.cc.

References Bjorken_x, DQMHelper::book1dHisto(), DeltaEcms, DeltaPx, DeltaPy, DeltaPz, genPtclNumber, genPtclStatus, genVrtxNumber, log10DeltaEcms, logPdfMax, logPdfMin, logPdfNbin, nEvt, otherPtclMomentum, otherPtclNumber, outVrtxPtclNumber, outVrtxStablePtclNumber, particles, parton1Id, parton2Id, partonNumber, partonpT, pdf_bbbar, pdf_ccbar, pdf_d, pdf_dbar, pdf_g, pdf_ssbar, pdf_u, pdf_ubar, scalePDF, MonitorElement::setBinLabel(), DQMStore::IBooker::setCurrentFolder(), stableChaNumber, stablePtclCharge, stablePtclEta, stablePtclNumber, stablePtclp, stablePtclPhi, stablePtclpT, status1ShortLived, unknownPDTNumber, vrtxRadius, and vrtxZ.

31  {
33  DQMHelper dqm(&i); i.setCurrentFolder("Generator/Particles");
34 
35  // Number of analyzed events
36  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1.);
37 
40  // quarks
41  particles.push_back(ParticleMonitor("u",1,i));
42  particles.push_back(ParticleMonitor("ubar",-1,i));
43  particles.push_back(ParticleMonitor("d",2,i));
44  particles.push_back(ParticleMonitor("dbar",-2,i));
45  particles.push_back(ParticleMonitor("s",3,i));
46  particles.push_back(ParticleMonitor("sbar",-3,i));
47  particles.push_back(ParticleMonitor("c",4,i));
48  particles.push_back(ParticleMonitor("cbar",-4,i));
49  particles.push_back(ParticleMonitor("b",5,i));
50  particles.push_back(ParticleMonitor("bbar",-5,i));
51  particles.push_back(ParticleMonitor("t",6,i));
52  particles.push_back(ParticleMonitor("tbar",-6,i));
53 
54  //leptons
55  particles.push_back(ParticleMonitor("eminus",11,i));
56  particles.push_back(ParticleMonitor("eplus",-11,i));
57  particles.push_back(ParticleMonitor("nue",12,i));
58  particles.push_back(ParticleMonitor("nuebar",-12,i));
59  particles.push_back(ParticleMonitor("muminus",13,i));
60  particles.push_back(ParticleMonitor("muplus",-13,i));
61  particles.push_back(ParticleMonitor("numu",14,i));
62  particles.push_back(ParticleMonitor("numubar",-14,i));
63  particles.push_back(ParticleMonitor("tauminus",15,i));
64  particles.push_back(ParticleMonitor("tauplus",-15,i));
65  particles.push_back(ParticleMonitor("nutau",16,i));
66  particles.push_back(ParticleMonitor("nutaubar",-16,i));
67 
68  //bosons
69  particles.push_back(ParticleMonitor("Wplus",24,i));
70  particles.push_back(ParticleMonitor("Wminus",-24,i));
71  particles.push_back(ParticleMonitor("Z",23,i));
72  particles.push_back(ParticleMonitor("gamma",22,i));
73  particles.push_back(ParticleMonitor("gluon",21,i));
74 
75  //mesons
76  particles.push_back(ParticleMonitor("piplus",211,i,true)); //log
77  particles.push_back(ParticleMonitor("piminus",-211,i,true)); //log
78  particles.push_back(ParticleMonitor("pizero",111,i,true)); //log
79  particles.push_back(ParticleMonitor("Kplus",321,i));
80  particles.push_back(ParticleMonitor("Kminus",-321,i));
81  particles.push_back(ParticleMonitor("Klzero",130,i));
82  particles.push_back(ParticleMonitor("Kszero",310,i));
83 
84  //baryons
85  particles.push_back(ParticleMonitor("p",2212,i,true)); //log
86  particles.push_back(ParticleMonitor("pbar",-2212,i,true)); //log
87  particles.push_back(ParticleMonitor("n",2112,i,true)); //log
88  particles.push_back(ParticleMonitor("nbar",-2112,i,true)); //log
89  particles.push_back(ParticleMonitor("lambda0",3122,i));
90  particles.push_back(ParticleMonitor("lambda0bar",-3122,i));
91 
92  //D mesons
93  particles.push_back(ParticleMonitor("Dplus",411,i));
94  particles.push_back(ParticleMonitor("Dminus",-411,i));
95  particles.push_back(ParticleMonitor("Dzero",421,i));
96  particles.push_back(ParticleMonitor("Dzerobar",-421,i));
97 
98  //B mesons
99  particles.push_back(ParticleMonitor("Bplus",521,i));
100  particles.push_back(ParticleMonitor("Bminus",-521,i));
101  particles.push_back(ParticleMonitor("Bzero",511,i));
102  particles.push_back(ParticleMonitor("Bzerobar",-511,i));
103  particles.push_back(ParticleMonitor("Bszero",531,i));
104  particles.push_back(ParticleMonitor("Bszerobar",-531,i));
105 
106  //
107  otherPtclNumber = dqm.book1dHisto("otherPtclNumber", "Log10(No. other ptcls)", 60, -1, 5,"log_{10}(No. other ptcls)","Number of Events"); //Log
108  otherPtclMomentum = dqm.book1dHisto("otherPtclMomentum", "Log10(p) other ptcls", 60, -2, 4,"log10(P^{other ptcls}) (log_{10}(GeV))","Number of Events");
109 
111  genPtclNumber = dqm.book1dHisto("genPtclNumber", "Log10(No. all particles)", 60, -1, 5,"log10(No. all particles)","Number of Events"); //Log
112  genVrtxNumber = dqm.book1dHisto("genVrtxNumber", "Log10(No. all vertexs)", 60, -1, 5,"log10(No. all vertexs)","Number of Events"); //Log
113  //
114  stablePtclNumber= dqm.book1dHisto("stablePtclNumber", "Log10(No. stable particles)", 50, 0, 5,"log10(No. stable particles)","Number of Events"); //Log
115  stablePtclPhi = dqm.book1dHisto("stablePtclPhi", "stable Ptcl Phi", 360, -180, 180,"#phi^{stable Ptcl} (rad)","Number of Events");
116  stablePtclEta = dqm.book1dHisto("stablePtclEta", "stable Ptcl Eta (pseudo rapidity)", 220, -11, 11,"#eta^{stable Ptcl}","Number of Events");
117  stablePtclCharge = dqm.book1dHisto("stablePtclCharge", "stablePtclCharge", 5, -2, 2,"Charge^{stable ptcls}","Number of Events");
118  stableChaNumber= dqm.book1dHisto("stableChaNumber", "Log10(No. stable charged particles)", 50, 0, 5,"log_{10}(No. stable charged particles)","Number of Events"); //Log
119  stablePtclp = dqm.book1dHisto("stablePtclp", "Log10(p) stable ptcl p", 80, -4, 4,"log_{10}(P^{stable ptcl}) (log_{10}(GeV))","Number of Events"); //Log
120  stablePtclpT = dqm.book1dHisto("stablePtclpT", "Log10(pT) stable ptcl pT", 80, -4, 4,"log_{10}(P_{t}^{stable ptcl}) (log_{10}(GeV))","Number of Events"); //Log
121  partonNumber = dqm.book1dHisto("partonNumber", "number of partons", 100, 0, 100,"number of partons","Number of Events");
122  partonpT = dqm.book1dHisto("partonpT", "Log10(pT) parton pT", 80, -4, 4,"Log10(P_{t}^{parton})","Number of Events"); //Log
123  outVrtxStablePtclNumber = dqm.book1dHisto("outVrtxStablePtclNumber", "No. outgoing stable ptcls from vrtx", 10, 0, 10,"No. outgoing stable ptcls from vrtx","Number of Events");
124  //
125  outVrtxPtclNumber = dqm.book1dHisto("outVrtxPtclNumber", "No. outgoing ptcls from vrtx", 30, 0, 30,"No. outgoing ptcls from vrtx","Number of Events");
126  vrtxZ = dqm.book1dHisto("VrtxZ", "VrtxZ", 50 , -250, 250,"Z_{Vtx}","Number of Events");
127  vrtxRadius = dqm.book1dHisto("vrtxRadius", "vrtxRadius", 50, 0, 50,"R_{vtx}","Number of Events");
128  //
129  unknownPDTNumber = dqm.book1dHisto("unknownPDTNumber", "Log10(No. unknown ptcls PDT)", 60, -1, 5,"log_{10}(No. unknown ptcls PDT)","Number of Events"); //Log
130  genPtclStatus = dqm.book1dHisto("genPtclStatus", "Status of genParticle", 200,0,200.,"","Number of Events");
131  //
132  Bjorken_x = dqm.book1dHisto("Bjorken_x", "Bjorken_x", 1000, 0.0, 1.0,"Bjorken_{x}","Number of Events");
133  pdf_u = dqm.book1dHisto("pdf_u", "Log10(PDF(u,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
134  pdf_ubar = dqm.book1dHisto("pdf_ubar", "Log10(PDF(ubar,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
135  pdf_d = dqm.book1dHisto("pdf_d", "Log10(PDF(d,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
136  pdf_dbar = dqm.book1dHisto("pdf_dbar", "Log10(PDF(dbar,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
137  pdf_ssbar = dqm.book1dHisto("pdf_ssbar", "Log10(PDF(ssbar,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
138  pdf_ccbar = dqm.book1dHisto("pdf_ccbar", "Log10(PDF(ccbar,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
139  pdf_bbbar = dqm.book1dHisto("pdf_bbbar", "Log10(PDF(bbbar,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
140  pdf_g = dqm.book1dHisto("pdf_g", "Log10(PDF(g,x,Q))", logPdfNbin, logPdfMin, logPdfMax,"log_{10}(x*f(x))","Number of Events"); //Log
141  scalePDF = dqm.book1dHisto("scalePDF", "Log10(Q-scale(GeV))", logQScaleNbin, logQScaleMin, logQScaleMax,"log_{10}(Q-scale(GeV))","Number of Events");
142  parton1Id = dqm.book1dHisto("parton1Id", "ID of parton 1", 45, -14.5, 30.5,"ID","Number of Events");
143  parton2Id = dqm.book1dHisto("parton2Id", "ID of parton 2", 45, -14.5, 30.5,"ID","Number of Events");
144  //
145  status1ShortLived = dqm.book1dHisto("status1ShortLived","Status 1 short lived", 11, 0, 11,"","Number of Events");
146  status1ShortLived->setBinLabel(1,"d/dbar");
147  status1ShortLived->setBinLabel(2,"u/ubar");
148  status1ShortLived->setBinLabel(3,"s/sbar");
149  status1ShortLived->setBinLabel(4,"c/cbar");
150  status1ShortLived->setBinLabel(5,"b/bbar");
151  status1ShortLived->setBinLabel(6,"t/tbar");
153  status1ShortLived->setBinLabel(8,"tau-/tau+");
155  status1ShortLived->setBinLabel(10,"W-/W+");
156  status1ShortLived->setBinLabel(11,"PDG = 7,8,17,25-99");
157 
158  log10DeltaEcms =dqm.book1dHisto("DeltaEcms1log10","log_{10} of deviation from nominal Ecms", 200,-1., 5.,"log_{10}(#DeltaE) (log_{10}(GeV))","Number of Events");
159  DeltaEcms = dqm.book1dHisto("DeltaEcms1","deviation from nominal Ecms", 200,-1., 1.,"#DeltaE (GeV)","Number of Events");
160  DeltaPx = dqm.book1dHisto("DeltaPx1","deviation from nominal Px", 200,-1., 1.,"#DeltaP_{x} (GeV)","Number of Events");
161  DeltaPy = dqm.book1dHisto("DeltaPy1","deviation from nominal Py", 200,-1., 1.,"#DeltaP_{y} (GeV)","Number of Events");
162  DeltaPz = dqm.book1dHisto("DeltaPz1","deviation from nominal Pz", 200,-1., 1.,"#DeltaP_{z} (GeV)","Number of Events");
163  return;
164 }
MonitorElement * genPtclStatus
MonitorElement * DeltaEcms
MonitorElement * genVrtxNumber
MonitorElement * unknownPDTNumber
MonitorElement * pdf_bbbar
MonitorElement * outVrtxPtclNumber
MonitorElement * pdf_ccbar
MonitorElement * Bjorken_x
MonitorElement * otherPtclNumber
other ME&#39;s
MonitorElement * parton2Id
MonitorElement * otherPtclMomentum
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * stablePtclCharge
MonitorElement * partonNumber
double logPdfMin
MonitorElement * vrtxRadius
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
int logQScaleNbin
MonitorElement * log10DeltaEcms
MonitorElement * stablePtclPhi
MonitorElement * stablePtclEta
int logPdfNbin
MonitorElement * stablePtclp
MonitorElement * status1ShortLived
MonitorElement * stablePtclpT
MonitorElement * stableChaNumber
double logPdfMax
MonitorElement * stablePtclNumber
double logQScaleMin
double logQScaleMax
MonitorElement * genPtclNumber
MonitorElement * outVrtxStablePtclNumber
std::vector< ParticleMonitor > particles
MonitorElement * pdf_ssbar
MonitorElement * parton1Id
void BasicHepMCValidation::dqmBeginRun ( const edm::Run r,
const edm::EventSetup c 
)
override

Definition at line 23 of file BasicHepMCValidation.cc.

References fPDGTable, and edm::EventSetup::getData().

23  {
24  c.getData( fPDGTable );
25 }
bool getData(T &iHolder) const
Definition: EventSetup.h:111
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.

Member Data Documentation

MonitorElement* BasicHepMCValidation::Bjorken_x
private

Definition at line 185 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaEcms
private

Definition at line 201 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaPx
private

Definition at line 202 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaPy
private

Definition at line 203 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaPz
private

Definition at line 204 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

edm::ESHandle<HepPDT::ParticleDataTable> BasicHepMCValidation::fPDGTable
private

PDT table.

Definition at line 49 of file BasicHepMCValidation.h.

Referenced by analyze(), and dqmBeginRun().

MonitorElement* BasicHepMCValidation::genPtclNumber
private

Definition at line 165 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::genPtclStatus
private

Definition at line 169 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::genVrtxNumber
private

Definition at line 166 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

edm::InputTag BasicHepMCValidation::hepmcCollection_
private

Definition at line 46 of file BasicHepMCValidation.h.

Referenced by BasicHepMCValidation().

edm::EDGetTokenT<edm::HepMCProduct> BasicHepMCValidation::hepmcCollectionToken_
private

Definition at line 206 of file BasicHepMCValidation.h.

Referenced by analyze(), and BasicHepMCValidation().

MonitorElement* BasicHepMCValidation::log10DeltaEcms
private

Definition at line 200 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::nEvt
private

Definition at line 159 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::otherPtclMomentum
private

Definition at line 164 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::otherPtclNumber
private

other ME's

Definition at line 163 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::outVrtxPtclNumber
private

Definition at line 168 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::outVrtxStablePtclNumber
private

Definition at line 180 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

std::vector<ParticleMonitor> BasicHepMCValidation::particles
private

Definition at line 160 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::parton1Id
private

Definition at line 195 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::parton2Id
private

Definition at line 196 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::partonNumber
private

Definition at line 178 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::partonpT
private

Definition at line 179 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_bbbar
private

Definition at line 192 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_ccbar
private

Definition at line 191 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_d
private

Definition at line 188 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_dbar
private

Definition at line 189 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_g
private

Definition at line 193 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_ssbar
private

Definition at line 190 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_u
private

Definition at line 186 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_ubar
private

Definition at line 187 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::scalePDF
private

Definition at line 194 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stableChaNumber
private

Definition at line 172 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclCharge
private

Definition at line 175 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclEta
private

Definition at line 174 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclNumber
private

Definition at line 171 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclp
private

Definition at line 176 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclPhi
private

Definition at line 173 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclpT
private

Definition at line 177 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::status1ShortLived
private

Definition at line 198 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::unknownPDTNumber
private

Definition at line 167 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::vrtxRadius
private

Definition at line 183 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::vrtxZ
private

Definition at line 182 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

WeightManager BasicHepMCValidation::wmanager_
private

Definition at line 45 of file BasicHepMCValidation.h.

Referenced by analyze().