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;
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  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 }
112 
113 void RivetAnalyzer::endLuminosityBlock(const edm::LuminosityBlock& iLumi, const edm::EventSetup& iSetup) { return; }
114 
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 }
171 
172 void RivetAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
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 }
184 
185 //from Rivet 2.X: Analysis.hh (cls 18Feb2014)
187 //const vector<AnalysisObjectPtr>& analysisObjects() const {
188 //return _analysisobjects;
189 //}
190 
192 
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 }
267 
void analyze(const edm::Event &, const edm::EventSetup &) override
double originalXWGTUP() const
T getParameter(std::string const &) const
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:280
void endRun(const edm::Run &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
bool getByToken(EDGetToken token, Handle< PROD > &result) const
const std::vector< std::string > & weightNames() const
headers_const_iterator headers_end() const
int nevent
Definition: AMPTWrapper.h:84
const std::vector< WGT > & weights() const
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:59
const edm::InputTag _lheLabel
Definition: RivetAnalyzer.h:64
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void beginLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override
headers_const_iterator headers_begin() const
void beginJob() override
void endLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override
void beginRun(const edm::Run &, const edm::EventSetup &) override
RivetAnalyzer(const edm::ParameterSet &)
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
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
double _xsection
Definition: RivetAnalyzer.h:65
DQMStore * dbe
Definition: RivetAnalyzer.h:69
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:639
std::vector< std::string > _lheWeightNames
Definition: RivetAnalyzer.h:67
std::string _outFileName
Definition: RivetAnalyzer.h:61
HLT enums.
void endJob() override
List of registered analysis data objects.
edm::EDGetTokenT< LHEEventProduct > _LHECollection
Definition: RivetAnalyzer.h:55
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::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< std::string > _weightNames
Definition: RivetAnalyzer.h:66
void normalizeTree()
Definition: Run.h:45
std::vector< MonitorElement * > _mes
Definition: RivetAnalyzer.h:70
~RivetAnalyzer() override