CMS 3D CMS Logo

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

#include <GenWeightValidation.h>

Inheritance diagram for GenWeightValidation:
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
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 
 GenWeightValidation (const edm::ParameterSet &)
 
 ~GenWeightValidation () override=default
 
- 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 Member Functions

void bookTemplates (DQMHelper &aDqmHelper, std::vector< MonitorElement *> &tmps, const std::string &name, const std::string &title, int nbin, float low, float high, const std::string &xtitle, const std::string &ytitle)
 
void fillTemplates (std::vector< MonitorElement *> &tmps, float obs)
 

Private Attributes

const edm::EDGetTokenT< reco::GenJetCollectiongenJetToken_
 
const edm::EDGetTokenT< reco::GenParticleCollectiongenParticleToken_
 
const int idxFSRdown_
 
const int idxFSRup_
 
const int idxGenEvtInfo_
 
const int idxISRdown_
 
const int idxISRup_
 
int idxMax_
 
const double jetEtaCut_
 
std::vector< MonitorElement * > jetMultTemp_
 
const double jetPtCut_
 
const int jetPtNbin_
 
const double jetPtRange_
 
std::vector< MonitorElement * > leadJetEtaTemp_
 
std::vector< MonitorElement * > leadJetPtTemp_
 
std::vector< MonitorElement * > leadLepEtaTemp_
 
const double leadLepPtCut_
 
const int leadLepPtNbin_
 
const double leadLepPtRange_
 
std::vector< MonitorElement * > leadLepPtTemp_
 
const double lepEtaCut_
 
MonitorElementnEvt_
 
const int nJetsNbin_
 
MonitorElementnlogWgt_
 
const int rapidityNbin_
 
const double rapidityRange_
 
double weight_
 
std::vector< std::vector< double > > weights_
 
MonitorElementwgtVal_
 
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 26 of file GenWeightValidation.h.

Constructor & Destructor Documentation

◆ GenWeightValidation()

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

Definition at line 8 of file GenWeightValidation.cc.

References idxFSRdown_, idxFSRup_, idxISRdown_, idxISRup_, idxMax_, and jetUpdater_cfi::sort.

9  : wmanager_(iPSet, consumesCollector()),
10  genParticleToken_(consumes<reco::GenParticleCollection>(iPSet.getParameter<edm::InputTag>("genParticles"))),
11  genJetToken_(consumes<reco::GenJetCollection>(iPSet.getParameter<edm::InputTag>("genJets"))),
12  idxGenEvtInfo_(iPSet.getParameter<int>("whichGenEventInfo")),
13  idxFSRup_(iPSet.getParameter<int>("idxFSRup")),
14  idxFSRdown_(iPSet.getParameter<int>("idxFSRdown")),
15  idxISRup_(iPSet.getParameter<int>("idxISRup")),
16  idxISRdown_(iPSet.getParameter<int>("idxISRdown")),
17  leadLepPtNbin_(iPSet.getParameter<int>("leadLepPtNbin")),
18  rapidityNbin_(iPSet.getParameter<int>("rapidityNbin")),
19  leadLepPtRange_(iPSet.getParameter<double>("leadLepPtRange")),
20  leadLepPtCut_(iPSet.getParameter<double>("leadLepPtCut")),
21  lepEtaCut_(iPSet.getParameter<double>("lepEtaCut")),
22  rapidityRange_(iPSet.getParameter<double>("rapidityRange")),
23  nJetsNbin_(iPSet.getParameter<int>("nJetsNbin")),
24  jetPtNbin_(iPSet.getParameter<int>("jetPtNbin")),
25  jetPtCut_(iPSet.getParameter<double>("jetPtCut")),
26  jetEtaCut_(iPSet.getParameter<double>("jetEtaCut")),
27  jetPtRange_(iPSet.getParameter<double>("jetPtRange")) {
28  std::vector<int> idxs = {idxFSRup_, idxFSRdown_, idxISRup_, idxISRdown_};
29  std::sort(idxs.begin(), idxs.end(), std::greater<int>());
30  idxMax_ = idxs.at(0);
31 }
const double leadLepPtRange_
const edm::EDGetTokenT< reco::GenParticleCollection > genParticleToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::EDGetTokenT< reco::GenJetCollection > genJetToken_

