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::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::EndLuminosityBlockProducer, edm::one::WatchLuminosityBlocks, edm::Accumulator > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

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
 
virtual void dqmBeginLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &)
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
virtual void dqmEndLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &)
 
virtual void dqmEndRun (edm::Run const &, edm::EventSetup const &)
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) final
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &, edm::EventSetup const &) final
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::one::EDProducer< edm::EndRunProducer, edm::one::WatchRuns, edm::EndLuminosityBlockProducer, edm::one::WatchLuminosityBlocks, edm::Accumulator >
 EDProducer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

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_
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

Definition at line 34 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_.

15  : wmanager_(iPSet, consumesCollector()), hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")) {
16  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
17 }
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
edm::InputTag hepmcCollection_
BasicHepMCValidation::~BasicHepMCValidation ( )
override

Definition at line 19 of file BasicHepMCValidation.cc.

19 {}

Member Function Documentation

void BasicHepMCValidation::analyze ( edm::Event const &  ,
edm::EventSetup const &   
)
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 274 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::Event::getByToken(), edm::HepMCProduct::GetEvent(), hepmcCollectionToken_, mps_fire::i, log10DeltaEcms, logPdfBinsize, logPdfMax, logPdfMin, logQScaleBinsize, 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_.

274  {
276  int partonNum = 0;
277  //
278  int outVrtxStablePtclNum = 0;
279  int stablePtclNum = 0;
280  int otherPtclNum = 0;
281  int unknownPDTNum = 0;
282  int stableChaNum = 0;
283  //
284  double bjorken = 0.;
285  double logPdf1 = 0.;
286  double logPdf2 = 0.;
287  double logQScale = 0.;
288  //
289  double etotal = 0.;
290  double pxtotal = 0.;
291  double pytotal = 0.;
292  double pztotal = 0.;
293 
296  iEvent.getByToken(hepmcCollectionToken_, evt);
297 
298  //Get EVENT
299  HepMC::GenEvent const *myGenEvent = evt->GetEvent();
300 
301  double weight = wmanager_.weight(iEvent);
302 
303  nEvt->Fill(0.5, weight);
304 
305  genPtclNumber->Fill(log10(myGenEvent->particles_size()), weight);
306  genVrtxNumber->Fill(log10(myGenEvent->vertices_size()), weight);
307 
309  HepMC::PdfInfo const *pdf = myGenEvent->pdf_info();
310  if (pdf) {
311  bjorken = ((pdf->x1()) / ((pdf->x1()) + (pdf->x2())));
312  logQScale = log10(pdf->scalePDF());
313  if (logQScale > logQScaleMax)
314  logQScale = logQScaleMax - 0.5 * logQScaleBinsize; // visualize overflow & underflow in the histograms
315  if (logQScale < logQScaleMin)
316  logQScale = logQScaleMin + 0.5 * logQScaleBinsize;
317  logPdf1 = log10(pdf->pdf1());
318  if (logPdf1 > logPdfMax)
319  logPdf1 = logPdfMax - 0.5 * logPdfBinsize;
320  if (logPdf1 < logPdfMin)
321  logPdf1 = logPdfMin + 0.5 * logPdfBinsize;
322  logPdf2 = log10(pdf->pdf2());
323  if (logPdf2 > logPdfMax)
324  logPdf2 = logPdfMax - 0.5 * logPdfBinsize;
325  if (logPdf2 < logPdfMin)
326  logPdf2 = logPdfMin + 0.5 * logPdfBinsize;
327  Bjorken_x->Fill(bjorken, weight);
328  scalePDF->Fill(logQScale, weight);
329  parton1Id->Fill((double)pdf->id1(), weight);
330  parton2Id->Fill((double)pdf->id2(), weight);
331  if (pdf->id1() == 2)
332  pdf_u->Fill(logPdf1, weight);
333  if (pdf->id2() == 2)
334  pdf_u->Fill(logPdf2, weight);
335  if (pdf->id1() == -2)
336  pdf_ubar->Fill(logPdf1, weight);
337  if (pdf->id2() == -2)
338  pdf_ubar->Fill(logPdf2, weight);
339  if (pdf->id1() == 1)
340  pdf_d->Fill(logPdf1, weight);
341  if (pdf->id2() == 1)
342  pdf_d->Fill(logPdf2, weight);
343  if (pdf->id1() == -1)
344  pdf_dbar->Fill(logPdf1, weight);
345  if (pdf->id2() == -1)
346  pdf_dbar->Fill(logPdf2, weight);
347  if (std::abs(pdf->id1()) == 3)
348  pdf_ssbar->Fill(logPdf1, weight);
349  if (std::abs(pdf->id2()) == 3)
350  pdf_ssbar->Fill(logPdf2, weight);
351  if (std::abs(pdf->id1()) == 4)
352  pdf_ccbar->Fill(logPdf1, weight);
353  if (std::abs(pdf->id2()) == 4)
354  pdf_ccbar->Fill(logPdf2, weight);
355  if (std::abs(pdf->id1()) == 5)
356  pdf_bbbar->Fill(logPdf1, weight);
357  if (std::abs(pdf->id2()) == 5)
358  pdf_bbbar->Fill(logPdf2, weight);
359  if (std::abs(pdf->id1()) == 21)
360  pdf_g->Fill(logPdf1, weight);
361  if (std::abs(pdf->id2()) == 21)
362  pdf_g->Fill(logPdf2, weight);
363  }
364 
365  //Looping through the VERTICES in the event
366  HepMC::GenEvent::vertex_const_iterator vrtxBegin = myGenEvent->vertices_begin();
367  HepMC::GenEvent::vertex_const_iterator vrtxEnd = myGenEvent->vertices_end();
368  unsigned int nvtx(0);
369  for (HepMC::GenEvent::vertex_const_iterator vrtxIt = vrtxBegin; vrtxIt != vrtxEnd; ++vrtxIt) {
371  HepMC::GenVertex const *vrtx = *vrtxIt;
372  outVrtxPtclNumber->Fill(vrtx->particles_out_size(), weight);
373  //std::cout << "all " << vrtx->particles_out_size() << '\n';
374 
375  if (nvtx == 0) {
376  vrtxZ->Fill(vrtx->point3d().z(), weight);
377  vrtxRadius->Fill(vrtx->point3d().perp(), weight);
378  }
380  HepMC::GenVertex::particles_out_const_iterator vrtxPtclBegin = vrtx->particles_out_const_begin();
381  HepMC::GenVertex::particles_out_const_iterator vrtxPtclEnd = vrtx->particles_out_const_end();
382  outVrtxStablePtclNum = 0;
383  for (HepMC::GenVertex::particles_out_const_iterator vrtxPtclIt = vrtxPtclBegin; vrtxPtclIt != vrtxPtclEnd;
384  ++vrtxPtclIt) {
385  HepMC::GenParticle const *vrtxPtcl = *vrtxPtclIt;
386  if (vrtxPtcl->status() == 1) {
387  ++outVrtxStablePtclNum;
388  //std::cout << "stable " << outVrtxStablePtclNum << '\n';
389  }
390  }
391  outVrtxStablePtclNumber->Fill(outVrtxStablePtclNum, weight);
392  nvtx++;
393  } //vertices
394 
396  HepMC::GenEvent::particle_const_iterator ptclBegin = myGenEvent->particles_begin();
397  HepMC::GenEvent::particle_const_iterator ptclEnd = myGenEvent->particles_end();
398  for (HepMC::GenEvent::particle_const_iterator ptclIt = ptclBegin; ptclIt != ptclEnd; ++ptclIt) {
400  HepMC::GenParticle const *ptcl = *ptclIt;
401  int Id = ptcl->pdg_id(); // std::cout << Id << '\n';
402  float Log_p = log10(ptcl->momentum().rho());
403  double charge = 999.; // for the charge it's needed a HepPDT method
404  int status = ptcl->status();
405  const HepPDT::ParticleData *PData = fPDGTable->particle(HepPDT::ParticleID(Id));
406  if (PData == nullptr) {
407  // std::cout << "Unknown id = " << Id << '\n';
408  ++unknownPDTNum;
409  } else
410  charge = PData->charge();
411 
413  genPtclStatus->Fill((float)status, weight);
414 
416  if (ptcl->status() == 1) {
417  ++stablePtclNum;
418  stablePtclPhi->Fill(ptcl->momentum().phi() / CLHEP::degree,
419  weight); //std::cout << ptcl->polarization().phi() << '\n';
420  stablePtclEta->Fill(ptcl->momentum().pseudoRapidity(), weight);
421  stablePtclCharge->Fill(charge, weight); // std::cout << ptclData.charge() << '\n';
422  stablePtclp->Fill(Log_p, weight);
423  stablePtclpT->Fill(log10(ptcl->momentum().perp()), weight);
424  if (charge != 0. && charge != 999.)
425  ++stableChaNum;
426  if (std::abs(Id) == 1)
427  status1ShortLived->Fill(1, weight);
428  if (std::abs(Id) == 2)
429  status1ShortLived->Fill(2, weight);
430  if (std::abs(Id) == 3)
431  status1ShortLived->Fill(3, weight);
432  if (std::abs(Id) == 4)
433  status1ShortLived->Fill(4, weight);
434  if (std::abs(Id) == 5)
435  status1ShortLived->Fill(5, weight);
436  if (std::abs(Id) == 6)
437  status1ShortLived->Fill(6, weight);
438  if (Id == 21)
439  status1ShortLived->Fill(7, weight);
440  if (std::abs(Id) == 15)
441  status1ShortLived->Fill(8, weight);
442  if (Id == 23)
443  status1ShortLived->Fill(9, weight);
444  if (std::abs(Id) == 24)
445  status1ShortLived->Fill(10, weight);
446  if (std::abs(Id) == 7 || std::abs(Id) == 8 || std::abs(Id) == 17 || (std::abs(Id) >= 25 && std::abs(Id) <= 99))
447  status1ShortLived->Fill(11, weight);
448  etotal += ptcl->momentum().e();
449  pxtotal += ptcl->momentum().px();
450  pytotal += ptcl->momentum().py();
451  pztotal += ptcl->momentum().pz();
452  }
453 
454  if (abs(Id) < 6 || abs(Id) == 22) {
455  ++partonNum;
456  partonpT->Fill(Log_p, weight);
457  }
458 
459  bool indentified = false;
460  for (unsigned int i = 0; i < particles.size(); i++) {
461  if (particles.at(i).Fill(ptcl, weight)) {
462  indentified = true;
463  break;
464  }
465  }
466  if (!indentified) {
467  ++otherPtclNum;
468  otherPtclMomentum->Fill(Log_p, weight);
469  }
470  } //event particles
471 
472  // set a default sqrt(s) and then check in the event
473  double ecms = 13000.;
474  if (myGenEvent->valid_beam_particles()) {
475  ecms = myGenEvent->beam_particles().first->momentum().e() + myGenEvent->beam_particles().second->momentum().e();
476  }
477  log10DeltaEcms->Fill(log10(fabs(etotal - ecms)), weight);
478  DeltaEcms->Fill(etotal - ecms, weight);
479  DeltaPx->Fill(pxtotal, weight);
480  DeltaPy->Fill(pytotal, weight);
481  DeltaPz->Fill(pztotal, weight);
482 
484  stablePtclNumber->Fill(log10(stablePtclNum + 0.1), weight);
485  stableChaNumber->Fill(log10(stableChaNum + 0.1), weight);
486  otherPtclNumber->Fill(log10(otherPtclNum + 0.1), weight);
487  unknownPDTNumber->Fill(log10(unknownPDTNum + 0.1), weight);
488  //
489  partonNumber->Fill(partonNum, weight);
490  for (unsigned int i = 0; i < particles.size(); i++) {
491  particles.at(i).FillCount(weight);
492  };
493 
494 } //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_
double logPdfMin
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
double logPdfBinsize
HepPDT::ParticleData ParticleData
MonitorElement * stablePtclp
MonitorElement * status1ShortLived
MonitorElement * stablePtclpT
MonitorElement * stableChaNumber
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
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 &   
)
overridevirtual

