CMS 3D CMS Logo

RivetAnalyzer.cc
Go to the documentation of this file.
2 
5 
8 
9 #include "Rivet/AnalysisHandler.hh"
10 #include "Rivet/Analysis.hh"
11 
12 using namespace Rivet;
13 using namespace edm;
14 
16 _analysisHandler(),
17 _isFirstEvent(true),
18 _outFileName(pset.getParameter<std::string>("OutputFile")),
19 //decide whether to finlaize tthe plots or not.
20 //deciding not to finalize them can be useful for further harvesting of many jobs
21 _doFinalize(pset.getParameter<bool>("DoFinalize")),
22 _produceDQM(pset.getParameter<bool>("ProduceDQMOutput")),
23 _xsection(-1.)
24 {
25  usesResource("Rivet");
26 
27  //retrive the analysis name from paarmeter set
28  std::vector<std::string> analysisNames = pset.getParameter<std::vector<std::string> >("AnalysisNames");
29 
30  _hepmcCollection = consumes<HepMCProduct>(pset.getParameter<edm::InputTag>("HepMCCollection"));
31 
32  _useExternalWeight = pset.getParameter<bool>("UseExternalWeight");
33  if (_useExternalWeight) {
34  if (!pset.exists("GenEventInfoCollection")){
35  throw cms::Exception("RivetAnalyzer") << "when using an external event weight you have to specify the GenEventInfoProduct collection from which the weight has to be taken " ;
36  }
37 
38  _genEventInfoCollection = consumes<GenEventInfoProduct>(pset.getParameter<edm::InputTag>("GenEventInfoCollection"));
39  _useGENweights = pset.getParameter<bool>("useGENweights");
40  _GENweightNumber = pset.getParameter<int>("GENweightNumber");
41  _LHECollection = consumes<LHEEventProduct>(pset.getParameter<edm::InputTag>("LHECollection"));
42  _useLHEweights = pset.getParameter<bool>("useLHEweights");
43  _LHEweightNumber = pset.getParameter<int>("LHEweightNumber");
44 
45  }
46 
47  //get the analyses
48  _analysisHandler.addAnalyses(analysisNames);
49 
50  //go through the analyses and check those that need the cross section
51  const std::set< AnaHandle, CmpAnaHandle > & analyses = _analysisHandler.analyses();
52 
53  std::set< AnaHandle, CmpAnaHandle >::const_iterator ibeg = analyses.begin();
54  std::set< AnaHandle, CmpAnaHandle >::const_iterator iend = analyses.end();
55  std::set< AnaHandle, CmpAnaHandle >::const_iterator iana;
56  _xsection = pset.getParameter<double>("CrossSection");
57  for (iana = ibeg; iana != iend; ++iana){
58  if ((*iana)->needsCrossSection())
59  (*iana)->setCrossSection(_xsection);
60  }
61  if (_produceDQM){
62  // book stuff needed for DQM
63  dbe = nullptr;
65  dbe->setVerbose(50);
66  }
67 
68 }
69 
71 }
72 
74  //set the environment, very ugly but rivet is monolithic when it comes to paths
75  char * cmsswbase = getenv("CMSSW_BASE");
76  char * cmsswrelease = getenv("CMSSW_RELEASE_BASE");
77  if ( !getenv("RIVET_REF_PATH") )
78  {
79  const std::string rivetref = "RIVET_REF_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
80  putenv(strdup(rivetref.c_str()));
81  }
82  if ( !getenv("RIVET_INFO_PATH") )
83  {
84  const std::string rivetinfo = "RIVET_INFO_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
85  putenv(strdup(rivetinfo.c_str()));
86  }
87 }
88 
89 void RivetAnalyzer::beginRun(const edm::Run& iRun,const edm::EventSetup& iSetup){
90  return;
91 }
92 
94 
95  //get the hepmc product from the event
97  iEvent.getByToken(_hepmcCollection, evt);
98 
99  // get HepMC GenEvent
100  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
101  std::unique_ptr<HepMC::GenEvent> tmpGenEvtPtr;
102  //if you want to use an external weight or set the cross section we have to clone the GenEvent and change the weight
103  if ( _useExternalWeight || _xsection > 0 ){
104  tmpGenEvtPtr = std::make_unique<HepMC::GenEvent>(*(evt->GetEvent()));
105 
106  if (_xsection > 0){
107  HepMC::GenCrossSection xsec;
108  xsec.set_cross_section(_xsection);
109  tmpGenEvtPtr->set_cross_section(xsec);
110  }
111 
112  if ( _useExternalWeight ){
113  if (tmpGenEvtPtr->weights().empty()) {
114  throw cms::Exception("RivetAnalyzer") << "Original weight container has 0 size ";
115  }
116  if (tmpGenEvtPtr->weights().size() > 1) {
117  edm::LogWarning("RivetAnalyzer") << "Original event weight size is " << tmpGenEvtPtr->weights().size() << ". Will change only the first one ";
118  }
119 
120  double weightForRivet = 1.;
121 
122  if(_useGENweights){
123  edm::Handle<GenEventInfoProduct> genEventInfoProduct;
124  iEvent.getByToken(_genEventInfoCollection, genEventInfoProduct);
125  weightForRivet *= genEventInfoProduct->weights().at(_GENweightNumber);
126  }
127  if(_useLHEweights){
128  edm::Handle<LHEEventProduct> lheEventHandle;
129  iEvent.getByToken(_LHECollection,lheEventHandle);
130  const LHEEventProduct::WGT& wgt = lheEventHandle->weights().at(_LHEweightNumber);
131  weightForRivet *= wgt.wgt;
132  }
133 
134  tmpGenEvtPtr->weights()[0] = weightForRivet;
135  }
136  myGenEvent = tmpGenEvtPtr.get();
137 
138  }
139 
140 
141  //aaply the beams initialization on the first event
142  if (_isFirstEvent){
143  _analysisHandler.init(*myGenEvent);
144  _isFirstEvent = false;
145  }
146 
147  //run the analysis
148  _analysisHandler.analyze(*myGenEvent);
149 
150 }
151 
152 
153 void RivetAnalyzer::endRun(const edm::Run& iRun,const edm::EventSetup& iSetup){
154  if (_doFinalize)
155  _analysisHandler.finalize();
156  else {
157  //if we don't finalize we just want to do the transformation from histograms to DPS
159  //normalizeTree();
160 
161  }
162  _analysisHandler.writeData(_outFileName);
163 
164  return;
165 }
166 
167 
168 
169 //from Rivet 2.X: Analysis.hh (cls 18Feb2014)
171 //const vector<AnalysisObjectPtr>& analysisObjects() const {
172 //return _analysisobjects;
173 //}
174 
175 
176 
178 }
179 
180 
182  using namespace YODA;
183  std::vector<string> analyses = _analysisHandler.analysisNames();
184 
185  //tree.ls(".", true);
186  const string tmpdir = "/RivetNormalizeTmp";
187  //tree.mkdir(tmpdir);
188  foreach (const string& analysis, analyses) {
189  if (_produceDQM){
190  dbe->setCurrentFolder("Rivet/"+analysis);
191  //global variables that are always present
192  //sumOfWeights
193  TH1F nevent("nEvt", "n analyzed Events", 1, 0., 1.);
194  nevent.SetBinContent(1,_analysisHandler.sumOfWeights());
195  _mes.push_back(dbe->book1D("nEvt",&nevent));
196  }
197  //cross section
198  //TH1F xsection("xSection", "Cross Section", 1, 0., 1.);
199  //xsection.SetBinContent(1,_analysisHandler.crossSection());
200  //_mes.push_back(dbe->book1D("xSection",&xsection));
201  //now loop over the histograms
202 
203  /*
204  const vector<string> paths = tree.listObjectNames("/"+analysis, true); // args set recursive listing
205  std::cout << "Number of objects in YODA tree for analysis " << analysis << " = " << paths.size() << std::endl;
206  foreach (const string& path, paths) {
207  IManagedObject* hobj = tree.find(path);
208  if (hobj) {
209  // Weird seg fault on SLC4 when trying to dyn cast an IProfile ptr to a IHistogram
210  // Fix by attempting to cast to IProfile first, only try IHistogram if it fails.
211  IHistogram1D* histo = 0;
212  IProfile1D* prof = dynamic_cast<IProfile1D*>(hobj);
213  if (!prof) histo = dynamic_cast<IHistogram1D*>(hobj);
214 
215  std::cout << "Converting histo " << path << " to DPS" << std::endl;
216  tree.mv(path, tmpdir);
217  const size_t lastslash = path.find_last_of("/");
218  const string basename = path.substr(lastslash+1, path.length() - (lastslash+1));
219  const string tmppath = tmpdir + "/" + basename;
220 
221  // If it's a normal histo:
222  if (histo) {
223  IHistogram1D* tmphisto = dynamic_cast<IHistogram1D*>(tree.find(tmppath));
224  if (tmphisto) {
225  _analysisHandler.datapointsetFactory().create(path, *tmphisto);
226  }
227  //now convert to root and then ME
228  //need aida2flat (from Rivet 1.X) & flat2root here
229  TH1F* h = aida2root<IHistogram1D, TH1F>(histo, basename);
230  if (_produceDQM)
231  _mes.push_back(dbe->book1D(h->GetName(), h));
232  delete h;
233  tree.rm(tmppath);
234  }
235  // If it's a profile histo:
236  else if (prof) {
237  IProfile1D* tmpprof = dynamic_cast<IProfile1D*>(tree.find(tmppath));
238  if (tmpprof) {
239  _analysisHandler.datapointsetFactory().create(path, *tmpprof);
240  }
241  //now convert to root and then ME
242  //need aida2flat (from Rivet 1.X) & flat2root here
243  TProfile* p = aida2root<IProfile1D, TProfile>(prof, basename);
244  if (_produceDQM)
245  _mes.push_back(dbe->bookProfile(p->GetName(), p));
246  delete p;
247  tree.rm(tmppath);
248  }
249  }
250  }
251  */
252  }
253  //tree.rmdir(tmpdir);
254 
255 }
256 
257 
258 
void analyze(const edm::Event &, const edm::EventSetup &) override
T getParameter(std::string const &) const
void endRun(const edm::Run &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
int nevent
Definition: AMPTWrapper.h:74
const std::vector< WGT > & weights() const
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:50
int iEvent
Definition: GenABIO.cc:230
void beginJob() override
void beginRun(const edm::Run &, const edm::EventSetup &) override
RivetAnalyzer(const edm::ParameterSet &)
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
Definition: RivetAnalyzer.h:42
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
double _xsection
Definition: RivetAnalyzer.h:55
DQMStore * dbe
Definition: RivetAnalyzer.h:57
std::string _outFileName
Definition: RivetAnalyzer.h:52
HLT enums.
void endJob() override
List of registered analysis data objects.
edm::EDGetTokenT< LHEEventProduct > _LHECollection
Definition: RivetAnalyzer.h:48
std::vector< double > & weights()
bool _useExternalWeight
Definition: RivetAnalyzer.h:43
edm::EDGetTokenT< GenEventInfoProduct > _genEventInfoCollection
Definition: RivetAnalyzer.h:49
void normalizeTree()
Definition: Run.h:44
std::vector< MonitorElement * > _mes
Definition: RivetAnalyzer.h:58
~RivetAnalyzer() override