◆ ~GenWeightValidation()

GenWeightValidation::~GenWeightValidation ( )
overridedefault

Member Function Documentation

◆ analyze()

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

Reimplemented from DQMEDAnalyzer.

Definition at line 122 of file GenWeightValidation.cc.

References funct::abs(), dqm::impl::MonitorElement::Fill(), fillTemplates(), genJetToken_, genParticleToken_, heavyIonCSV_trainingSettings::idx, idxGenEvtInfo_, idxMax_, iEvent, jetEtaCut_, jetMultTemp_, jetPtCut_, leadJetEtaTemp_, leadJetPtTemp_, leadLepEtaTemp_, leadLepPtCut_, leadLepPtTemp_, lepEtaCut_, HLT_2022v11_cff::leptons, nEvt_, nlogWgt_, jetUpdater_cfi::sort, weight_, weights_, WeightManager::weightsCollection(), wgtVal_, and wmanager_.

122  {
124 
125  unsigned weightsSize = weights_.at(idxGenEvtInfo_).size();
126  if (weightsSize < 2)
127  return; // no baseline weight in GenEventInfo
128 
130  nEvt_->Fill(0.5, weight_);
131  nlogWgt_->Fill(std::log10(weightsSize), weight_);
132 
133  for (unsigned idx = 0; idx < weightsSize; idx++)
135 
136  if ((int)weightsSize <= idxMax_)
137  return; // no PS weights in GenEventInfo
138 
140  iEvent.getByToken(genParticleToken_, ptcls);
142  iEvent.getByToken(genJetToken_, genjets);
143 
144  std::vector<reco::GenParticleRef> leptons;
145 
146  for (unsigned iptc = 0; iptc < ptcls->size(); iptc++) {
147  reco::GenParticleRef ptc(ptcls, iptc);
148  if (ptc->status() == 1 && (std::abs(ptc->pdgId()) == 11 || std::abs(ptc->pdgId()) == 13)) {
149  if (ptc->pt() > leadLepPtCut_ && std::abs(ptc->eta()) < lepEtaCut_)
150  leptons.push_back(ptc);
151  }
152  }
153 
154  std::sort(leptons.begin(), leptons.end(), HepMCValidationHelper::sortByPtRef<reco::GenParticleRef>);
155 
156  if (!leptons.empty()) {
157  reco::GenParticleRef leadLep = leptons.at(0);
158  fillTemplates(leadLepPtTemp_, leadLep->pt());
159  fillTemplates(leadLepEtaTemp_, leadLep->eta());
160  }
161 
162  std::vector<reco::GenJetRef> genjetVec;
163 
164  for (unsigned igj = 0; igj < genjets->size(); igj++) {
165  reco::GenJetRef genjet(genjets, igj);
166 
167  if (genjet->pt() > jetPtCut_ && std::abs(genjet->eta()) < jetEtaCut_)
168  genjetVec.push_back(genjet);
169  }
170 
171  fillTemplates(jetMultTemp_, (float)genjetVec.size());
172 
173  if (!genjetVec.empty()) {
174  std::sort(genjetVec.begin(), genjetVec.end(), HepMCValidationHelper::sortByPtRef<reco::GenJetRef>);
175 
176  auto leadJet = genjetVec.at(0);
177  fillTemplates(leadJetPtTemp_, leadJet->pt());
178  fillTemplates(leadJetEtaTemp_, leadJet->eta());
179  }
180 } //analyze
std::vector< MonitorElement * > leadLepEtaTemp_
const edm::EDGetTokenT< reco::GenParticleCollection > genParticleToken_
MonitorElement * nEvt_
void fillTemplates(std::vector< MonitorElement *> &tmps, float obs)
std::vector< MonitorElement * > jetMultTemp_
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
MonitorElement * nlogWgt_
std::vector< MonitorElement * > leadJetEtaTemp_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< MonitorElement * > leadJetPtTemp_
std::vector< MonitorElement * > leadLepPtTemp_
MonitorElement * wgtVal_
const edm::EDGetTokenT< reco::GenJetCollection > genJetToken_
std::vector< std::vector< double > > weightsCollection(const edm::Event &)
std::vector< std::vector< double > > weights_

