CMS 3D CMS Logo

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

#include <RivetAnalyzer.h>

Inheritance diagram for RivetAnalyzer:
edm::one::EDAnalyzer< edm::one::WatchRuns, edm::one::WatchLuminosityBlocks, edm::one::SharedResources > edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Types

typedef dqm::legacy::DQMStore DQMStore
 
typedef dqm::legacy::MonitorElement MonitorElement
 
- Public Types inherited from edm::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void beginJob () override
 
void beginLuminosityBlock (const edm::LuminosityBlock &, const edm::EventSetup &) override
 
void beginRun (const edm::Run &, const edm::EventSetup &) override
 
void endJob () override
 List of registered analysis data objects. More...
 
void endLuminosityBlock (const edm::LuminosityBlock &, const edm::EventSetup &) override
 
void endRun (const edm::Run &, const edm::EventSetup &) override
 
 RivetAnalyzer (const edm::ParameterSet &)
 
 ~RivetAnalyzer () override
 
- Public Member Functions inherited from edm::one::EDAnalyzer< edm::one::WatchRuns, edm::one::WatchLuminosityBlocks, edm::one::SharedResources >
 EDAnalyzer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDAnalyzerBase () 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 Member Functions

void normalizeTree ()
 

Private Attributes

Rivet::AnalysisHandler _analysisHandler
 
bool _doFinalize
 
edm::EDGetTokenT< GenEventInfoProduct_genEventInfoCollection
 
edm::EDGetTokenT< GenLumiInfoHeader_genLumiInfoToken
 
int _GENweightNumber
 
edm::EDGetTokenT< edm::HepMCProduct_hepmcCollection
 
bool _isFirstEvent
 
edm::EDGetTokenT< LHEEventProduct_LHECollection
 
const edm::InputTag _lheLabel
 
edm::EDGetTokenT< LHERunInfoProduct_lheRunInfoToken
 
std::vector< std::string > _lheWeightNames
 
int _LHEweightNumber
 
std::vector< MonitorElement * > _mes
 
std::string _outFileName
 
bool _produceDQM
 
bool _useExternalWeight
 
bool _useGENweights
 
bool _useLHEweights
 
std::vector< std::string > _weightNames
 
double _xsection
 
DQMStoredbe
 

Additional Inherited Members

- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- 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)
 

Detailed Description

Definition at line 22 of file RivetAnalyzer.h.

Member Typedef Documentation

Definition at line 25 of file RivetAnalyzer.h.

Definition at line 26 of file RivetAnalyzer.h.

Constructor & Destructor Documentation

RivetAnalyzer::RivetAnalyzer ( const edm::ParameterSet pset)

Definition at line 19 of file RivetAnalyzer.cc.

References _analysisHandler, _genLumiInfoToken, _hepmcCollection, _LHECollection, _lheLabel, _lheRunInfoToken, _produceDQM, _useLHEweights, _xsection, dbe, edm::ParameterSet::getParameter(), and Utilities::operator.

20  : _analysisHandler(),
21  _isFirstEvent(true),
22  _outFileName(pset.getParameter<std::string>("OutputFile")),
23  //decide whether to finalize the plots or not.
24  //deciding not to finalize them can be useful for further harvesting of many jobs
25  _doFinalize(pset.getParameter<bool>("DoFinalize")),
26  _produceDQM(pset.getParameter<bool>("ProduceDQMOutput")),
27  _lheLabel(pset.getParameter<edm::InputTag>("LHECollection")),
28  _xsection(-1.) {
29  usesResource("Rivet");
30 
31  //retrive the analysis name from parameter set
32  std::vector<std::string> analysisNames = pset.getParameter<std::vector<std::string> >("AnalysisNames");
33 
34  _hepmcCollection = consumes<HepMCProduct>(pset.getParameter<edm::InputTag>("HepMCCollection"));
35  _genLumiInfoToken = consumes<GenLumiInfoHeader, edm::InLumi>(pset.getParameter<edm::InputTag>("genLumiInfo"));
36 
37  _useLHEweights = pset.getParameter<bool>("useLHEweights");
38  if (_useLHEweights) {
39  _lheRunInfoToken = consumes<LHERunInfoProduct, edm::InRun>(_lheLabel);
40  _LHECollection = consumes<LHEEventProduct>(_lheLabel);
41  }
42 
43  //get the analyses
44  _analysisHandler.addAnalyses(analysisNames);
45 
46  //set user cross section if needed
47  _xsection = pset.getParameter<double>("CrossSection");
48 
49  if (_produceDQM) {
50  // book stuff needed for DQM
51  dbe = nullptr;
53  }
54 }
T getParameter(std::string const &) const
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:59
const edm::InputTag _lheLabel
Definition: RivetAnalyzer.h:64
edm::EDGetTokenT< GenLumiInfoHeader > _genLumiInfoToken
Definition: RivetAnalyzer.h:57
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
Definition: RivetAnalyzer.h:49
edm::EDGetTokenT< LHERunInfoProduct > _lheRunInfoToken
Definition: RivetAnalyzer.h:58
double _xsection
Definition: RivetAnalyzer.h:65
DQMStore * dbe
Definition: RivetAnalyzer.h:69
std::string _outFileName
Definition: RivetAnalyzer.h:61
edm::EDGetTokenT< LHEEventProduct > _LHECollection
Definition: RivetAnalyzer.h:55
RivetAnalyzer::~RivetAnalyzer ( )
override

