CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BasicHepMCValidation.cc
Go to the documentation of this file.
1 /*class BasicHepMCValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  *
5  */
6 
8 
9 #include "CLHEP/Units/defs.h"
10 #include "CLHEP/Units/PhysicalConstants.h"
12 using namespace edm;
13 
15  : wmanager_(iPSet, consumesCollector()), hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")) {
16  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
17  fPDGTableToken = esConsumes<edm::Transition::BeginRun>();
18 }
19 
21 
24 }
25 
26 namespace {
27  // Set upper bound & lower bound for PDF & Scale related histograms
28  constexpr double logPdfMax = 3.0;
29  constexpr double logPdfMin = -3.0;
30  constexpr int logPdfNbin = 150;
31  constexpr double logPdfBinsize = (logPdfMax - logPdfMin) / (double)logPdfNbin;
32  constexpr double logQScaleMax = 4.0;
33  constexpr double logQScaleMin = -1.0;
34  constexpr int logQScaleNbin = 500;
35  constexpr double logQScaleBinsize = (logQScaleMax - logQScaleMin) / (double)logQScaleNbin;
36 } // namespace
37 
40  DQMHelper dqm(&i);
41  i.setCurrentFolder("Generator/Particles");
42 
43  // Number of analyzed events
44  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1.);
45 
48  // quarks
49  particles.push_back(ParticleMonitor("u", 1, i));
50  particles.push_back(ParticleMonitor("ubar", -1, i));
51  particles.push_back(ParticleMonitor("d", 2, i));
52  particles.push_back(ParticleMonitor("dbar", -2, i));
53  particles.push_back(ParticleMonitor("s", 3, i));
54  particles.push_back(ParticleMonitor("sbar", -3, i));
55  particles.push_back(ParticleMonitor("c", 4, i));
56  particles.push_back(ParticleMonitor("cbar", -4, i));
57  particles.push_back(ParticleMonitor("b", 5, i));
58  particles.push_back(ParticleMonitor("bbar", -5, i));
59  particles.push_back(ParticleMonitor("t", 6, i));
60  particles.push_back(ParticleMonitor("tbar", -6, i));
61 
62  //leptons
63  particles.push_back(ParticleMonitor("eminus", 11, i));
64  particles.push_back(ParticleMonitor("eplus", -11, i));
65  particles.push_back(ParticleMonitor("nue", 12, i));
66  particles.push_back(ParticleMonitor("nuebar", -12, i));
67  particles.push_back(ParticleMonitor("muminus", 13, i));
68  particles.push_back(ParticleMonitor("muplus", -13, i));
69  particles.push_back(ParticleMonitor("numu", 14, i));
70  particles.push_back(ParticleMonitor("numubar", -14, i));
71  particles.push_back(ParticleMonitor("tauminus", 15, i));
72  particles.push_back(ParticleMonitor("tauplus", -15, i));
73  particles.push_back(ParticleMonitor("nutau", 16, i));
74  particles.push_back(ParticleMonitor("nutaubar", -16, i));
75 
76  //bosons
77  particles.push_back(ParticleMonitor("Wplus", 24, i));
78  particles.push_back(ParticleMonitor("Wminus", -24, i));
79  particles.push_back(ParticleMonitor("Z", 23, i));
80  particles.push_back(ParticleMonitor("gamma", 22, i));
81  particles.push_back(ParticleMonitor("gluon", 21, i));
82 
83  //mesons
84  particles.push_back(ParticleMonitor("piplus", 211, i, true)); //log
85  particles.push_back(ParticleMonitor("piminus", -211, i, true)); //log
86  particles.push_back(ParticleMonitor("pizero", 111, i, true)); //log
87  particles.push_back(ParticleMonitor("Kplus", 321, i));
88  particles.push_back(ParticleMonitor("Kminus", -321, i));
89  particles.push_back(ParticleMonitor("Klzero", 130, i));
90  particles.push_back(ParticleMonitor("Kszero", 310, i));
91 
92  //baryons
93  particles.push_back(ParticleMonitor("p", 2212, i, true)); //log
94  particles.push_back(ParticleMonitor("pbar", -2212, i, true)); //log
95  particles.push_back(ParticleMonitor("n", 2112, i, true)); //log
96  particles.push_back(ParticleMonitor("nbar", -2112, i, true)); //log
97  particles.push_back(ParticleMonitor("lambda0", 3122, i));
98  particles.push_back(ParticleMonitor("lambda0bar", -3122, i));
99 
100  //D mesons
101  particles.push_back(ParticleMonitor("Dplus", 411, i));
102  particles.push_back(ParticleMonitor("Dminus", -411, i));
103  particles.push_back(ParticleMonitor("Dzero", 421, i));
104  particles.push_back(ParticleMonitor("Dzerobar", -421, i));
105 
106  //B mesons
107  particles.push_back(ParticleMonitor("Bplus", 521, i));
108  particles.push_back(ParticleMonitor("Bminus", -521, i));
109  particles.push_back(ParticleMonitor("Bzero", 511, i));
110  particles.push_back(ParticleMonitor("Bzerobar", -511, i));
111  particles.push_back(ParticleMonitor("Bszero", 531, i));
112  particles.push_back(ParticleMonitor("Bszerobar", -531, i));
113 
114  //
116  "otherPtclNumber", "Log10(No. other ptcls)", 60, -1, 5, "log_{10}(No. other ptcls)", "Number of Events"); //Log
117  otherPtclMomentum = dqm.book1dHisto("otherPtclMomentum",
118  "Log10(p) other ptcls",
119  60,
120  -2,
121  4,
122  "log10(P^{other ptcls}) (log_{10}(GeV))",
123  "Number of Events");
124 
127  "genPtclNumber", "Log10(No. all particles)", 60, -1, 5, "log10(No. all particles)", "Number of Events"); //Log
129  "genVrtxNumber", "Log10(No. all vertexs)", 60, -1, 5, "log10(No. all vertexs)", "Number of Events"); //Log
130  //
131  stablePtclNumber = dqm.book1dHisto("stablePtclNumber",
132  "Log10(No. stable particles)",
133  50,
134  0,
135  5,
136  "log10(No. stable particles)",
137  "Number of Events"); //Log
139  "stablePtclPhi", "stable Ptcl Phi", 360, -180, 180, "#phi^{stable Ptcl} (rad)", "Number of Events");
141  "stablePtclEta", "stable Ptcl Eta (pseudo rapidity)", 220, -11, 11, "#eta^{stable Ptcl}", "Number of Events");
143  dqm.book1dHisto("stablePtclCharge", "stablePtclCharge", 5, -2, 2, "Charge^{stable ptcls}", "Number of Events");
144  stableChaNumber = dqm.book1dHisto("stableChaNumber",
145  "Log10(No. stable charged particles)",
146  50,
147  0,
148  5,
149  "log_{10}(No. stable charged particles)",
150  "Number of Events"); //Log
151  stablePtclp = dqm.book1dHisto("stablePtclp",
152  "Log10(p) stable ptcl p",
153  80,
154  -4,
155  4,
156  "log_{10}(P^{stable ptcl}) (log_{10}(GeV))",
157  "Number of Events"); //Log
158  stablePtclpT = dqm.book1dHisto("stablePtclpT",
159  "Log10(pT) stable ptcl pT",
160  80,
161  -4,
162  4,
163  "log_{10}(P_{t}^{stable ptcl}) (log_{10}(GeV))",
164  "Number of Events"); //Log
165  partonNumber =
166  dqm.book1dHisto("partonNumber", "number of partons", 100, 0, 100, "number of partons", "Number of Events");
167  partonpT =
168  dqm.book1dHisto("partonpT", "Log10(pT) parton pT", 80, -4, 4, "Log10(P_{t}^{parton})", "Number of Events"); //Log
169  outVrtxStablePtclNumber = dqm.book1dHisto("outVrtxStablePtclNumber",
170  "No. outgoing stable ptcls from vrtx",
171  10,
172  0,
173  10,
174  "No. outgoing stable ptcls from vrtx",
175  "Number of Events");
176  //
177  outVrtxPtclNumber = dqm.book1dHisto("outVrtxPtclNumber",
178  "No. outgoing ptcls from vrtx",
179  30,
180  0,
181  30,
182  "No. outgoing ptcls from vrtx",
183  "Number of Events");
184  vrtxZ = dqm.book1dHisto("VrtxZ", "VrtxZ", 50, -250, 250, "Z_{Vtx}", "Number of Events");
185  vrtxRadius = dqm.book1dHisto("vrtxRadius", "vrtxRadius", 50, 0, 50, "R_{vtx}", "Number of Events");
186  //
187  unknownPDTNumber = dqm.book1dHisto("unknownPDTNumber",
188  "Log10(No. unknown ptcls PDT)",
189  60,
190  -1,
191  5,
192  "log_{10}(No. unknown ptcls PDT)",
193  "Number of Events"); //Log
194  genPtclStatus = dqm.book1dHisto("genPtclStatus", "Status of genParticle", 200, 0, 200., "", "Number of Events");
195  //
196  Bjorken_x = dqm.book1dHisto("Bjorken_x", "Bjorken_x", 1000, 0.0, 1.0, "Bjorken_{x}", "Number of Events");
197  pdf_u = dqm.book1dHisto(
198  "pdf_u", "Log10(PDF(u,x,Q))", logPdfNbin, logPdfMin, logPdfMax, "log_{10}(x*f(x))", "Number of Events"); //Log
199  pdf_ubar = dqm.book1dHisto("pdf_ubar",
200  "Log10(PDF(ubar,x,Q))",
201  logPdfNbin,
202  logPdfMin,
203  logPdfMax,
204  "log_{10}(x*f(x))",
205  "Number of Events"); //Log
206  pdf_d = dqm.book1dHisto(
207  "pdf_d", "Log10(PDF(d,x,Q))", logPdfNbin, logPdfMin, logPdfMax, "log_{10}(x*f(x))", "Number of Events"); //Log
208  pdf_dbar = dqm.book1dHisto("pdf_dbar",
209  "Log10(PDF(dbar,x,Q))",
210  logPdfNbin,
211  logPdfMin,
212  logPdfMax,
213  "log_{10}(x*f(x))",
214  "Number of Events"); //Log
215  pdf_ssbar = dqm.book1dHisto("pdf_ssbar",
216  "Log10(PDF(ssbar,x,Q))",
217  logPdfNbin,
218  logPdfMin,
219  logPdfMax,
220  "log_{10}(x*f(x))",
221  "Number of Events"); //Log
222  pdf_ccbar = dqm.book1dHisto("pdf_ccbar",
223  "Log10(PDF(ccbar,x,Q))",
224  logPdfNbin,
225  logPdfMin,
226  logPdfMax,
227  "log_{10}(x*f(x))",
228  "Number of Events"); //Log
229  pdf_bbbar = dqm.book1dHisto("pdf_bbbar",
230  "Log10(PDF(bbbar,x,Q))",
231  logPdfNbin,
232  logPdfMin,
233  logPdfMax,
234  "log_{10}(x*f(x))",
235  "Number of Events"); //Log
236  pdf_g = dqm.book1dHisto(
237  "pdf_g", "Log10(PDF(g,x,Q))", logPdfNbin, logPdfMin, logPdfMax, "log_{10}(x*f(x))", "Number of Events"); //Log
238  scalePDF = dqm.book1dHisto("scalePDF",
239  "Log10(Q-scale(GeV))",
240  logQScaleNbin,
241  logQScaleMin,
242  logQScaleMax,
243  "log_{10}(Q-scale(GeV))",
244  "Number of Events");
245  parton1Id = dqm.book1dHisto("parton1Id", "ID of parton 1", 45, -14.5, 30.5, "ID", "Number of Events");
246  parton2Id = dqm.book1dHisto("parton2Id", "ID of parton 2", 45, -14.5, 30.5, "ID", "Number of Events");
247  //
248  status1ShortLived = dqm.book1dHisto("status1ShortLived", "Status 1 short lived", 11, 0, 11, "", "Number of Events");
249  status1ShortLived->setBinLabel(1, "d/dbar");
250  status1ShortLived->setBinLabel(2, "u/ubar");
251  status1ShortLived->setBinLabel(3, "s/sbar");
252  status1ShortLived->setBinLabel(4, "c/cbar");
253  status1ShortLived->setBinLabel(5, "b/bbar");
254  status1ShortLived->setBinLabel(6, "t/tbar");
256  status1ShortLived->setBinLabel(8, "tau-/tau+");
257  status1ShortLived->setBinLabel(9, "Z0");
258  status1ShortLived->setBinLabel(10, "W-/W+");
259  status1ShortLived->setBinLabel(11, "PDG = 7,8,17,25-99");
260 
261  log10DeltaEcms = dqm.book1dHisto("DeltaEcms1log10",
262  "log_{10} of deviation from nominal Ecms",
263  200,
264  -1.,
265  5.,
266  "log_{10}(#DeltaE) (log_{10}(GeV))",
267  "Number of Events");
268  DeltaEcms =
269  dqm.book1dHisto("DeltaEcms1", "deviation from nominal Ecms", 200, -1., 1., "#DeltaE (GeV)", "Number of Events");
270  DeltaPx =
271  dqm.book1dHisto("DeltaPx1", "deviation from nominal Px", 200, -1., 1., "#DeltaP_{x} (GeV)", "Number of Events");
272  DeltaPy =
273  dqm.book1dHisto("DeltaPy1", "deviation from nominal Py", 200, -1., 1., "#DeltaP_{y} (GeV)", "Number of Events");
274  DeltaPz =
275  dqm.book1dHisto("DeltaPz1", "deviation from nominal Pz", 200, -1., 1., "#DeltaP_{z} (GeV)", "Number of Events");
276  return;
277 }
278 
281  int partonNum = 0;
282  //
283  int outVrtxStablePtclNum = 0;
284  int stablePtclNum = 0;
285  int otherPtclNum = 0;
286  int unknownPDTNum = 0;
287  int stableChaNum = 0;
288  //
289  double bjorken = 0.;
290  double logPdf1 = 0.;
291  double logPdf2 = 0.;
292  double logQScale = 0.;
293  //
294  double etotal = 0.;
295  double pxtotal = 0.;
296  double pytotal = 0.;
297  double pztotal = 0.;
298 
301  iEvent.getByToken(hepmcCollectionToken_, evt);
302 
303  //Get EVENT
304  HepMC::GenEvent const *myGenEvent = evt->GetEvent();
305 
306  double weight = wmanager_.weight(iEvent);
307 
308  nEvt->Fill(0.5, weight);
309 
310  genPtclNumber->Fill(log10(myGenEvent->particles_size()), weight);
311  genVrtxNumber->Fill(log10(myGenEvent->vertices_size()), weight);
312 
314  HepMC::PdfInfo const *pdf = myGenEvent->pdf_info();
315  if (pdf) {
316  bjorken = ((pdf->x1()) / ((pdf->x1()) + (pdf->x2())));
317  logQScale = log10(pdf->scalePDF());
318  if (logQScale > logQScaleMax)
319  logQScale = logQScaleMax - 0.5 * logQScaleBinsize; // visualize overflow & underflow in the histograms
320  if (logQScale < logQScaleMin)
321  logQScale = logQScaleMin + 0.5 * logQScaleBinsize;
322  logPdf1 = log10(pdf->pdf1());
323  if (logPdf1 > logPdfMax)
324  logPdf1 = logPdfMax - 0.5 * logPdfBinsize;
325  if (logPdf1 < logPdfMin)
326  logPdf1 = logPdfMin + 0.5 * logPdfBinsize;
327  logPdf2 = log10(pdf->pdf2());
328  if (logPdf2 > logPdfMax)
329  logPdf2 = logPdfMax - 0.5 * logPdfBinsize;
330  if (logPdf2 < logPdfMin)
331  logPdf2 = logPdfMin + 0.5 * logPdfBinsize;
332  Bjorken_x->Fill(bjorken, weight);
333  scalePDF->Fill(logQScale, weight);
334  parton1Id->Fill((double)pdf->id1(), weight);
335  parton2Id->Fill((double)pdf->id2(), weight);
336  if (pdf->id1() == 2)
337  pdf_u->Fill(logPdf1, weight);
338  if (pdf->id2() == 2)
339  pdf_u->Fill(logPdf2, weight);
340  if (pdf->id1() == -2)
341  pdf_ubar->Fill(logPdf1, weight);
342  if (pdf->id2() == -2)
343  pdf_ubar->Fill(logPdf2, weight);
344  if (pdf->id1() == 1)
345  pdf_d->Fill(logPdf1, weight);
346  if (pdf->id2() == 1)
347  pdf_d->Fill(logPdf2, weight);
348  if (pdf->id1() == -1)
349  pdf_dbar->Fill(logPdf1, weight);
350  if (pdf->id2() == -1)
351  pdf_dbar->Fill(logPdf2, weight);
352  if (std::abs(pdf->id1()) == 3)
353  pdf_ssbar->Fill(logPdf1, weight);
354  if (std::abs(pdf->id2()) == 3)
355  pdf_ssbar->Fill(logPdf2, weight);
356  if (std::abs(pdf->id1()) == 4)
357  pdf_ccbar->Fill(logPdf1, weight);
358  if (std::abs(pdf->id2()) == 4)
359  pdf_ccbar->Fill(logPdf2, weight);
360  if (std::abs(pdf->id1()) == 5)
361  pdf_bbbar->Fill(logPdf1, weight);
362  if (std::abs(pdf->id2()) == 5)
363  pdf_bbbar->Fill(logPdf2, weight);
364  if (std::abs(pdf->id1()) == 21)
365  pdf_g->Fill(logPdf1, weight);
366  if (std::abs(pdf->id2()) == 21)
367  pdf_g->Fill(logPdf2, weight);
368  }
369 
370  //Looping through the VERTICES in the event
371  HepMC::GenEvent::vertex_const_iterator vrtxBegin = myGenEvent->vertices_begin();
372  HepMC::GenEvent::vertex_const_iterator vrtxEnd = myGenEvent->vertices_end();
373  unsigned int nvtx(0);
374  for (HepMC::GenEvent::vertex_const_iterator vrtxIt = vrtxBegin; vrtxIt != vrtxEnd; ++vrtxIt) {
376  HepMC::GenVertex const *vrtx = *vrtxIt;
377  outVrtxPtclNumber->Fill(vrtx->particles_out_size(), weight);
378  //std::cout << "all " << vrtx->particles_out_size() << '\n';
379 
380  if (nvtx == 0) {
381  vrtxZ->Fill(vrtx->point3d().z(), weight);
382  vrtxRadius->Fill(vrtx->point3d().perp(), weight);
383  }
385  HepMC::GenVertex::particles_out_const_iterator vrtxPtclBegin = vrtx->particles_out_const_begin();
386  HepMC::GenVertex::particles_out_const_iterator vrtxPtclEnd = vrtx->particles_out_const_end();
387  outVrtxStablePtclNum = 0;
388  for (HepMC::GenVertex::particles_out_const_iterator vrtxPtclIt = vrtxPtclBegin; vrtxPtclIt != vrtxPtclEnd;
389  ++vrtxPtclIt) {
390  HepMC::GenParticle const *vrtxPtcl = *vrtxPtclIt;
391  if (vrtxPtcl->status() == 1) {
392  ++outVrtxStablePtclNum;
393  //std::cout << "stable " << outVrtxStablePtclNum << '\n';
394  }
395  }
396  outVrtxStablePtclNumber->Fill(outVrtxStablePtclNum, weight);
397  nvtx++;
398  } //vertices
399 
401  HepMC::GenEvent::particle_const_iterator ptclBegin = myGenEvent->particles_begin();
402  HepMC::GenEvent::particle_const_iterator ptclEnd = myGenEvent->particles_end();
403  for (HepMC::GenEvent::particle_const_iterator ptclIt = ptclBegin; ptclIt != ptclEnd; ++ptclIt) {
405  HepMC::GenParticle const *ptcl = *ptclIt;
406  int Id = ptcl->pdg_id(); // std::cout << Id << '\n';
407  float Log_p = log10(ptcl->momentum().rho());
408  double charge = 999.; // for the charge it's needed a HepPDT method
409  int status = ptcl->status();
410  const HepPDT::ParticleData *PData = fPDGTable->particle(HepPDT::ParticleID(Id));
411  if (PData == nullptr) {
412  // std::cout << "Unknown id = " << Id << '\n';
413  ++unknownPDTNum;
414  } else
415  charge = PData->charge();
416 
418  genPtclStatus->Fill((float)status, weight);
419 
421  if (ptcl->status() == 1) {
422  ++stablePtclNum;
423  stablePtclPhi->Fill(ptcl->momentum().phi() / CLHEP::degree,
424  weight); //std::cout << ptcl->polarization().phi() << '\n';
425  stablePtclEta->Fill(ptcl->momentum().pseudoRapidity(), weight);
426  stablePtclCharge->Fill(charge, weight); // std::cout << ptclData.charge() << '\n';
427  stablePtclp->Fill(Log_p, weight);
428  stablePtclpT->Fill(log10(ptcl->momentum().perp()), weight);
429  if (charge != 0. && charge != 999.)
430  ++stableChaNum;
431  if (std::abs(Id) == 1)
432  status1ShortLived->Fill(1, weight);
433  if (std::abs(Id) == 2)
434  status1ShortLived->Fill(2, weight);
435  if (std::abs(Id) == 3)
436  status1ShortLived->Fill(3, weight);
437  if (std::abs(Id) == 4)
438  status1ShortLived->Fill(4, weight);
439  if (std::abs(Id) == 5)
440  status1ShortLived->Fill(5, weight);
441  if (std::abs(Id) == 6)
442  status1ShortLived->Fill(6, weight);
443  if (Id == 21)
444  status1ShortLived->Fill(7, weight);
445  if (std::abs(Id) == 15)
446  status1ShortLived->Fill(8, weight);
447  if (Id == 23)
448  status1ShortLived->Fill(9, weight);
449  if (std::abs(Id) == 24)
450  status1ShortLived->Fill(10, weight);
451  if (std::abs(Id) == 7 || std::abs(Id) == 8 || std::abs(Id) == 17 || (std::abs(Id) >= 25 && std::abs(Id) <= 99))
452  status1ShortLived->Fill(11, weight);
453  etotal += ptcl->momentum().e();
454  pxtotal += ptcl->momentum().px();
455  pytotal += ptcl->momentum().py();
456  pztotal += ptcl->momentum().pz();
457  }
458 
459  if (abs(Id) < 6 || abs(Id) == 22) {
460  ++partonNum;
461  partonpT->Fill(Log_p, weight);
462  }
463 
464  bool indentified = false;
465  for (unsigned int i = 0; i < particles.size(); i++) {
466  if (particles.at(i).Fill(ptcl, weight)) {
467  indentified = true;
468  break;
469  }
470  }
471  if (!indentified) {
472  ++otherPtclNum;
473  otherPtclMomentum->Fill(Log_p, weight);
474  }
475  } //event particles
476 
477  // set a default sqrt(s) and then check in the event
478  double ecms = 13000.;
479  if (myGenEvent->valid_beam_particles()) {
480  ecms = myGenEvent->beam_particles().first->momentum().e() + myGenEvent->beam_particles().second->momentum().e();
481  }
482  log10DeltaEcms->Fill(log10(fabs(etotal - ecms)), weight);
483  DeltaEcms->Fill(etotal - ecms, weight);
484  DeltaPx->Fill(pxtotal, weight);
485  DeltaPy->Fill(pytotal, weight);
486  DeltaPz->Fill(pztotal, weight);
487 
489  stablePtclNumber->Fill(log10(stablePtclNum + 0.1), weight);
490  stableChaNumber->Fill(log10(stableChaNum + 0.1), weight);
491  otherPtclNumber->Fill(log10(otherPtclNum + 0.1), weight);
492  unknownPDTNumber->Fill(log10(unknownPDTNum + 0.1), weight);
493  //
494  partonNumber->Fill(partonNum, weight);
495  for (unsigned int i = 0; i < particles.size(); i++) {
496  particles.at(i).FillCount(weight);
497  };
498 
499 } //analyze
MonitorElement * genPtclStatus
MonitorElement * DeltaEcms
MonitorElement * genVrtxNumber
const edm::EventSetup & c
MonitorElement * unknownPDTNumber
MonitorElement * pdf_bbbar
MonitorElement * outVrtxPtclNumber
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
MonitorElement * pdf_ccbar
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
MonitorElement * Bjorken_x
MonitorElement * otherPtclNumber
other ME&#39;s
void dqmBeginRun(const edm::Run &r, const edm::EventSetup &c) override
BasicHepMCValidation(const edm::ParameterSet &)
MonitorElement * parton2Id
list status
Definition: mps_update.py:107
MonitorElement * otherPtclMomentum
MonitorElement * stablePtclCharge
MonitorElement * partonNumber
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
MonitorElement * vrtxRadius
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
MonitorElement * log10DeltaEcms
MonitorElement * stablePtclPhi
MonitorElement * stablePtclEta
virtual 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)
HepPDT::ParticleData ParticleData
void analyze(edm::Event const &, edm::EventSetup const &) override
MonitorElement * stablePtclp
MonitorElement * book1dHisto(const std::string &name, const std::string &title, int n, double xmin, double xmax, const std::string &xaxis, const std::string &yaxis)
Definition: DQMHelper.cc:7
MonitorElement * status1ShortLived
MonitorElement * stablePtclpT
MonitorElement * stableChaNumber
MonitorElement * stablePtclNumber
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > fPDGTableToken
MonitorElement * genPtclNumber
edm::InputTag hepmcCollection_
MonitorElement * outVrtxStablePtclNumber
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
int weight
Definition: histoStyle.py:51
double weight(const edm::Event &)
std::vector< ParticleMonitor > particles
MonitorElement * pdf_ssbar
Definition: Run.h:45
MonitorElement * parton1Id