◆ bookHistograms()

void GenWeightValidation::bookHistograms ( DQMStore::IBooker iBook,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Setting the DQM top directories

Implements DQMEDAnalyzer.

Definition at line 33 of file GenWeightValidation.cc.

References DQMHelper::book1dHisto(), bookTemplates(), ALCARECODTCalibSynchCosmicsDQM_cff::folderName, jetMultTemp_, jetPtNbin_, jetPtRange_, leadJetEtaTemp_, leadJetPtTemp_, leadLepEtaTemp_, leadLepPtNbin_, leadLepPtRange_, leadLepPtTemp_, nEvt_, nJetsNbin_, nlogWgt_, rapidityNbin_, rapidityRange_, dqm::implementation::NavigatorBase::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, and wgtVal_.

33  {
35  std::string folderName = "Generator/GenWeight";
36  DQMHelper aDqmHelper(&iBook);
38 
39  // Number of analyzed events
40  nEvt_ = aDqmHelper.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "bin", "Number of Events");
41  nlogWgt_ = aDqmHelper.book1dHisto("nlogWgt", "Log10(n weights)", 100, 0., 3., "log_{10}(nWgts)", "Number of Events");
42  wgtVal_ = aDqmHelper.book1dHisto("wgtVal", "weights", 100, -1.5, 3., "weight", "Number of Weights");
43  bookTemplates(aDqmHelper,
45  "leadLepPt",
46  "leading lepton Pt",
48  0.,
50  "Pt_{l} (GeV)",
51  "Number of Events");
52  bookTemplates(aDqmHelper,
54  "leadLepEta",
55  "leading lepton #eta",
59  "#eta_{l}",
60  "Number of Events");
61  bookTemplates(aDqmHelper,
63  "JetMultiplicity",
64  "Gen jet multiplicity",
65  nJetsNbin_,
66  0,
67  nJetsNbin_,
68  "n",
69  "Number of Events");
70  bookTemplates(aDqmHelper,
72  "leadJetPt",
73  "leading Gen jet Pt",
74  jetPtNbin_,
75  0.,
77  "Pt_{j} (GeV)",
78  "Number of Events");
79  bookTemplates(aDqmHelper,
81  "leadJetEta",
82  "leading Gen jet #eta",
86  "#eta_{j}",
87  "Number of Events");
88 
89  return;
90 }
const double leadLepPtRange_
std::vector< MonitorElement * > leadLepEtaTemp_
MonitorElement * nEvt_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
std::vector< MonitorElement * > jetMultTemp_
MonitorElement * nlogWgt_
std::vector< MonitorElement * > leadJetEtaTemp_
std::vector< MonitorElement * > leadJetPtTemp_
std::vector< MonitorElement * > leadLepPtTemp_
MonitorElement * wgtVal_
void bookTemplates(DQMHelper &aDqmHelper, std::vector< MonitorElement *> &tmps, const std::string &name, const std::string &title, int nbin, float low, float high, const std::string &xtitle, const std::string &ytitle)

◆ bookTemplates()

void GenWeightValidation::bookTemplates ( DQMHelper aDqmHelper,
std::vector< MonitorElement *> &  tmps,
const std::string &  name,
const std::string &  title,
int  nbin,
float  low,
float  high,
const std::string &  xtitle,
const std::string &  ytitle 
)
private

Definition at line 92 of file GenWeightValidation.cc.

References DQMHelper::book1dHisto(), LaserClient_cfi::high, LaserClient_cfi::low, Skims_PA_cff::name, runGCPTkAlMap::title, compareTotals::xtitle, and compareTotals::ytitle.

Referenced by bookHistograms().