Setting the DQM top directories

Booking the ME's multiplicity

other

Implements DQMEDAnalyzer.

Definition at line 33 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, dqm::impl::MonitorElement::setBinLabel(), dqm::dqmstoreimpl::DQMStore::IBooker::setCurrentFolder(), stableChaNumber, stablePtclCharge, stablePtclEta, stablePtclNumber, stablePtclp, stablePtclPhi, stablePtclpT, status1ShortLived, unknownPDTNumber, vrtxRadius, and vrtxZ.

33  {
35  DQMHelper dqm(&i);
36  i.setCurrentFolder("Generator/Particles");
37 
38  // Number of analyzed events
39  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1.);
40 
43  // quarks
44  particles.push_back(ParticleMonitor("u", 1, i));
45  particles.push_back(ParticleMonitor("ubar", -1, i));
46  particles.push_back(ParticleMonitor("d", 2, i));
47  particles.push_back(ParticleMonitor("dbar", -2, i));
48  particles.push_back(ParticleMonitor("s", 3, i));
49  particles.push_back(ParticleMonitor("sbar", -3, i));
50  particles.push_back(ParticleMonitor("c", 4, i));
51  particles.push_back(ParticleMonitor("cbar", -4, i));
52  particles.push_back(ParticleMonitor("b", 5, i));
53  particles.push_back(ParticleMonitor("bbar", -5, i));
54  particles.push_back(ParticleMonitor("t", 6, i));
55  particles.push_back(ParticleMonitor("tbar", -6, i));
56 
57  //leptons
58  particles.push_back(ParticleMonitor("eminus", 11, i));
59  particles.push_back(ParticleMonitor("eplus", -11, i));
60  particles.push_back(ParticleMonitor("nue", 12, i));
61  particles.push_back(ParticleMonitor("nuebar", -12, i));
62  particles.push_back(ParticleMonitor("muminus", 13, i));
63  particles.push_back(ParticleMonitor("muplus", -13, i));
64  particles.push_back(ParticleMonitor("numu", 14, i));
65  particles.push_back(ParticleMonitor("numubar", -14, i));
66  particles.push_back(ParticleMonitor("tauminus", 15, i));
67  particles.push_back(ParticleMonitor("tauplus", -15, i));
68  particles.push_back(ParticleMonitor("nutau", 16, i));
69  particles.push_back(ParticleMonitor("nutaubar", -16, i));
70 
71  //bosons
72  particles.push_back(ParticleMonitor("Wplus", 24, i));
73  particles.push_back(ParticleMonitor("Wminus", -24, i));
74  particles.push_back(ParticleMonitor("Z", 23, i));
75  particles.push_back(ParticleMonitor("gamma", 22, i));
76  particles.push_back(ParticleMonitor("gluon", 21, i));
77 
78  //mesons
79  particles.push_back(ParticleMonitor("piplus", 211, i, true)); //log
80  particles.push_back(ParticleMonitor("piminus", -211, i, true)); //log
81  particles.push_back(ParticleMonitor("pizero", 111, i, true)); //log
82  particles.push_back(ParticleMonitor("Kplus", 321, i));
83  particles.push_back(ParticleMonitor("Kminus", -321, i));
84  particles.push_back(ParticleMonitor("Klzero", 130, i));
85  particles.push_back(ParticleMonitor("Kszero", 310, i));
86 
87  //baryons
88  particles.push_back(ParticleMonitor("p", 2212, i, true)); //log
89  particles.push_back(ParticleMonitor("pbar", -2212, i, true)); //log
90  particles.push_back(ParticleMonitor("n", 2112, i, true)); //log
91  particles.push_back(ParticleMonitor("nbar", -2112, i, true)); //log
92  particles.push_back(ParticleMonitor("lambda0", 3122, i));
93  particles.push_back(ParticleMonitor("lambda0bar", -3122, i));
94 
95  //D mesons
96  particles.push_back(ParticleMonitor("Dplus", 411, i));
97  particles.push_back(ParticleMonitor("Dminus", -411, i));
98  particles.push_back(ParticleMonitor("Dzero", 421, i));
99  particles.push_back(ParticleMonitor("Dzerobar", -421, i));
100 
101  //B mesons
102  particles.push_back(ParticleMonitor("Bplus", 521, i));
103  particles.push_back(ParticleMonitor("Bminus", -521, i));
104  particles.push_back(ParticleMonitor("Bzero", 511, i));
105  particles.push_back(ParticleMonitor("Bzerobar", -511, i));
106  particles.push_back(ParticleMonitor("Bszero", 531, i));
107  particles.push_back(ParticleMonitor("Bszerobar", -531, i));
108 
109  //
110  otherPtclNumber = dqm.book1dHisto(
111  "otherPtclNumber", "Log10(No. other ptcls)", 60, -1, 5, "log_{10}(No. other ptcls)", "Number of Events"); //Log
112  otherPtclMomentum = dqm.book1dHisto("otherPtclMomentum",
113  "Log10(p) other ptcls",
114  60,
115  -2,
116  4,
117  "log10(P^{other ptcls}) (log_{10}(GeV))",
118  "Number of Events");
119 
121  genPtclNumber = dqm.book1dHisto(
122  "genPtclNumber", "Log10(No. all particles)", 60, -1, 5, "log10(No. all particles)", "Number of Events"); //Log
123  genVrtxNumber = dqm.book1dHisto(
124  "genVrtxNumber", "Log10(No. all vertexs)", 60, -1, 5, "log10(No. all vertexs)", "Number of Events"); //Log
125  //
126  stablePtclNumber = dqm.book1dHisto("stablePtclNumber",
127  "Log10(No. stable particles)",
128  50,
129  0,
130  5,
131  "log10(No. stable particles)",
132  "Number of Events"); //Log
133  stablePtclPhi = dqm.book1dHisto(
134  "stablePtclPhi", "stable Ptcl Phi", 360, -180, 180, "#phi^{stable Ptcl} (rad)", "Number of Events");
135  stablePtclEta = dqm.book1dHisto(
136  "stablePtclEta", "stable Ptcl Eta (pseudo rapidity)", 220, -11, 11, "#eta^{stable Ptcl}", "Number of Events");
138  dqm.book1dHisto("stablePtclCharge", "stablePtclCharge", 5, -2, 2, "Charge^{stable ptcls}", "Number of Events");
139  stableChaNumber = dqm.book1dHisto("stableChaNumber",
140  "Log10(No. stable charged particles)",
141  50,
142  0,
143  5,
144  "log_{10}(No. stable charged particles)",
145  "Number of Events"); //Log
146  stablePtclp = dqm.book1dHisto("stablePtclp",
147  "Log10(p) stable ptcl p",
148  80,
149  -4,
150  4,
151  "log_{10}(P^{stable ptcl}) (log_{10}(GeV))",
152  "Number of Events"); //Log
153  stablePtclpT = dqm.book1dHisto("stablePtclpT",
154  "Log10(pT) stable ptcl pT",
155  80,
156  -4,
157  4,
158  "log_{10}(P_{t}^{stable ptcl}) (log_{10}(GeV))",
159  "Number of Events"); //Log
160  partonNumber =
161  dqm.book1dHisto("partonNumber", "number of partons", 100, 0, 100, "number of partons", "Number of Events");
162  partonpT =
163  dqm.book1dHisto("partonpT", "Log10(pT) parton pT", 80, -4, 4, "Log10(P_{t}^{parton})", "Number of Events"); //Log
164  outVrtxStablePtclNumber = dqm.book1dHisto("outVrtxStablePtclNumber",
165  "No. outgoing stable ptcls from vrtx",
166  10,
167  0,
168  10,
169  "No. outgoing stable ptcls from vrtx",
170  "Number of Events");
171  //
172  outVrtxPtclNumber = dqm.book1dHisto("outVrtxPtclNumber",
173  "No. outgoing ptcls from vrtx",
174  30,
175  0,
176  30,
177  "No. outgoing ptcls from vrtx",
178  "Number of Events");
179  vrtxZ = dqm.book1dHisto("VrtxZ", "VrtxZ", 50, -250, 250, "Z_{Vtx}", "Number of Events");
180  vrtxRadius = dqm.book1dHisto("vrtxRadius", "vrtxRadius", 50, 0, 50, "R_{vtx}", "Number of Events");
181  //
182  unknownPDTNumber = dqm.book1dHisto("unknownPDTNumber",
183  "Log10(No. unknown ptcls PDT)",
184  60,
185  -1,
186  5,
187  "log_{10}(No. unknown ptcls PDT)",
188  "Number of Events"); //Log
189  genPtclStatus = dqm.book1dHisto("genPtclStatus", "Status of genParticle", 200, 0, 200., "", "Number of Events");
190  //
191  Bjorken_x = dqm.book1dHisto("Bjorken_x", "Bjorken_x", 1000, 0.0, 1.0, "Bjorken_{x}", "Number of Events");
192  pdf_u = dqm.book1dHisto(
193  "pdf_u", "Log10(PDF(u,x,Q))", logPdfNbin, logPdfMin, logPdfMax, "log_{10}(x*f(x))", "Number of Events"); //Log
194  pdf_ubar = dqm.book1dHisto("pdf_ubar",
195  "Log10(PDF(ubar,x,Q))",
196  logPdfNbin,
197  logPdfMin,
198  logPdfMax,
199  "log_{10}(x*f(x))",
200  "Number of Events"); //Log
201  pdf_d = dqm.book1dHisto(
202  "pdf_d", "Log10(PDF(d,x,Q))", logPdfNbin, logPdfMin, logPdfMax, "log_{10}(x*f(x))", "Number of Events"); //Log
203  pdf_dbar = dqm.book1dHisto("pdf_dbar",
204  "Log10(PDF(dbar,x,Q))",
205  logPdfNbin,
206  logPdfMin,
207  logPdfMax,
208  "log_{10}(x*f(x))",
209  "Number of Events"); //Log
210  pdf_ssbar = dqm.book1dHisto("pdf_ssbar",
211  "Log10(PDF(ssbar,x,Q))",
212  logPdfNbin,
213  logPdfMin,
214  logPdfMax,
215  "log_{10}(x*f(x))",
216  "Number of Events"); //Log
217  pdf_ccbar = dqm.book1dHisto("pdf_ccbar",
218  "Log10(PDF(ccbar,x,Q))",
219  logPdfNbin,
220  logPdfMin,
221  logPdfMax,
222  "log_{10}(x*f(x))",
223  "Number of Events"); //Log
224  pdf_bbbar = dqm.book1dHisto("pdf_bbbar",
225  "Log10(PDF(bbbar,x,Q))",
226  logPdfNbin,
227  logPdfMin,
228  logPdfMax,
229  "log_{10}(x*f(x))",
230  "Number of Events"); //Log
231  pdf_g = dqm.book1dHisto(
232  "pdf_g", "Log10(PDF(g,x,Q))", logPdfNbin, logPdfMin, logPdfMax, "log_{10}(x*f(x))", "Number of Events"); //Log
233  scalePDF = dqm.book1dHisto("scalePDF",
234  "Log10(Q-scale(GeV))",
236  logQScaleMin,
237  logQScaleMax,
238  "log_{10}(Q-scale(GeV))",
239  "Number of Events");
240  parton1Id = dqm.book1dHisto("parton1Id", "ID of parton 1", 45, -14.5, 30.5, "ID", "Number of Events");
241  parton2Id = dqm.book1dHisto("parton2Id", "ID of parton 2", 45, -14.5, 30.5, "ID", "Number of Events");
242  //
243  status1ShortLived = dqm.book1dHisto("status1ShortLived", "Status 1 short lived", 11, 0, 11, "", "Number of Events");
244  status1ShortLived->setBinLabel(1, "d/dbar");
245  status1ShortLived->setBinLabel(2, "u/ubar");
246  status1ShortLived->setBinLabel(3, "s/sbar");
247  status1ShortLived->setBinLabel(4, "c/cbar");
248  status1ShortLived->setBinLabel(5, "b/bbar");
249  status1ShortLived->setBinLabel(6, "t/tbar");
251  status1ShortLived->setBinLabel(8, "tau-/tau+");
252  status1ShortLived->setBinLabel(9, "Z0");
253  status1ShortLived->setBinLabel(10, "W-/W+");
254  status1ShortLived->setBinLabel(11, "PDG = 7,8,17,25-99");
255 
256  log10DeltaEcms = dqm.book1dHisto("DeltaEcms1log10",
257  "log_{10} of deviation from nominal Ecms",
258  200,
259  -1.,
260  5.,
261  "log_{10}(#DeltaE) (log_{10}(GeV))",
262  "Number of Events");
263  DeltaEcms =
264  dqm.book1dHisto("DeltaEcms1", "deviation from nominal Ecms", 200, -1., 1., "#DeltaE (GeV)", "Number of Events");
265  DeltaPx =
266  dqm.book1dHisto("DeltaPx1", "deviation from nominal Px", 200, -1., 1., "#DeltaP_{x} (GeV)", "Number of Events");
267  DeltaPy =
268  dqm.book1dHisto("DeltaPy1", "deviation from nominal Py", 200, -1., 1., "#DeltaP_{y} (GeV)", "Number of Events");
269  DeltaPz =
270  dqm.book1dHisto("DeltaPz1", "deviation from nominal Pz", 200, -1., 1., "#DeltaP_{z} (GeV)", "Number of Events");
271  return;
272 }
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
double logPdfMin
MonitorElement * vrtxRadius
int logQScaleNbin
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)
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 21 of file BasicHepMCValidation.cc.

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

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