Definition at line 56 of file RivetAnalyzer.cc.

56 {}

Member Function Documentation

void RivetAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 115 of file RivetAnalyzer.cc.

References _analysisHandler, _hepmcCollection, _isFirstEvent, _LHECollection, _lheWeightNames, _useLHEweights, _weightNames, _xsection, edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), mps_fire::i, LHEEventProduct::originalXWGTUP(), AlCaHLTBitMon_QueryRunRegistry::string, and LHEEventProduct::weights().

115  {
116  //get the hepmc product from the event
118  iEvent.getByToken(_hepmcCollection, evt);
119 
120  // get HepMC GenEvent
121  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
122  std::unique_ptr<HepMC::GenEvent> tmpGenEvtPtr;
123  //if you want to use an external weight or set the cross section we have to clone the GenEvent and change the weight
124  if (_useLHEweights || _xsection > 0) {
125  tmpGenEvtPtr = std::make_unique<HepMC::GenEvent>(*(evt->GetEvent()));
126 
127  if (_xsection > 0) {
128  HepMC::GenCrossSection xsec;
129  xsec.set_cross_section(_xsection);
130  tmpGenEvtPtr->set_cross_section(xsec);
131  }
132 
133  if (_useLHEweights) {
134  std::vector<double> mergedWeights;
135  for (unsigned int i = 0; i < tmpGenEvtPtr->weights().size(); i++) {
136  mergedWeights.push_back(tmpGenEvtPtr->weights()[i]);
137  }
138 
139  edm::Handle<LHEEventProduct> lheEventHandle;
140  iEvent.getByToken(_LHECollection, lheEventHandle);
141  for (unsigned int i = 0; i < _lheWeightNames.size(); i++) {
142  mergedWeights.push_back(tmpGenEvtPtr->weights()[0] * lheEventHandle->weights().at(i).wgt /
143  lheEventHandle->originalXWGTUP());
144  }
145 
146  tmpGenEvtPtr->weights() = mergedWeights;
147  }
148  myGenEvent = tmpGenEvtPtr.get();
149  }
150 
151  //apply the beams initialization on the first event
152  if (_isFirstEvent) {
153  if (_useLHEweights) {
154  _weightNames.insert(_weightNames.end(), _lheWeightNames.begin(), _lheWeightNames.end());
155  }
156  // clean weight names to be accepted by Rivet plotting
157  std::vector<std::string> cleanedWeightNames;
158  for (std::string wn : _weightNames) {
159  cleanedWeightNames.push_back(std::regex_replace(wn, std::regex("[^A-Za-z\\d\\._=]"), "_"));
160  }
161  _analysisHandler.init(*myGenEvent, cleanedWeightNames);
162  const HepMC::GenCrossSection* xs = myGenEvent->cross_section();
163  _analysisHandler.setCrossSection(make_pair(xs->cross_section(), xs->cross_section_error()));
164 
165  _isFirstEvent = false;
166  }
167 
168  //run the analysis
169  _analysisHandler.analyze(*myGenEvent);
170 }
double originalXWGTUP() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const std::vector< WGT > & weights() const
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:59
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
Definition: RivetAnalyzer.h:49
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
double _xsection
Definition: RivetAnalyzer.h:65
std::vector< std::string > _lheWeightNames
Definition: RivetAnalyzer.h:67
edm::EDGetTokenT< LHEEventProduct > _LHECollection
Definition: RivetAnalyzer.h:55
std::vector< std::string > _weightNames
Definition: RivetAnalyzer.h:66
void RivetAnalyzer::beginJob ( void  )
overridevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 58 of file RivetAnalyzer.cc.

