CMS 3D CMS Logo

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

#include <BasicGenParticleValidation.h>

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

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
 BasicGenParticleValidation (const edm::ParameterSet &)
 
void bookHistograms (DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
 
bool matchParticles (const HepMC::GenParticle *&, const reco::GenParticle *&)
 
 ~BasicGenParticleValidation () 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
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 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

MonitorElementgenJetCentral
 
edm::InputTag genjetCollection_
 
edm::EDGetTokenT< reco::GenJetCollectiongenjetCollectionToken_
 
MonitorElementgenJetDeltaEtaMin
 
MonitorElementgenJetEnergy
 
MonitorElementgenJetEta
 
MonitorElementgenJetMult
 
MonitorElementgenJetPhi
 
MonitorElementgenJetPt
 
MonitorElementgenJetPto1
 
MonitorElementgenJetPto10
 
MonitorElementgenJetPto100
 
MonitorElementgenJetTotPt
 
MonitorElementgenMatched
 
edm::InputTag genparticleCollection_
 
edm::EDGetTokenT< reco::GenParticleCollectiongenparticleCollectionToken_
 
MonitorElementgenPMultiplicity
 
edm::InputTag hepmcCollection_
 
edm::EDGetTokenT< edm::HepMCProducthepmcCollectionToken_
 
MonitorElementmatchedResolution
 
double matchPr_
 
MonitorElementmultipleMatching
 
MonitorElementnEvt
 
bool signalParticlesOnly_
 
unsigned int verbosity_
 
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 36 of file BasicGenParticleValidation.h.

Constructor & Destructor Documentation

◆ BasicGenParticleValidation()

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

Definition at line 15 of file BasicGenParticleValidation.cc.

References genjetCollection_, genjetCollectionToken_, genparticleCollection_, genparticleCollectionToken_, hepmcCollection_, and hepmcCollectionToken_.

16  : wmanager_(iPSet, consumesCollector()),
17  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
18  genparticleCollection_(iPSet.getParameter<edm::InputTag>("genparticleCollection")),
19  genjetCollection_(iPSet.getParameter<edm::InputTag>("genjetsCollection")),
20  matchPr_(iPSet.getParameter<double>("matchingPrecision")),
21  verbosity_(iPSet.getUntrackedParameter<unsigned int>("verbosity", 0)),
22  signalParticlesOnly_(iPSet.getParameter<bool>("signalParticlesOnly")) {
23  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
24  genparticleCollectionToken_ = consumes<reco::GenParticleCollection>(genparticleCollection_);
25  genjetCollectionToken_ = consumes<reco::GenJetCollection>(genjetCollection_);
26 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< reco::GenJetCollection > genjetCollectionToken_
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_

◆ ~BasicGenParticleValidation()

BasicGenParticleValidation::~BasicGenParticleValidation ( )
override

Definition at line 28 of file BasicGenParticleValidation.cc.

28 {}

Member Function Documentation

◆ analyze()

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

Gathering the HepMCProduct information

Reimplemented from DQMEDAnalyzer.

Definition at line 98 of file BasicGenParticleValidation.cc.

References gather_cfg::cout, spr::deltaEta, PVValHelper::eta, dqm::impl::MonitorElement::Fill(), alignBH_cfg::fixed, genJetCentral, genjetCollectionToken_, genJetDeltaEtaMin, genJetEnergy, genJetEta, genJetMult, genJetPhi, genJetPt, genJetPto1, genJetPto10, genJetPto100, l1tCaloJetHTTProducer_cfi::genJets, genJetTotPt, genMatched, genparticleCollectionToken_, AJJGenJetFilter_cfi::genParticles, genPMultiplicity, edm::HepMCProduct::GetEvent(), hepmcCollectionToken_, mps_fire::i, iEvent, dqmiolumiharvest::j, reco::btau::jetEta, matchedResolution, matchParticles(), SiStripPI::min, multipleMatching, nEvt, l1ctJetFileWriter_cfi::nJets, ecalTrigSettings_cff::particles, DiDispStaMuonMonitor_cfi::pt, rho, signalParticlesOnly_, verbosity_, WeightManager::weight(), mps_merge::weight, and wmanager_.

98  {
99  unsigned int initSize = 1000;
100 
103  iEvent.getByToken(hepmcCollectionToken_, evt);
104 
105  //Get HepMC EVENT
106  HepMC::GenEvent* myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
107 
108  double weight = wmanager_.weight(iEvent);
109 
110  nEvt->Fill(0.5, weight);
111 
112  std::vector<const HepMC::GenParticle*> hepmcGPCollection;
113  std::vector<int> barcodeList;
114  hepmcGPCollection.reserve(initSize);
115  barcodeList.reserve(initSize);
116 
117  //Looping through HepMC::GenParticle collection to search for status 1 particles
118  for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
119  iter != myGenEvent->particles_end();
120  ++iter) {
121  if ((*iter)->status() == 1) {
122  hepmcGPCollection.push_back(*iter);
123  barcodeList.push_back((*iter)->barcode());
124  if (verbosity_ > 0) {
125  std::cout << "HepMC " << std::setw(14) << std::fixed << (*iter)->pdg_id() << std::setw(14) << std::fixed
126  << (*iter)->momentum().px() << std::setw(14) << std::fixed << (*iter)->momentum().py()
127  << std::setw(14) << std::fixed << (*iter)->momentum().pz() << std::endl;
128  }
129  }
130  }
131 
132  // Gather information on the reco::GenParticle collection
135 
136  std::vector<const reco::GenParticle*> particles;
137  particles.reserve(initSize);
138  for (reco::GenParticleCollection::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
139  if ((*iter).status() == 1) {
140  particles.push_back(&*iter);
141  if (verbosity_ > 0) {
142  std::cout << "reco " << std::setw(14) << std::fixed << (*iter).pdgId() << std::setw(14) << std::fixed
143  << (*iter).px() << std::setw(14) << std::fixed << (*iter).py() << std::setw(14) << std::fixed
144  << (*iter).pz() << std::endl;
145  }
146  }
147  }
148 
149  unsigned int nReco = particles.size();
150  unsigned int nHepMC = hepmcGPCollection.size();
151 
152  if (signalParticlesOnly_) {
153  for (unsigned int i = 0; i < particles.size(); ++i) {
154  if (particles[i]->collisionId() != 0)
155  nReco -= 1;
156  }
157  }
158 
159  genPMultiplicity->Fill(std::log10(nReco), weight);
160 
161  // Define vector containing index of hepmc corresponding to the reco::GenParticle
162  std::vector<int> hepmcMatchIndex;
163  hepmcMatchIndex.reserve(initSize);
164 
165  // Matching procedure
166 
167  // Check array size consistency
168 
169  if (nReco != nHepMC) {
170  edm::LogWarning("CollectionSizeInconsistency")
171  << "Collection size inconsistency: HepMC::GenParticle = " << nHepMC << " reco::GenParticle = " << nReco;
172  }
173 
174  // Match each HepMC with a reco
175 
176  for (unsigned int i = 0; i < nReco; ++i) {
177  for (unsigned int j = 0; j < nHepMC; ++j) {
178  if (matchParticles(hepmcGPCollection[j], particles[i])) {
179  hepmcMatchIndex.push_back((int)j);
180  if (hepmcGPCollection[j]->momentum().rho() != 0.) {
181  double reso = 1. - particles[i]->p() / hepmcGPCollection[j]->momentum().rho();
182  if (verbosity_ > 0) {
183  std::cout << "Matching momentum: reco = " << particles[i]->p()
184  << " HepMC = " << hepmcGPCollection[j]->momentum().rho() << " resoultion = " << reso << std::endl;
185  }
186  matchedResolution->Fill(std::log10(std::fabs(reso)), weight);
187  }
188  continue;
189  }
190  }
191  }
192 
193  // Check unicity of matching
194 
195  unsigned int nMatched = hepmcMatchIndex.size();
196 
197  if (nMatched != nReco) {
198  edm::LogWarning("IncorrectMatching") << "Incorrect number of matched indexes: GenParticle = " << nReco
199  << " matched indexes = " << nMatched;
200  }
201  genMatched->Fill(int(nReco - nMatched), weight);
202 
203  unsigned int nWrMatch = 0;
204 
205  for (unsigned int i = 0; i < nMatched; ++i) {
206  for (unsigned int j = i + 1; j < nMatched; ++j) {
207  if (hepmcMatchIndex[i] == hepmcMatchIndex[j]) {
208  int theIndex = hepmcMatchIndex[i];
209  edm::LogWarning("DuplicatedMatching")
210  << "Multiple matching occurencies for GenParticle barcode = " << barcodeList[theIndex];
211  nWrMatch++;
212  }
213  }
214  }
215  multipleMatching->Fill(int(nWrMatch), weight);
216 
217  // Gather information in the GenJet collection
220 
221  int nJets = 0;
222  int nJetso1 = 0;
223  int nJetso10 = 0;
224  int nJetso100 = 0;
225  int nJetsCentral = 0;
226  double totPt = 0.;
227 
228  std::vector<double> jetEta;
229  jetEta.reserve(initSize);
230 
231  for (reco::GenJetCollection::const_iterator iter = genJets->begin(); iter != genJets->end(); ++iter) {
232  nJets++;
233  double pt = (*iter).pt();
234  totPt += pt;
235  if (pt > 1.)
236  nJetso1++;
237  if (pt > 10.)
238  nJetso10++;
239  if (pt > 100.)
240  nJetso100++;
241  double eta = (*iter).eta();
242  if (std::fabs(eta) < 2.5)
243  nJetsCentral++;
244  jetEta.push_back(eta);
245 
246  genJetEnergy->Fill(std::log10((*iter).energy()), weight);
247  genJetPt->Fill(std::log10(pt), weight);
249  genJetPhi->Fill((*iter).phi() / CLHEP::degree, weight);
250  }
251 
253  genJetPto1->Fill(nJetso1, weight);
254  genJetPto10->Fill(nJetso10, weight);
255  genJetPto100->Fill(nJetso100, weight);
256  genJetCentral->Fill(nJetsCentral, weight);
257 
258  genJetTotPt->Fill(std::log10(totPt), weight);
259 
260  double deltaEta = 999.;
261  if (jetEta.size() > 1) {
262  for (unsigned int i = 0; i < jetEta.size(); i++) {
263  for (unsigned int j = i + 1; j < jetEta.size(); j++) {
264  deltaEta = std::min(deltaEta, std::fabs(jetEta[i] - jetEta[j]));
265  }
266  }
267  }
268 
270 
271  delete myGenEvent;
272 } //analyze
Definition: weight.py:1
bool matchParticles(const HepMC::GenParticle *&, const reco::GenParticle *&)
static const double deltaEta
Definition: CaloConstants.h:8
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::GenJetCollection > genjetCollectionToken_
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
Log< level::Warning, false > LogWarning
double weight(const edm::Event &)
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_

◆ bookHistograms()

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

Setting the DQM top directories

Booking the ME's

multiplicity

Implements DQMEDAnalyzer.

Definition at line 30 of file BasicGenParticleValidation.cc.

References genJetCentral, genJetDeltaEtaMin, genJetEnergy, genJetEta, genJetMult, genJetPhi, genJetPt, genJetPto1, genJetPto10, genJetPto100, genJetTotPt, genMatched, genPMultiplicity, mps_fire::i, matchedResolution, multipleMatching, and nEvt.

30  {
32  DQMHelper dqm(&i);
33  i.setCurrentFolder("Generator/GenParticles");
35 
36  // Number of analyzed events
37  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "", "Number of Events");
38 
40  genPMultiplicity = dqm.book1dHisto("genPMultiplicty",
41  "Log(No. all GenParticles)",
42  50,
43  -1,
44  5,
45  "log_{10}(N_{All GenParticles})",
46  "Number of Events"); //Log
47  //difference in HepMC and reco multiplicity
48  genMatched = dqm.book1dHisto(
49  "genMatched", "Difference reco - matched", 50, -25, 25, "N_{All GenParticles}-N_{Matched}", "Number of Events");
50  //multiple matching
51  multipleMatching = dqm.book1dHisto("multipleMatching",
52  "multiple reco HepMC matching",
53  50,
54  0,
55  50,
56  "N_{multiple reco HepMC matching}",
57  "Number of Events");
58  //momentum difference of matched particles
59  matchedResolution = dqm.book1dHisto("matchedResolution",
60  "log10(momentum difference of matched particles)",
61  70,
62  -10.,
63  -3.,
64  "log_{10}(#DeltaP_{matched Particles})",
65  "Number of Events");
66 
67  // GenJet general distributions
68  genJetMult = dqm.book1dHisto("genJetMult", "GenJet multiplicity", 50, 0, 50, "N_{gen-jets}", "Number of Events");
69  genJetEnergy = dqm.book1dHisto(
70  "genJetEnergy", "Log10(GenJet energy)", 60, -1, 5, "log_{10}(E^{gen-jets}) (log_{10}(GeV))", "Number of Events");
71  genJetPt = dqm.book1dHisto(
72  "genJetPt", "Log10(GenJet pt)", 60, -1, 5, "log_{10}(P_{t}^{gen-jets}) (log_{10}(GeV))", "Number of Events");
73  genJetEta = dqm.book1dHisto("genJetEta", "GenJet eta", 220, -11, 11, "#eta^{gen-jets}", "Number of Events");
74  genJetPhi = dqm.book1dHisto("genJetPhi", "GenJet phi", 360, -180, 180, "#phi^{gen-jets} (rad)", "Number of Events");
75  genJetDeltaEtaMin = dqm.book1dHisto(
76  "genJetDeltaEtaMin", "GenJet minimum rapidity gap", 30, 0, 30, "#delta#eta_{min}^{gen-jets}", "Number of Events");
77 
78  genJetPto1 = dqm.book1dHisto(
79  "genJetPto1", "GenJet multiplicity above 1 GeV", 50, 0, 50, "N_{gen-jets P_{t}>1GeV}", "Number of Events");
80  genJetPto10 = dqm.book1dHisto(
81  "genJetPto10", "GenJet multiplicity above 10 GeV", 50, 0, 50, "N_{gen-jets P_{t}>10GeV}", "Number of Events");
82  genJetPto100 = dqm.book1dHisto(
83  "genJetPto100", "GenJet multiplicity above 100 GeV", 50, 0, 50, "N_{gen-jets P_{t}>100GeV}", "Number of Events");
84  genJetCentral = dqm.book1dHisto(
85  "genJetCentral", "GenJet multiplicity |eta|.lt.2.5", 50, 0, 50, "N_{gen-jets |#eta|#leq2.5}", "Number of Events");
86 
87  genJetTotPt = dqm.book1dHisto("genJetTotPt",
88  "Log10(GenJet total pt)",
89  100,
90  -5,
91  5,
92  "log_{10}(#SigmaP_{t}^{gen-jets}) (log_{10}(GeV))",
93  "Number of Events");
94 
95  return;
96 }
Definition: DQMStore.h:18