Member Data Documentation

MonitorElement* BasicHepMCValidation::Bjorken_x
private

Definition at line 220 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaEcms
private

Definition at line 236 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaPx
private

Definition at line 237 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaPy
private

Definition at line 238 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::DeltaPz
private

Definition at line 239 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

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

PDT table.

Definition at line 48 of file BasicHepMCValidation.h.

Referenced by analyze(), and dqmBeginRun().

MonitorElement* BasicHepMCValidation::genPtclNumber
private

Definition at line 200 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::genPtclStatus
private

Definition at line 204 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::genVrtxNumber
private

Definition at line 201 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

edm::InputTag BasicHepMCValidation::hepmcCollection_
private

Definition at line 45 of file BasicHepMCValidation.h.

Referenced by BasicHepMCValidation().

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

Definition at line 241 of file BasicHepMCValidation.h.

Referenced by analyze(), and BasicHepMCValidation().

MonitorElement* BasicHepMCValidation::log10DeltaEcms
private

Definition at line 235 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::nEvt
private

Definition at line 194 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::otherPtclMomentum
private

Definition at line 199 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::otherPtclNumber
private

other ME's

Definition at line 198 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::outVrtxPtclNumber
private

Definition at line 203 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::outVrtxStablePtclNumber
private

Definition at line 215 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

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