100  {
101  tmps.push_back(aDqmHelper.book1dHisto(name, title, nbin, low, high, xtitle, ytitle));
102  tmps.push_back(aDqmHelper.book1dHisto(name + "FSRup", title + " FSR up", nbin, low, high, xtitle, ytitle));
103  tmps.push_back(aDqmHelper.book1dHisto(name + "FSRdn", title + " FSR down", nbin, low, high, xtitle, ytitle));
104  tmps.push_back(aDqmHelper.book1dHisto(
105  name + "FSRup_ratio", "Ratio of " + title + " FSR up / Nominal", nbin, low, high, xtitle, ytitle));
106  tmps.at(3)->setEfficiencyFlag();
107  tmps.push_back(aDqmHelper.book1dHisto(
108  name + "FSRdn_ratio", "Ratio of " + title + " FSR down / Nominal", nbin, low, high, xtitle, ytitle));
109  tmps.at(4)->setEfficiencyFlag();
110  tmps.push_back(aDqmHelper.book1dHisto(name + "ISRup", title + " ISR up", nbin, low, high, xtitle, ytitle));
111  tmps.push_back(aDqmHelper.book1dHisto(name + "ISRdn", title + " ISR down", nbin, low, high, xtitle, ytitle));
112  tmps.push_back(aDqmHelper.book1dHisto(
113  name + "ISRup_ratio", "Ratio of " + title + " ISR up / Nominal", nbin, low, high, xtitle, ytitle));
114  tmps.at(7)->setEfficiencyFlag();
115  tmps.push_back(aDqmHelper.book1dHisto(
116  name + "ISRdn_ratio", "Ratio of " + title + " ISR down / Nominal", nbin, low, high, xtitle, ytitle));
117  tmps.at(8)->setEfficiencyFlag();
118 } // to get ratio plots correctly - need to modify PostProcessor_cff.py as well!
MonitorElement * book1dHisto(const std::string &name, const std::string &title, int n, double xmin, double xmax, const std::string &xaxis, const std::string &yaxis)
Definition: DQMHelper.cc:7

◆ dqmBeginRun()