◆ matchParticles()

bool BasicGenParticleValidation::matchParticles ( const HepMC::GenParticle *&  hepmcP,
const reco::GenParticle *&  recoP 
)

Definition at line 274 of file BasicGenParticleValidation.cc.

References matchPr_, reco::LeafCandidate::pdgId(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), and reco::LeafCandidate::pz().

Referenced by analyze().

274  {
275  bool state = false;
276 
277  if (hepmcP->pdg_id() != recoP->pdgId())
278  return state;
279  if (std::fabs(hepmcP->momentum().px() - recoP->px()) < std::fabs(matchPr_ * hepmcP->momentum().px()) &&
280  std::fabs(hepmcP->momentum().py() - recoP->py()) < std::fabs(matchPr_ * hepmcP->momentum().py()) &&
281  std::fabs(hepmcP->momentum().pz() - recoP->pz()) < std::fabs(matchPr_ * hepmcP->momentum().pz())) {
282  state = true;
283  }
284 
285  return state;
286 }
double pz() const final
z coordinate of momentum vector
int pdgId() const final
PDG identifier.
double px() const final
x coordinate of momentum vector
double py() const final
y coordinate of momentum vector

Member Data Documentation

◆ genJetCentral

MonitorElement* BasicGenParticleValidation::genJetCentral
private

