CMS 3D CMS Logo

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

#include <RivetAnalyzer.h>

Inheritance diagram for RivetAnalyzer:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
virtual void beginJob () override
 
virtual void beginRun (const edm::Run &, const edm::EventSetup &) override
 
virtual void endJob () override
 List of registered analysis data objects. More...
 
virtual void endRun (const edm::Run &, const edm::EventSetup &) override
 
 RivetAnalyzer (const edm::ParameterSet &)
 
virtual ~RivetAnalyzer ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
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
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void normalizeTree ()
 

Private Attributes

Rivet::AnalysisHandler _analysisHandler
 
bool _doFinalize
 
edm::EDGetTokenT< GenEventInfoProduct_genEventInfoCollection
 
edm::EDGetTokenT< edm::HepMCProduct_hepmcCollection
 
bool _isFirstEvent
 
edm::EDGetTokenT< LHEEventProduct_LHECollection
 
int _LHEweightNumber
 
std::vector< MonitorElement * > _mes
 
std::string _outFileName
 
bool _produceDQM
 
bool _useExternalWeight
 
bool _useLHEweights
 
double _xsection
 
DQMStoredbe
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- 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 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 21 of file RivetAnalyzer.h.

Constructor & Destructor Documentation

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

Definition at line 15 of file RivetAnalyzer.cc.

References _analysisHandler, _genEventInfoCollection, _hepmcCollection, _LHECollection, _LHEweightNumber, _produceDQM, _useExternalWeight, _useLHEweights, _xsection, dbe, Exception, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), and Utilities::operator.

15  :
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  //retrive the analysis name from paarmeter set
26  std::vector<std::string> analysisNames = pset.getParameter<std::vector<std::string> >("AnalysisNames");
27 
28  _hepmcCollection = consumes<HepMCProduct>(pset.getParameter<edm::InputTag>("HepMCCollection"));
29 
30  _useExternalWeight = pset.getParameter<bool>("UseExternalWeight");
31  if (_useExternalWeight) {
32  if (!pset.exists("GenEventInfoCollection")){
33  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 " ;
34  }
35  _LHECollection = consumes<LHEEventProduct>(pset.getParameter<edm::InputTag>("LHECollection"));
36  _useLHEweights = pset.getParameter<bool>("useLHEweights");
37  _LHEweightNumber = pset.getParameter<int>("LHEweightNumber");
38 
39  _genEventInfoCollection = consumes<GenEventInfoProduct>(pset.getParameter<edm::InputTag>("GenEventInfoCollection"));
40  _LHECollection = consumes<LHEEventProduct>(pset.getParameter<edm::InputTag>("LHECollection"));
41  _useLHEweights = pset.getParameter<bool>("useLHEweights");
42  _LHEweightNumber = pset.getParameter<int>("LHEweightNumber");
43 
44  }
45 
46  //get the analyses
47  _analysisHandler.addAnalyses(analysisNames);
48 
49  //go through the analyses and check those that need the cross section
50  const std::set< AnaHandle, CmpAnaHandle > & analyses = _analysisHandler.analyses();
51 
52  std::set< AnaHandle, CmpAnaHandle >::const_iterator ibeg = analyses.begin();
53  std::set< AnaHandle, CmpAnaHandle >::const_iterator iend = analyses.end();
54  std::set< AnaHandle, CmpAnaHandle >::const_iterator iana;
55  _xsection = pset.getParameter<double>("CrossSection");
56  for (iana = ibeg; iana != iend; ++iana){
57  if ((*iana)->needsCrossSection())
58  (*iana)->setCrossSection(_xsection);
59  }
60  if (_produceDQM){
61  // book stuff needed for DQM
62  dbe = 0;
64  dbe->setVerbose(50);
65  }
66 
67 }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:48
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
Definition: RivetAnalyzer.h:42
double _xsection
Definition: RivetAnalyzer.h:53
DQMStore * dbe
Definition: RivetAnalyzer.h:55
std::string _outFileName
Definition: RivetAnalyzer.h:50
edm::EDGetTokenT< LHEEventProduct > _LHECollection
Definition: RivetAnalyzer.h:46
bool _useExternalWeight
Definition: RivetAnalyzer.h:43
edm::EDGetTokenT< GenEventInfoProduct > _genEventInfoCollection
Definition: RivetAnalyzer.h:47
RivetAnalyzer::~RivetAnalyzer ( )
virtual

Definition at line 69 of file RivetAnalyzer.cc.

69  {
70 }

Member Function Documentation

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

Definition at line 92 of file RivetAnalyzer.cc.

References _analysisHandler, _genEventInfoCollection, _hepmcCollection, _isFirstEvent, _LHECollection, _LHEweightNumber, _useExternalWeight, _useLHEweights, _xsection, Exception, edm::Event::getByToken(), edm::HepMCProduct::GetEvent(), GenEventInfoProduct::weight(), LHEEventProduct::weights(), and gen::WeightsInfo::wgt.