void GenWeightValidation::dqmBeginRun ( const edm::Run ,
const edm::EventSetup  
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 120 of file GenWeightValidation.cc.

120 {}

◆ fillTemplates()

void GenWeightValidation::fillTemplates ( std::vector< MonitorElement *> &  tmps,
float  obs 
)
private

Definition at line 182 of file GenWeightValidation.cc.

References idxFSRdown_, idxFSRup_, idxGenEvtInfo_, idxISRdown_, idxISRup_, weight_, and weights_.

Referenced by analyze().

182  {
183  tmps.at(0)->Fill(obs, weight_);
184  tmps.at(1)->Fill(obs, weights_.at(idxGenEvtInfo_)[idxFSRup_] / weights_.at(idxGenEvtInfo_)[1]);
185  tmps.at(2)->Fill(obs, weights_.at(idxGenEvtInfo_)[idxFSRdown_] / weights_.at(idxGenEvtInfo_)[1]);
186  tmps.at(5)->Fill(obs, weights_.at(idxGenEvtInfo_)[idxISRup_] / weights_.at(idxGenEvtInfo_)[1]);
187  tmps.at(6)->Fill(obs, weights_.at(idxGenEvtInfo_)[idxISRdown_] / weights_.at(idxGenEvtInfo_)[1]);
188 }
std::vector< std::vector< double > > weights_

Member Data Documentation

◆ genJetToken_

const edm::EDGetTokenT<reco::GenJetCollection> GenWeightValidation::genJetToken_
private

Definition at line 60 of file GenWeightValidation.h.

Referenced by analyze().

◆ genParticleToken_

const edm::EDGetTokenT<reco::GenParticleCollection> GenWeightValidation::genParticleToken_
private

Definition at line 59 of file GenWeightValidation.h.

Referenced by analyze().

◆ idxFSRdown_

const int GenWeightValidation::idxFSRdown_
private

Definition at line 62 of file GenWeightValidation.h.

Referenced by fillTemplates(), and GenWeightValidation().

◆ idxFSRup_

const int GenWeightValidation::idxFSRup_
private

Definition at line 62 of file GenWeightValidation.h.

Referenced by fillTemplates(), and GenWeightValidation().

◆ idxGenEvtInfo_

const int GenWeightValidation::idxGenEvtInfo_
private

Definition at line 62 of file GenWeightValidation.h.

Referenced by analyze(), and fillTemplates().

◆ idxISRdown_

const int GenWeightValidation::idxISRdown_
private

Definition at line 62 of file GenWeightValidation.h.

Referenced by fillTemplates(), and GenWeightValidation().

◆ idxISRup_

const int GenWeightValidation::idxISRup_
private

Definition at line 62 of file GenWeightValidation.h.

Referenced by fillTemplates(), and GenWeightValidation().

◆ idxMax_

int GenWeightValidation::idxMax_
private

Definition at line 66 of file GenWeightValidation.h.

Referenced by analyze(), and GenWeightValidation().

◆ jetEtaCut_

const double GenWeightValidation::jetEtaCut_
private

Definition at line 65 of file GenWeightValidation.h.

Referenced by analyze().

◆ jetMultTemp_

std::vector<MonitorElement *> GenWeightValidation::jetMultTemp_
private

Definition at line 55 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ jetPtCut_

const double GenWeightValidation::jetPtCut_
private

Definition at line 65 of file GenWeightValidation.h.

Referenced by analyze().

◆ jetPtNbin_

const int GenWeightValidation::jetPtNbin_
private

Definition at line 64 of file GenWeightValidation.h.

Referenced by bookHistograms().

◆ jetPtRange_

const double GenWeightValidation::jetPtRange_
private

Definition at line 65 of file GenWeightValidation.h.

Referenced by bookHistograms().

◆ leadJetEtaTemp_

std::vector<MonitorElement *> GenWeightValidation::leadJetEtaTemp_
private

Definition at line 57 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ leadJetPtTemp_

std::vector<MonitorElement *> GenWeightValidation::leadJetPtTemp_
private

Definition at line 56 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ leadLepEtaTemp_

std::vector<MonitorElement *> GenWeightValidation::leadLepEtaTemp_
private

Definition at line 54 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ leadLepPtCut_

const double GenWeightValidation::leadLepPtCut_
private

Definition at line 63 of file GenWeightValidation.h.

Referenced by analyze().

◆ leadLepPtNbin_

const int GenWeightValidation::leadLepPtNbin_
private

Definition at line 62 of file GenWeightValidation.h.

Referenced by bookHistograms().

◆ leadLepPtRange_

const double GenWeightValidation::leadLepPtRange_
private

Definition at line 63 of file GenWeightValidation.h.

Referenced by bookHistograms().

◆ leadLepPtTemp_

std::vector<MonitorElement *> GenWeightValidation::leadLepPtTemp_
private

Definition at line 53 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ lepEtaCut_

const double GenWeightValidation::lepEtaCut_
private

Definition at line 63 of file GenWeightValidation.h.

Referenced by analyze().

◆ nEvt_

MonitorElement* GenWeightValidation::nEvt_
private

Definition at line 50 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ nJetsNbin_

const int GenWeightValidation::nJetsNbin_
private

Definition at line 64 of file GenWeightValidation.h.

Referenced by bookHistograms().

◆ nlogWgt_

MonitorElement* GenWeightValidation::nlogWgt_
private

Definition at line 51 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ rapidityNbin_

const int GenWeightValidation::rapidityNbin_
private

Definition at line 62 of file GenWeightValidation.h.

Referenced by bookHistograms().

◆ rapidityRange_

const double GenWeightValidation::rapidityRange_
private

Definition at line 63 of file GenWeightValidation.h.

Referenced by bookHistograms().

◆ weight_

double GenWeightValidation::weight_
private

Definition at line 47 of file GenWeightValidation.h.

Referenced by analyze(), and fillTemplates().

◆ weights_

std::vector<std::vector<double> > GenWeightValidation::weights_
private

Definition at line 48 of file GenWeightValidation.h.

Referenced by analyze(), and fillTemplates().

◆ wgtVal_

MonitorElement* GenWeightValidation::wgtVal_
private

Definition at line 52 of file GenWeightValidation.h.

Referenced by analyze(), and bookHistograms().

◆ wmanager_

WeightManager GenWeightValidation::wmanager_
private

Definition at line 45 of file GenWeightValidation.h.

Referenced by analyze().