References SiStripCommissioningSource_FromEDM_cfg::cmsswbase, and AlCaHLTBitMon_QueryRunRegistry::string.

58  {
59  //set the environment, very ugly but rivet is monolithic when it comes to paths
60  char* cmsswbase = std::getenv("CMSSW_BASE");
61  char* cmsswrelease = std::getenv("CMSSW_RELEASE_BASE");
62  if (!std::getenv("RIVET_REF_PATH")) {
63  const std::string rivetref = "RIVET_REF_PATH=" + string(cmsswbase) +
64  "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) +
65  "/src/GeneratorInterface/RivetInterface/data";
66  char* rivetrefCstr = strdup(rivetref.c_str());
67  putenv(rivetrefCstr);
68  free(rivetrefCstr);
69  }
70  if (!std::getenv("RIVET_INFO_PATH")) {
71  const std::string rivetinfo = "RIVET_INFO_PATH=" + string(cmsswbase) +
72  "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) +
73  "/src/GeneratorInterface/RivetInterface/data";
74  char* rivetinfoCstr = strdup(rivetinfo.c_str());
75  putenv(rivetinfoCstr);
76  free(rivetinfoCstr);
77  }
78 }
void RivetAnalyzer::beginLuminosityBlock ( const edm::LuminosityBlock iLumi,
const edm::EventSetup iSetup 
)
override

Definition at line 101 of file RivetAnalyzer.cc.

References _genLumiInfoToken, _weightNames, edm::LuminosityBlock::getByToken(), and GenLumiInfoHeader::weightNames().