Definition at line 77 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genjetCollection_

edm::InputTag BasicGenParticleValidation::genjetCollection_
private

Definition at line 50 of file BasicGenParticleValidation.h.

Referenced by BasicGenParticleValidation().

◆ genjetCollectionToken_

edm::EDGetTokenT<reco::GenJetCollection> BasicGenParticleValidation::genjetCollectionToken_
private

Definition at line 83 of file BasicGenParticleValidation.h.

Referenced by analyze(), and BasicGenParticleValidation().

◆ genJetDeltaEtaMin

MonitorElement* BasicGenParticleValidation::genJetDeltaEtaMin
private

Definition at line 72 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetEnergy

MonitorElement* BasicGenParticleValidation::genJetEnergy
private

Definition at line 68 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetEta

MonitorElement* BasicGenParticleValidation::genJetEta
private

Definition at line 70 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetMult

MonitorElement* BasicGenParticleValidation::genJetMult
private

Definition at line 67 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetPhi

MonitorElement* BasicGenParticleValidation::genJetPhi
private

Definition at line 71 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetPt

MonitorElement* BasicGenParticleValidation::genJetPt
private

Definition at line 69 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetPto1

MonitorElement* BasicGenParticleValidation::genJetPto1
private

Definition at line 74 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetPto10

