CMS 3D CMS Logo

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

#include <BasicHepMCValidation.h>

Inheritance diagram for BasicHepMCValidation:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

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 DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Attributes

MonitorElementBjorken_x
 
MonitorElementDeltaEcms
 
MonitorElementDeltaPx
 
MonitorElementDeltaPy
 
MonitorElementDeltaPz
 
edm::ESHandle< HepPDT::ParticleDataTablefPDGTable
 PDT table. More...
 
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecordfPDGTableToken
 
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_
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

Definition at line 34 of file BasicHepMCValidation.h.

Constructor & Destructor Documentation

◆ BasicHepMCValidation()

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

Definition at line 14 of file BasicHepMCValidation.cc.

References fPDGTableToken, hepmcCollection_, and hepmcCollectionToken_.

15  : wmanager_(iPSet, consumesCollector()), hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")) {
16  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
17  fPDGTableToken = esConsumes<edm::Transition::BeginRun>();
18 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > fPDGTableToken
edm::InputTag hepmcCollection_

◆ ~BasicHepMCValidation()

BasicHepMCValidation::~BasicHepMCValidation ( )
override

Definition at line 20 of file BasicHepMCValidation.cc.

20 {}

Member Function Documentation

◆ analyze()

void BasicHepMCValidation::analyze ( edm::Event const &  iEvent,
edm::EventSetup const &  iSetup 
)
overridevirtual

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

Reimplemented from DQMEDAnalyzer.

Definition at line 279 of file BasicHepMCValidation.cc.

References funct::abs(), Bjorken_x, ALCARECOTkAlJpsiMuMu_cff::charge, DeltaEcms, DeltaPx, DeltaPy, DeltaPz, dqm::impl::MonitorElement::Fill(), fPDGTable, GenParticle::GenParticle, genPtclNumber, genPtclStatus, genVrtxNumber, edm::HepMCProduct::GetEvent(), hepmcCollectionToken_, mps_fire::i, iEvent, log10DeltaEcms, nEvt, otherPtclMomentum, otherPtclNumber, outVrtxPtclNumber, outVrtxStablePtclNumber, LHEGenericFilter_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_.

279  {
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)
433  if (std::abs(Id) == 2)
435  if (std::abs(Id) == 3)
437  if (std::abs(Id) == 4)
439  if (std::abs(Id) == 5)
441  if (std::abs(Id) == 6)
443  if (Id == 21)
445  if (std::abs(Id) == 15)
447  if (Id == 23)
449  if (std::abs(Id) == 24)
451  if (std::abs(Id) == 7 || std::abs(Id) == 8 || std::abs(Id) == 17 || (std::abs(Id) >= 25 && std::abs(Id) <= 99))
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
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)
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
HepPDT::ParticleData ParticleData
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
MonitorElement * stablePtclp
MonitorElement * status1ShortLived
MonitorElement * stablePtclpT
MonitorElement * stableChaNumber
MonitorElement * stablePtclNumber
MonitorElement * genPtclNumber
MonitorElement * outVrtxStablePtclNumber
double weight(const edm::Event &)
std::vector< ParticleMonitor > particles
MonitorElement * pdf_ssbar
MonitorElement * parton1Id

◆ bookHistograms()