101  {
102  edm::Handle<GenLumiInfoHeader> genLumiInfoHandle;
103  iLumi.getByToken(_genLumiInfoToken, genLumiInfoHandle);
104 
105  _weightNames = genLumiInfoHandle->weightNames();
106 
107  // need to reset the default weight name (or plotting will fail)
108  if (!_weightNames.empty()) {
109  _weightNames[0] = "";
110  }
111 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
const std::vector< std::string > & weightNames() const
edm::EDGetTokenT< GenLumiInfoHeader > _genLumiInfoToken
Definition: RivetAnalyzer.h:57
std::vector< std::string > _weightNames
Definition: RivetAnalyzer.h:66
void RivetAnalyzer::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
override

Definition at line 80 of file RivetAnalyzer.cc.

References _lheLabel, _lheWeightNames, _useLHEweights, edm::Run::getByLabel(), LHERunInfoProduct::headers_begin(), LHERunInfoProduct::headers_end(), groupFilesInBlocks::lines, and match().

80  {
81  if (_useLHEweights) {
82  edm::Handle<LHERunInfoProduct> lheRunInfoHandle;
83  iRun.getByLabel(_lheLabel, lheRunInfoHandle);
84  typedef std::vector<LHERunInfoProduct::Header>::const_iterator headers_const_iterator;
85 
86  std::regex reg("<weight.*> ?(.*?) ?<\\/weight>");
87  for (headers_const_iterator iter = lheRunInfoHandle->headers_begin(); iter != lheRunInfoHandle->headers_end();
88  iter++) {
89  std::vector<std::string> lines = iter->lines();
90  for (unsigned int iLine = 0; iLine < lines.size(); iLine++) {
91  std::smatch match;
92  std::regex_search(lines.at(iLine), match, reg);
93  if (!match.empty()) {
94  _lheWeightNames.push_back(match[1]);
95  }
96  }
97  }
98  }
99 }
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:280
headers_const_iterator headers_end() const
const edm::InputTag _lheLabel
Definition: RivetAnalyzer.h:64
headers_const_iterator headers_begin() const
std::vector< std::string > _lheWeightNames
Definition: RivetAnalyzer.h:67
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
void RivetAnalyzer::endJob ( void  )
overridevirtual

List of registered analysis data objects.

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 191 of file RivetAnalyzer.cc.

191 {}
void RivetAnalyzer::endLuminosityBlock ( const edm::LuminosityBlock iLumi,
const edm::EventSetup iSetup 
)
override

Definition at line 113 of file RivetAnalyzer.cc.

113 { return; }
void RivetAnalyzer::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
override

Definition at line 172 of file RivetAnalyzer.cc.

References _analysisHandler, _doFinalize, and _outFileName.

172  {
173  if (_doFinalize)
174  _analysisHandler.finalize();
175  else {
176  //if we don't finalize we just want to do the transformation from histograms to DPS
178  //normalizeTree();
179  }
180  _analysisHandler.writeData(_outFileName);
181 
182  return;
183 }
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:59
std::string _outFileName
Definition: RivetAnalyzer.h:61
void RivetAnalyzer::normalizeTree ( )
private

Definition at line 193 of file RivetAnalyzer.cc.

References _analysisHandler, _mes, _produceDQM, hltExoticaValidator_cfi::analyses, dqm::dqmstoreimpl::DQMStore::book1D(), dbe, DEFINE_FWK_MODULE, nevent, and dqm::dqmstoreimpl::DQMStore::setCurrentFolder().

193  {
194  using namespace YODA;
195  std::vector<string> analyses = _analysisHandler.analysisNames();
196 
197  //tree.ls(".", true);
198  const string tmpdir = "/RivetNormalizeTmp";
199  //tree.mkdir(tmpdir);
200  for (const string& analysis : analyses) {
201  if (_produceDQM) {
202  dbe->setCurrentFolder("Rivet/" + analysis);
203  //global variables that are always present
204  //sumOfWeights
205  TH1F nevent("nEvt", "n analyzed Events", 1, 0., 1.);
206  nevent.SetBinContent(1, _analysisHandler.sumW());
207  _mes.push_back(dbe->book1D("nEvt", &nevent));
208  }
209  //cross section
210  //TH1F xsection("xSection", "Cross Section", 1, 0., 1.);
211  //xsection.SetBinContent(1,_analysisHandler.crossSection());
212  //_mes.push_back(dbe->book1D("xSection",&xsection));
213  //now loop over the histograms
214 
215  /*
216  const vector<string> paths = tree.listObjectNames("/"+analysis, true); // args set recursive listing
217  std::cout << "Number of objects in YODA tree for analysis " << analysis << " = " << paths.size() << std::endl;
218  foreach (const string& path, paths) {
219  IManagedObject* hobj = tree.find(path);
220  if (hobj) {
221  // Weird seg fault on SLC4 when trying to dyn cast an IProfile ptr to a IHistogram
222  // Fix by attempting to cast to IProfile first, only try IHistogram if it fails.
223  IHistogram1D* histo = 0;
224  IProfile1D* prof = dynamic_cast<IProfile1D*>(hobj);
225  if (!prof) histo = dynamic_cast<IHistogram1D*>(hobj);
226 
227  std::cout << "Converting histo " << path << " to DPS" << std::endl;
228  tree.mv(path, tmpdir);
229  const size_t lastslash = path.find_last_of("/");
230  const string basename = path.substr(lastslash+1, path.length() - (lastslash+1));
231  const string tmppath = tmpdir + "/" + basename;
232 
233  // If it's a normal histo:
234  if (histo) {
235  IHistogram1D* tmphisto = dynamic_cast<IHistogram1D*>(tree.find(tmppath));
236  if (tmphisto) {
237  _analysisHandler.datapointsetFactory().create(path, *tmphisto);
238  }
239  //now convert to root and then ME
240  //need aida2flat (from Rivet 1.X) & flat2root here
241  TH1F* h = aida2root<IHistogram1D, TH1F>(histo, basename);
242  if (_produceDQM)
243  _mes.push_back(dbe->book1D(h->GetName(), h));
244  delete h;
245  tree.rm(tmppath);
246  }
247  // If it's a profile histo:
248  else if (prof) {
249  IProfile1D* tmpprof = dynamic_cast<IProfile1D*>(tree.find(tmppath));
250  if (tmpprof) {
251  _analysisHandler.datapointsetFactory().create(path, *tmpprof);
252  }
253  //now convert to root and then ME
254  //need aida2flat (from Rivet 1.X) & flat2root here
255  TProfile* p = aida2root<IProfile1D, TProfile>(prof, basename);
256  if (_produceDQM)
257  _mes.push_back(dbe->bookProfile(p->GetName(), p));
258  delete p;
259  tree.rm(tmppath);
260  }
261  }
262  }
263  */
264  }
265  //tree.rmdir(tmpdir);
266 }
int nevent
Definition: AMPTWrapper.h:84
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:59
DQMStore * dbe
Definition: RivetAnalyzer.h:69
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:639
MonitorElement * book1D(char_string const &name, char_string const &title, int const nchX, double const lowX, double const highX)
Book 1D histogram.
Definition: DQMStore.cc:1121
std::vector< MonitorElement * > _mes
Definition: RivetAnalyzer.h:70

Member Data Documentation

Rivet::AnalysisHandler RivetAnalyzer::_analysisHandler
private

Definition at line 59 of file RivetAnalyzer.h.

Referenced by analyze(), endRun(), normalizeTree(), and RivetAnalyzer().

bool RivetAnalyzer::_doFinalize
private

Definition at line 62 of file RivetAnalyzer.h.

Referenced by endRun().

edm::EDGetTokenT<GenEventInfoProduct> RivetAnalyzer::_genEventInfoCollection
private

Definition at line 56 of file RivetAnalyzer.h.

edm::EDGetTokenT<GenLumiInfoHeader> RivetAnalyzer::_genLumiInfoToken
private

Definition at line 57 of file RivetAnalyzer.h.

Referenced by beginLuminosityBlock(), and RivetAnalyzer().

int RivetAnalyzer::_GENweightNumber
private

Definition at line 54 of file RivetAnalyzer.h.

edm::EDGetTokenT<edm::HepMCProduct> RivetAnalyzer::_hepmcCollection
private

Definition at line 49 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

bool RivetAnalyzer::_isFirstEvent
private

Definition at line 60 of file RivetAnalyzer.h.

Referenced by analyze().

edm::EDGetTokenT<LHEEventProduct> RivetAnalyzer::_LHECollection
private

Definition at line 55 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

const edm::InputTag RivetAnalyzer::_lheLabel
private

Definition at line 64 of file RivetAnalyzer.h.

Referenced by beginRun(), and RivetAnalyzer().

edm::EDGetTokenT<LHERunInfoProduct> RivetAnalyzer::_lheRunInfoToken
private

Definition at line 58 of file RivetAnalyzer.h.

Referenced by RivetAnalyzer().

std::vector<std::string> RivetAnalyzer::_lheWeightNames
private

Definition at line 67 of file RivetAnalyzer.h.

Referenced by analyze(), and beginRun().

int RivetAnalyzer::_LHEweightNumber
private

Definition at line 52 of file RivetAnalyzer.h.

std::vector<MonitorElement *> RivetAnalyzer::_mes
private

Definition at line 70 of file RivetAnalyzer.h.

Referenced by normalizeTree().

std::string RivetAnalyzer::_outFileName
private

Definition at line 61 of file RivetAnalyzer.h.

Referenced by endRun().

bool RivetAnalyzer::_produceDQM
private

Definition at line 63 of file RivetAnalyzer.h.

Referenced by normalizeTree(), and RivetAnalyzer().

bool RivetAnalyzer::_useExternalWeight
private

Definition at line 50 of file RivetAnalyzer.h.

bool RivetAnalyzer::_useGENweights
private

Definition at line 53 of file RivetAnalyzer.h.

bool RivetAnalyzer::_useLHEweights
private

Definition at line 51 of file RivetAnalyzer.h.

Referenced by analyze(), beginRun(), and RivetAnalyzer().

std::vector<std::string> RivetAnalyzer::_weightNames
private

Definition at line 66 of file RivetAnalyzer.h.

Referenced by analyze(), and beginLuminosityBlock().

double RivetAnalyzer::_xsection
private

Definition at line 65 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

DQMStore* RivetAnalyzer::dbe
private

Definition at line 69 of file RivetAnalyzer.h.

Referenced by normalizeTree(), and RivetAnalyzer().