MonitorElement* BasicGenParticleValidation::genJetPto10
private

Definition at line 75 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetPto100

MonitorElement* BasicGenParticleValidation::genJetPto100
private

Definition at line 76 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genJetTotPt

MonitorElement* BasicGenParticleValidation::genJetTotPt
private

Definition at line 79 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genMatched

MonitorElement* BasicGenParticleValidation::genMatched
private

Definition at line 61 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ genparticleCollection_

edm::InputTag BasicGenParticleValidation::genparticleCollection_
private

Definition at line 49 of file BasicGenParticleValidation.h.

Referenced by BasicGenParticleValidation().

◆ genparticleCollectionToken_

edm::EDGetTokenT<reco::GenParticleCollection> BasicGenParticleValidation::genparticleCollectionToken_
private

Definition at line 82 of file BasicGenParticleValidation.h.

Referenced by analyze(), and BasicGenParticleValidation().

◆ genPMultiplicity

MonitorElement* BasicGenParticleValidation::genPMultiplicity
private

Definition at line 60 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ hepmcCollection_

edm::InputTag BasicGenParticleValidation::hepmcCollection_
private

Definition at line 48 of file BasicGenParticleValidation.h.

Referenced by BasicGenParticleValidation().

◆ hepmcCollectionToken_

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

Definition at line 81 of file BasicGenParticleValidation.h.

Referenced by analyze(), and BasicGenParticleValidation().

◆ matchedResolution

MonitorElement* BasicGenParticleValidation::matchedResolution
private

Definition at line 63 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ matchPr_

double BasicGenParticleValidation::matchPr_
private

Definition at line 51 of file BasicGenParticleValidation.h.

Referenced by matchParticles().

◆ multipleMatching

MonitorElement* BasicGenParticleValidation::multipleMatching
private

Definition at line 62 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ nEvt

MonitorElement* BasicGenParticleValidation::nEvt
private

Definition at line 56 of file BasicGenParticleValidation.h.

Referenced by analyze(), and bookHistograms().

◆ signalParticlesOnly_

bool BasicGenParticleValidation::signalParticlesOnly_
private

Definition at line 54 of file BasicGenParticleValidation.h.

Referenced by analyze().

◆ verbosity_

unsigned int BasicGenParticleValidation::verbosity_
private

Definition at line 53 of file BasicGenParticleValidation.h.

Referenced by analyze().

◆ wmanager_

WeightManager BasicGenParticleValidation::wmanager_
private

Definition at line 47 of file BasicGenParticleValidation.h.

Referenced by analyze().