92  {
93 
94  //get the hepmc product from the event
96  iEvent.getByToken(_hepmcCollection, evt);
97 
98  // get HepMC GenEvent
99  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
100  std::unique_ptr<HepMC::GenEvent> tmpGenEvtPtr;
101  //if you want to use an external weight or set the cross section we have to clone the GenEvent and change the weight
102  if ( _useExternalWeight || _xsection > 0 ){
103  tmpGenEvtPtr = std::make_unique<HepMC::GenEvent>(*(evt->GetEvent()));
104 
105  if (_xsection > 0){
106  HepMC::GenCrossSection xsec;
107  xsec.set_cross_section(_xsection);
108  tmpGenEvtPtr->set_cross_section(xsec);
109  }
110 
111  if ( _useExternalWeight ){
112  if (tmpGenEvtPtr->weights().size() == 0) {
113  throw cms::Exception("RivetAnalyzer") << "Original weight container has 0 size ";
114  }
115  if (tmpGenEvtPtr->weights().size() > 1) {
116  edm::LogWarning("RivetAnalyzer") << "Original event weight size is " << tmpGenEvtPtr->weights().size() << ". Will change only the first one ";
117  }
118 
119  if(!_useLHEweights){
120  edm::Handle<GenEventInfoProduct> genEventInfoProduct;
121  iEvent.getByToken(_genEventInfoCollection, genEventInfoProduct);
122  tmpGenEvtPtr->weights()[0] = genEventInfoProduct->weight();
123  }else{
124  edm::Handle<LHEEventProduct> lheEventHandle;
125  iEvent.getByToken(_LHECollection,lheEventHandle);
126  const LHEEventProduct::WGT& wgt = lheEventHandle->weights().at(_LHEweightNumber);
127  tmpGenEvtPtr->weights()[0] = wgt.wgt;
128  }
129  }
130  myGenEvent = tmpGenEvtPtr.get();
131 
132  }
133 
134 
135  //aaply the beams initialization on the first event
136  if (_isFirstEvent){
137  _analysisHandler.init(*myGenEvent);
138  _isFirstEvent = false;
139  }
140 
141  //run the analysis
142  _analysisHandler.analyze(*myGenEvent);
143 
144 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
double weight() const
const std::vector< WGT > & weights() const
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:48
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
Definition: RivetAnalyzer.h:42
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
double _xsection
Definition: RivetAnalyzer.h:53
edm::EDGetTokenT< LHEEventProduct > _LHECollection
Definition: RivetAnalyzer.h:46
bool _useExternalWeight
Definition: RivetAnalyzer.h:43
edm::EDGetTokenT< GenEventInfoProduct > _genEventInfoCollection
Definition: RivetAnalyzer.h:47
void RivetAnalyzer::beginJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 72 of file RivetAnalyzer.cc.

References AlCaHLTBitMon_QueryRunRegistry::string.

72  {
73  //set the environment, very ugly but rivet is monolithic when it comes to paths
74  char * cmsswbase = getenv("CMSSW_BASE");
75  char * cmsswrelease = getenv("CMSSW_RELEASE_BASE");
76  if ( !getenv("RIVET_REF_PATH") )
77  {
78  const std::string rivetref = "RIVET_REF_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
79  putenv(strdup(rivetref.c_str()));
80  }
81  if ( !getenv("RIVET_INFO_PATH") )
82  {
83  const std::string rivetinfo = "RIVET_INFO_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
84  putenv(strdup(rivetinfo.c_str()));
85  }
86 }
void RivetAnalyzer::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
overridevirtual

Definition at line 88 of file RivetAnalyzer.cc.

88  {
89  return;
90 }
void RivetAnalyzer::endJob ( void  )
overridevirtual

List of registered analysis data objects.

Reimplemented from edm::EDAnalyzer.

Definition at line 171 of file RivetAnalyzer.cc.

171  {
172 }
void RivetAnalyzer::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
overridevirtual

Definition at line 147 of file RivetAnalyzer.cc.

References _analysisHandler, _doFinalize, and _outFileName.

147  {
148  if (_doFinalize)
149  _analysisHandler.finalize();
150  else {
151  //if we don't finalize we just want to do the transformation from histograms to DPS
153  //normalizeTree();
154 
155  }
156  _analysisHandler.writeData(_outFileName);
157 
158  return;
159 }
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:48
std::string _outFileName
Definition: RivetAnalyzer.h:50
void RivetAnalyzer::normalizeTree ( )
private

Definition at line 175 of file RivetAnalyzer.cc.

References _analysisHandler, _mes, _produceDQM, dbe, DEFINE_FWK_MODULE, and nevent.

175  {
176  using namespace YODA;
177  std::vector<string> analyses = _analysisHandler.analysisNames();
178 
179  //tree.ls(".", true);
180  const string tmpdir = "/RivetNormalizeTmp";
181  //tree.mkdir(tmpdir);
182  foreach (const string& analysis, analyses) {
183  if (_produceDQM){
184  dbe->setCurrentFolder(("Rivet/"+analysis).c_str());
185  //global variables that are always present
186  //sumOfWeights
187  TH1F nevent("nEvt", "n analyzed Events", 1, 0., 1.);
188  nevent.SetBinContent(1,_analysisHandler.sumOfWeights());
189  _mes.push_back(dbe->book1D("nEvt",&nevent));
190  }
191  //cross section
192  //TH1F xsection("xSection", "Cross Section", 1, 0., 1.);
193  //xsection.SetBinContent(1,_analysisHandler.crossSection());
194  //_mes.push_back(dbe->book1D("xSection",&xsection));
195  //now loop over the histograms
196 
197  /*
198  const vector<string> paths = tree.listObjectNames("/"+analysis, true); // args set recursive listing
199  std::cout << "Number of objects in YODA tree for analysis " << analysis << " = " << paths.size() << std::endl;
200  foreach (const string& path, paths) {
201  IManagedObject* hobj = tree.find(path);
202  if (hobj) {
203  // Weird seg fault on SLC4 when trying to dyn cast an IProfile ptr to a IHistogram
204  // Fix by attempting to cast to IProfile first, only try IHistogram if it fails.
205  IHistogram1D* histo = 0;
206  IProfile1D* prof = dynamic_cast<IProfile1D*>(hobj);
207  if (!prof) histo = dynamic_cast<IHistogram1D*>(hobj);
208 
209  std::cout << "Converting histo " << path << " to DPS" << std::endl;
210  tree.mv(path, tmpdir);
211  const size_t lastslash = path.find_last_of("/");
212  const string basename = path.substr(lastslash+1, path.length() - (lastslash+1));
213  const string tmppath = tmpdir + "/" + basename;
214 
215  // If it's a normal histo:
216  if (histo) {
217  IHistogram1D* tmphisto = dynamic_cast<IHistogram1D*>(tree.find(tmppath));
218  if (tmphisto) {
219  _analysisHandler.datapointsetFactory().create(path, *tmphisto);
220  }
221  //now convert to root and then ME
222  //need aida2flat (from Rivet 1.X) & flat2root here
223  TH1F* h = aida2root<IHistogram1D, TH1F>(histo, basename);
224  if (_produceDQM)
225  _mes.push_back(dbe->book1D(h->GetName(), h));
226  delete h;
227  tree.rm(tmppath);
228  }
229  // If it's a profile histo:
230  else if (prof) {
231  IProfile1D* tmpprof = dynamic_cast<IProfile1D*>(tree.find(tmppath));
232  if (tmpprof) {
233  _analysisHandler.datapointsetFactory().create(path, *tmpprof);
234  }
235  //now convert to root and then ME
236  //need aida2flat (from Rivet 1.X) & flat2root here
237  TProfile* p = aida2root<IProfile1D, TProfile>(prof, basename);
238  if (_produceDQM)
239  _mes.push_back(dbe->bookProfile(p->GetName(), p));
240  delete p;
241  tree.rm(tmppath);
242  }
243  }
244  }
245  */
246  }
247  //tree.rmdir(tmpdir);
248 
249 }
int nevent
Definition: AMPTWrapper.h:74
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:48
DQMStore * dbe
Definition: RivetAnalyzer.h:55
std::vector< MonitorElement * > _mes
Definition: RivetAnalyzer.h:56

Member Data Documentation

Rivet::AnalysisHandler RivetAnalyzer::_analysisHandler
private

Definition at line 48 of file RivetAnalyzer.h.

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

bool RivetAnalyzer::_doFinalize
private

Definition at line 51 of file RivetAnalyzer.h.

Referenced by endRun().

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

Definition at line 47 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

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

Definition at line 42 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

bool RivetAnalyzer::_isFirstEvent
private

Definition at line 49 of file RivetAnalyzer.h.

Referenced by analyze().

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

Definition at line 46 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

int RivetAnalyzer::_LHEweightNumber
private

Definition at line 45 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

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

Definition at line 56 of file RivetAnalyzer.h.

Referenced by normalizeTree().

std::string RivetAnalyzer::_outFileName
private

Definition at line 50 of file RivetAnalyzer.h.

Referenced by endRun().

bool RivetAnalyzer::_produceDQM
private

Definition at line 52 of file RivetAnalyzer.h.

Referenced by normalizeTree(), and RivetAnalyzer().

bool RivetAnalyzer::_useExternalWeight
private

Definition at line 43 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

bool RivetAnalyzer::_useLHEweights
private

Definition at line 44 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

double RivetAnalyzer::_xsection
private

Definition at line 53 of file RivetAnalyzer.h.

Referenced by analyze(), and RivetAnalyzer().

DQMStore* RivetAnalyzer::dbe
private

Definition at line 55 of file RivetAnalyzer.h.

Referenced by normalizeTree(), and RivetAnalyzer().