Definition at line 195 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::parton1Id
private

Definition at line 230 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::parton2Id
private

Definition at line 231 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::partonNumber
private

Definition at line 213 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::partonpT
private

Definition at line 214 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_bbbar
private

Definition at line 227 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_ccbar
private

Definition at line 226 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_d
private

Definition at line 223 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_dbar
private

Definition at line 224 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_g
private

Definition at line 228 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_ssbar
private

Definition at line 225 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_u
private

Definition at line 221 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::pdf_ubar
private

Definition at line 222 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::scalePDF
private

Definition at line 229 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stableChaNumber
private

Definition at line 207 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclCharge
private

Definition at line 210 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclEta
private

Definition at line 209 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclNumber
private

Definition at line 206 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclp
private

Definition at line 211 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclPhi
private

Definition at line 208 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::stablePtclpT
private

Definition at line 212 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::status1ShortLived
private

Definition at line 233 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::unknownPDTNumber
private

Definition at line 202 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::vrtxRadius
private

Definition at line 218 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

MonitorElement* BasicHepMCValidation::vrtxZ
private

Definition at line 217 of file BasicHepMCValidation.h.

Referenced by analyze(), and bookHistograms().

WeightManager BasicHepMCValidation::wmanager_
private

Definition at line 44 of file BasicHepMCValidation.h.

Referenced by analyze().