CMS 3D CMS Logo

RivetAnalyzer.cc
Go to the documentation of this file.
2 
6 
9 
10 #include "Rivet/Run.hh"
11 #include "Rivet/AnalysisHandler.hh"
12 #include "Rivet/Analysis.hh"
13 
14 #include <regex>
15 
16 using namespace Rivet;
17 using namespace edm;
18 
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;
52  dbe = edm::Service<DQMStore>().operator->();
53  }
54 }
55 
57 
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 }
79 
80 void RivetAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
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 }
100 
102  edm::Handle<GenLumiInfoHeader> genLumiInfoHandle;
103  if (iLumi.getByToken(_genLumiInfoToken, genLumiInfoHandle)) {
104  _weightNames = genLumiInfoHandle->weightNames();
105  }
106 
107  // need to reset the default weight name (or plotting will fail)
108  if (!_weightNames.empty()) {
109  _weightNames[0] = "";
110  } else { // Summer16 samples have 1 weight stored in HepMC but no weightNames
111  _weightNames.push_back("");
112  }
113 }
114 
115 void RivetAnalyzer::endLuminosityBlock(const edm::LuminosityBlock& iLumi, const edm::EventSetup& iSetup) { return; }
116 
118  //get the hepmc product from the event
120  iEvent.getByToken(_hepmcCollection, evt);
121 
122  // get HepMC GenEvent
123  const HepMC::GenEvent* myGenEvent = evt->GetEvent();
124  std::unique_ptr<HepMC::GenEvent> tmpGenEvtPtr;
125  //if you want to use an external weight or set the cross section we have to clone the GenEvent and change the weight
126  if (_useLHEweights || _xsection > 0) {
127  tmpGenEvtPtr = std::make_unique<HepMC::GenEvent>(*(evt->GetEvent()));
128 
129  if (_xsection > 0) {
130  HepMC::GenCrossSection xsec;
131  xsec.set_cross_section(_xsection);
132  tmpGenEvtPtr->set_cross_section(xsec);
133  }
134 
135  if (_useLHEweights) {
136  std::vector<double> mergedWeights;
137  for (unsigned int i = 0; i < tmpGenEvtPtr->weights().size(); i++) {
138  mergedWeights.push_back(tmpGenEvtPtr->weights()[i]);
139  }
140 
141  edm::Handle<LHEEventProduct> lheEventHandle;
142  iEvent.getByToken(_LHECollection, lheEventHandle);
143  for (unsigned int i = 0; i < _lheWeightNames.size(); i++) {
144  mergedWeights.push_back(tmpGenEvtPtr->weights()[0] * lheEventHandle->weights().at(i).wgt /
145  lheEventHandle->originalXWGTUP());
146  }
147 
148  tmpGenEvtPtr->weights() = mergedWeights;
149  }
150  myGenEvent = tmpGenEvtPtr.get();
151  }
152 
153  //apply the beams initialization on the first event
154  if (_isFirstEvent) {
155  if (_useLHEweights) {
156  _weightNames.insert(_weightNames.end(), _lheWeightNames.begin(), _lheWeightNames.end());
157  }
158  // clean weight names to be accepted by Rivet plotting
159  std::vector<std::string> cleanedWeightNames;
160  for (std::string wn : _weightNames) {
161  cleanedWeightNames.push_back(std::regex_replace(wn, std::regex("[^A-Za-z\\d\\._=]"), "_"));
162  }
163  _analysisHandler.init(*myGenEvent, cleanedWeightNames);
164  const HepMC::GenCrossSection* xs = myGenEvent->cross_section();
165  _analysisHandler.setCrossSection(make_pair(xs->cross_section(), xs->cross_section_error()));
166 
167  _isFirstEvent = false;
168  }
169 
170  //run the analysis
171  _analysisHandler.analyze(*myGenEvent);
172 }
173 
174 void RivetAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
175  if (_doFinalize)
176  _analysisHandler.finalize();
177  else {
178  //if we don't finalize we just want to do the transformation from histograms to DPS
180  //normalizeTree();
181  }
182  _analysisHandler.writeData(_outFileName);
183 
184  return;
185 }
186 
187 //from Rivet 2.X: Analysis.hh (cls 18Feb2014)
189 //const vector<AnalysisObjectPtr>& analysisObjects() const {
190 //return _analysisobjects;
191 //}
192 
194 
196  using namespace YODA;
197  std::vector<string> analyses = _analysisHandler.analysisNames();
198 
199  //tree.ls(".", true);
200  const string tmpdir = "/RivetNormalizeTmp";
201  //tree.mkdir(tmpdir);
202  for (const string& analysis : analyses) {
203  if (_produceDQM) {
204  dbe->setCurrentFolder("Rivet/" + analysis);
205  //global variables that are always present
206  //sumOfWeights
207  TH1F nevent("nEvt", "n analyzed Events", 1, 0., 1.);
208  nevent.SetBinContent(1, _analysisHandler.sumW());
209  _mes.push_back(dbe->book1D("nEvt", &nevent));
210  }
211  //cross section
212  //TH1F xsection("xSection", "Cross Section", 1, 0., 1.);
213  //xsection.SetBinContent(1,_analysisHandler.crossSection());
214  //_mes.push_back(dbe->book1D("xSection",&xsection));
215  //now loop over the histograms
216 
217  /*
218  const vector<string> paths = tree.listObjectNames("/"+analysis, true); // args set recursive listing
219  std::cout << "Number of objects in YODA tree for analysis " << analysis << " = " << paths.size() << std::endl;
220  foreach (const string& path, paths) {
221  IManagedObject* hobj = tree.find(path);
222  if (hobj) {
223  // Weird seg fault on SLC4 when trying to dyn cast an IProfile ptr to a IHistogram
224  // Fix by attempting to cast to IProfile first, only try IHistogram if it fails.
225  IHistogram1D* histo = 0;
226  IProfile1D* prof = dynamic_cast<IProfile1D*>(hobj);
227  if (!prof) histo = dynamic_cast<IHistogram1D*>(hobj);
228 
229  std::cout << "Converting histo " << path << " to DPS" << std::endl;
230  tree.mv(path, tmpdir);
231  const size_t lastslash = path.find_last_of("/");
232  const string basename = path.substr(lastslash+1, path.length() - (lastslash+1));
233  const string tmppath = tmpdir + "/" + basename;
234 
235  // If it's a normal histo:
236  if (histo) {
237  IHistogram1D* tmphisto = dynamic_cast<IHistogram1D*>(tree.find(tmppath));
238  if (tmphisto) {
239  _analysisHandler.datapointsetFactory().create(path, *tmphisto);
240  }
241  //now convert to root and then ME
242  //need aida2flat (from Rivet 1.X) & flat2root here
243  TH1F* h = aida2root<IHistogram1D, TH1F>(histo, basename);
244  if (_produceDQM)
245  _mes.push_back(dbe->book1D(h->GetName(), h));
246  delete h;
247  tree.rm(tmppath);
248  }
249  // If it's a profile histo:
250  else if (prof) {
251  IProfile1D* tmpprof = dynamic_cast<IProfile1D*>(tree.find(tmppath));
252  if (tmpprof) {
253  _analysisHandler.datapointsetFactory().create(path, *tmpprof);
254  }
255  //now convert to root and then ME
256  //need aida2flat (from Rivet 1.X) & flat2root here
257  TProfile* p = aida2root<IProfile1D, TProfile>(prof, basename);
258  if (_produceDQM)
259  _mes.push_back(dbe->bookProfile(p->GetName(), p));
260  delete p;
261  tree.rm(tmppath);
262  }
263  }
264  }
265  */
266  }
267  //tree.rmdir(tmpdir);
268 }
269 
RivetAnalyzer::normalizeTree
void normalizeTree()
Definition: RivetAnalyzer.cc:195
nevent
int nevent
Definition: AMPTWrapper.h:84
Handle.h
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
RivetAnalyzer::_lheRunInfoToken
edm::EDGetTokenT< LHERunInfoProduct > _lheRunInfoToken
Definition: RivetAnalyzer.h:58
RivetAnalyzer::beginLuminosityBlock
void beginLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override
Definition: RivetAnalyzer.cc:101
RivetAnalyzer::_xsection
double _xsection
Definition: RivetAnalyzer.h:65
RivetAnalyzer::_useLHEweights
bool _useLHEweights
Definition: RivetAnalyzer.h:51
RivetAnalyzer::beginJob
void beginJob() override
Definition: RivetAnalyzer.cc:58
edm::LuminosityBlock
Definition: LuminosityBlock.h:50
edm::Run
Definition: Run.h:45
LHEEventProduct::originalXWGTUP
double originalXWGTUP() const
Definition: LHEEventProduct.h:34
edm
HLT enums.
Definition: AlignableModifier.h:19
Rivet
Definition: RivetAnalysis.h:21
RivetAnalyzer::_doFinalize
bool _doFinalize
Definition: RivetAnalyzer.h:62
RivetAnalyzer::_genLumiInfoToken
edm::EDGetTokenT< GenLumiInfoHeader > _genLumiInfoToken
Definition: RivetAnalyzer.h:57
RivetAnalyzer::_isFirstEvent
bool _isFirstEvent
Definition: RivetAnalyzer.h:60
LHERunInfoProduct::headers_end
headers_const_iterator headers_end() const
Definition: LHERunInfoProduct.h:59
RivetAnalyzer::_lheWeightNames
std::vector< std::string > _lheWeightNames
Definition: RivetAnalyzer.h:67
edm::Handle
Definition: AssociativeIterator.h:50
RivetAnalyzer::beginRun
void beginRun(const edm::Run &, const edm::EventSetup &) override
Definition: RivetAnalyzer.cc:80
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
LHERunInfoProduct::headers_begin
headers_const_iterator headers_begin() const
Definition: LHERunInfoProduct.h:58
MakerMacros.h
edm::LuminosityBlock::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: LuminosityBlock.h:318
GenLumiInfoHeader::weightNames
const std::vector< std::string > & weightNames() const
Definition: GenLumiInfoHeader.h:25
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Service.h
Run.h
dqm::implementation::DQMStore::setCurrentFolder
void setCurrentFolder(std::string const &fullpath) override
Definition: DQMStore.h:569
RivetAnalyzer::RivetAnalyzer
RivetAnalyzer(const edm::ParameterSet &)
Definition: RivetAnalyzer.cc:19
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
funct::true
true
Definition: Factorize.h:173
RivetAnalyzer::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: RivetAnalyzer.cc:117
RivetAnalyzer::endJob
void endJob() override
List of registered analysis data objects.
Definition: RivetAnalyzer.cc:193
edm::ParameterSet
Definition: ParameterSet.h:36
groupFilesInBlocks.lines
lines
Definition: groupFilesInBlocks.py:95
Event.h
RivetAnalyzer::_analysisHandler
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:59
RivetAnalyzer::_hepmcCollection
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
Definition: RivetAnalyzer.h:49
RivetAnalyzer::_mes
std::vector< MonitorElement * > _mes
Definition: RivetAnalyzer.h:70
match
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
LHEEventProduct::weights
const std::vector< WGT > & weights() const
Definition: LHEEventProduct.h:35
edm::Service
Definition: Service.h:30
iEvent
int iEvent
Definition: GenABIO.cc:224
hltExoticaValidator_cfi.analyses
analyses
Definition: hltExoticaValidator_cfi.py:54
edm::EventSetup
Definition: EventSetup.h:57
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
RivetAnalyzer::_lheLabel
const edm::InputTag _lheLabel
Definition: RivetAnalyzer.h:64
RivetAnalyzer::dbe
DQMStore * dbe
Definition: RivetAnalyzer.h:69
RivetAnalyzer
Definition: RivetAnalyzer.h:22
std
Definition: JetResolutionObject.h:76
RivetAnalyzer::_weightNames
std::vector< std::string > _weightNames
Definition: RivetAnalyzer.h:66
RivetAnalyzer::_produceDQM
bool _produceDQM
Definition: RivetAnalyzer.h:63
RivetAnalyzer::~RivetAnalyzer
~RivetAnalyzer() override
Definition: RivetAnalyzer.cc:56
RivetAnalyzer::endRun
void endRun(const edm::Run &, const edm::EventSetup &) override
Definition: RivetAnalyzer.cc:174
RivetAnalyzer::_LHECollection
edm::EDGetTokenT< LHEEventProduct > _LHECollection
Definition: RivetAnalyzer.h:55
RivetAnalyzer::endLuminosityBlock
void endLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override
Definition: RivetAnalyzer.cc:115
edm::Run::getByLabel
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:280
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
RivetAnalyzer::_outFileName
std::string _outFileName
Definition: RivetAnalyzer.h:61
RivetAnalyzer.h
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27
dqm::implementation::IBooker::book1D
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
SiStripCommissioningSource_FromEDM_cfg.cmsswbase
cmsswbase
Definition: SiStripCommissioningSource_FromEDM_cfg.py:7