void BasicHepMCValidation::bookHistograms ( DQMStore::IBooker i,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Setting the DQM top directories

Booking the ME's multiplicity
other

Implements DQMEDAnalyzer.

Definition at line 38 of file BasicHepMCValidation.cc.

References Bjorken_x, DeltaEcms, DeltaPx, DeltaPy, DeltaPz, genPtclNumber, genPtclStatus, genVrtxNumber, mps_fire::i, log10DeltaEcms, 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, dqm::impl::MonitorElement::setBinLabel(), stableChaNumber, stablePtclCharge, stablePtclEta, stablePtclNumber, stablePtclp, stablePtclPhi, stablePtclpT, status1ShortLived, unknownPDTNumber, vrtxRadius, and vrtxZ.

38  {
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  //
115  otherPtclNumber = dqm.book1dHisto(
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 
126  genPtclNumber = dqm.book1dHisto(
127  "genPtclNumber", "Log10(No. all particles)", 60, -1, 5, "log10(No. all particles)", "Number of Events"); //Log
128  genVrtxNumber = dqm.book1dHisto(
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
138  stablePtclPhi = dqm.book1dHisto(
139  "stablePtclPhi", "stable Ptcl Phi", 360, -180, 180, "#phi^{stable Ptcl} (rad)", "Number of Events");
140  stablePtclEta = dqm.book1dHisto(
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 }
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
MonitorElement * stablePtclCharge
MonitorElement * partonNumber
MonitorElement * vrtxRadius
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)
MonitorElement * stablePtclp
MonitorElement * status1ShortLived
MonitorElement * stablePtclpT
MonitorElement * stableChaNumber
MonitorElement * stablePtclNumber
MonitorElement * genPtclNumber
MonitorElement * outVrtxStablePtclNumber
Definition: DQMStore.h:18
std::vector< ParticleMonitor > particles
MonitorElement * pdf_ssbar
MonitorElement * parton1Id

◆ dqmBeginRun()

void BasicHepMCValidation::dqmBeginRun ( const edm::Run r,
const edm::EventSetup c 
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 22 of file BasicHepMCValidation.cc.

References HltBtagPostValidation_cff::c, fPDGTable, and fPDGTableToken.

22  {
23  fPDGTable = c.getHandle(fPDGTableToken);
24 }
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
PDT table.
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > fPDGTableToken

Member Data Documentation

◆ Bjorken_x

MonitorElement* BasicHepMCValidation::Bjorken_x
private

Definition at line 221 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ DeltaEcms

MonitorElement* BasicHepMCValidation::DeltaEcms
private

Definition at line 237 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ DeltaPx

MonitorElement* BasicHepMCValidation::DeltaPx
private

Definition at line 238 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ DeltaPy

MonitorElement* BasicHepMCValidation::DeltaPy
private

Definition at line 239 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ DeltaPz

MonitorElement* BasicHepMCValidation::DeltaPz
private

Definition at line 240 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ fPDGTable

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

PDT table.

Definition at line 48 of file BasicHepMCValidation.h.

Referenced by analyze(), and dqmBeginRun().

◆ fPDGTableToken

edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> BasicHepMCValidation::fPDGTableToken
private

Definition at line 49 of file BasicHepMCValidation.h.

Referenced by BasicHepMCValidation(), and dqmBeginRun().

◆ genPtclNumber

MonitorElement* BasicHepMCValidation::genPtclNumber
private

Definition at line 201 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genPtclStatus

MonitorElement* BasicHepMCValidation::genPtclStatus
private

Definition at line 205 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genVrtxNumber

MonitorElement* BasicHepMCValidation::genVrtxNumber
private

Definition at line 202 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ hepmcCollection_

edm::InputTag BasicHepMCValidation::hepmcCollection_
private

Definition at line 45 of file BasicHepMCValidation.h.

Referenced by BasicHepMCValidation().

◆ hepmcCollectionToken_

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

Definition at line 242 of file BasicHepMCValidation.h.

Referenced by analyze(), and BasicHepMCValidation().

◆ log10DeltaEcms

MonitorElement* BasicHepMCValidation::log10DeltaEcms
private

Definition at line 236 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ nEvt

MonitorElement* BasicHepMCValidation::nEvt
private

Definition at line 195 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ otherPtclMomentum

MonitorElement* BasicHepMCValidation::otherPtclMomentum
private

Definition at line 200 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ otherPtclNumber

MonitorElement* BasicHepMCValidation::otherPtclNumber
private

other ME's

Definition at line 199 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ outVrtxPtclNumber

MonitorElement* BasicHepMCValidation::outVrtxPtclNumber
private

Definition at line 204 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ outVrtxStablePtclNumber

MonitorElement* BasicHepMCValidation::outVrtxStablePtclNumber
private

Definition at line 216 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ particles

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

Definition at line 196 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ parton1Id

MonitorElement* BasicHepMCValidation::parton1Id
private

Definition at line 231 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ parton2Id

MonitorElement* BasicHepMCValidation::parton2Id
private

Definition at line 232 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ partonNumber

MonitorElement* BasicHepMCValidation::partonNumber
private

Definition at line 214 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ partonpT

MonitorElement* BasicHepMCValidation::partonpT
private

Definition at line 215 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_bbbar

MonitorElement* BasicHepMCValidation::pdf_bbbar
private

Definition at line 228 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_ccbar

MonitorElement* BasicHepMCValidation::pdf_ccbar
private

Definition at line 227 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_d

MonitorElement* BasicHepMCValidation::pdf_d
private

Definition at line 224 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_dbar

MonitorElement* BasicHepMCValidation::pdf_dbar
private

Definition at line 225 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_g

MonitorElement* BasicHepMCValidation::pdf_g
private

Definition at line 229 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_ssbar

MonitorElement* BasicHepMCValidation::pdf_ssbar
private

Definition at line 226 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_u

MonitorElement* BasicHepMCValidation::pdf_u
private

Definition at line 222 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ pdf_ubar

MonitorElement* BasicHepMCValidation::pdf_ubar
private

Definition at line 223 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ scalePDF

MonitorElement* BasicHepMCValidation::scalePDF
private

Definition at line 230 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ stableChaNumber

MonitorElement* BasicHepMCValidation::stableChaNumber
private

Definition at line 208 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ stablePtclCharge

MonitorElement* BasicHepMCValidation::stablePtclCharge
private

Definition at line 211 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ stablePtclEta

MonitorElement* BasicHepMCValidation::stablePtclEta
private

Definition at line 210 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ stablePtclNumber

MonitorElement* BasicHepMCValidation::stablePtclNumber
private

Definition at line 207 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ stablePtclp

MonitorElement* BasicHepMCValidation::stablePtclp
private

Definition at line 212 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ stablePtclPhi

MonitorElement* BasicHepMCValidation::stablePtclPhi
private

Definition at line 209 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ stablePtclpT

MonitorElement* BasicHepMCValidation::stablePtclpT
private

Definition at line 213 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ status1ShortLived

MonitorElement* BasicHepMCValidation::status1ShortLived
private

Definition at line 234 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ unknownPDTNumber

MonitorElement* BasicHepMCValidation::unknownPDTNumber
private

Definition at line 203 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ vrtxRadius

MonitorElement* BasicHepMCValidation::vrtxRadius
private

Definition at line 219 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ vrtxZ

MonitorElement* BasicHepMCValidation::vrtxZ
private

Definition at line 218 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

◆ wmanager_

WeightManager BasicHepMCValidation::wmanager_
private

Definition at line 44 of file BasicHepMCValidation.h.

Referenced by analyze().