CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (std::string const &iProcessName, std::string const &iModuleLabel, bool iPrint, std::vector< char const * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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

65  {
66 }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 83 of file RivetAnalyzer.cc.

References _analysisHandler, _genEventInfoCollection, _hepmcCollection, _isFirstEvent, _LHECollection, _LHEweightNumber, _useExternalWeight, _useLHEweights, Exception, edm::Event::getByToken(), and gen::WeightsInfo::wgt.

83  {
84 
85  //get the hepmc product from the event
87  iEvent.getByToken(_hepmcCollection, evt);
88 
89  // get HepMC GenEvent
90  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
91 
92  //if you want to use an external weight we have to clone the GenEvent and change the weight
93  if ( _useExternalWeight ){
94 
95  HepMC::GenEvent * tmpGenEvtPtr = new HepMC::GenEvent( *(evt->GetEvent()) );
96  if (tmpGenEvtPtr->weights().size() == 0) {
97  throw cms::Exception("RivetAnalyzer") << "Original weight container has 0 size ";
98  }
99  if (tmpGenEvtPtr->weights().size() > 1) {
100  edm::LogWarning("RivetAnalyzer") << "Original event weight size is " << tmpGenEvtPtr->weights().size() << ". Will change only the first one ";
101  }
102 
103  if(!_useLHEweights){
104  edm::Handle<GenEventInfoProduct> genEventInfoProduct;
105  iEvent.getByToken(_genEventInfoCollection, genEventInfoProduct);
106  tmpGenEvtPtr->weights()[0] = genEventInfoProduct->weight();
107  }else{
108  edm::Handle<LHEEventProduct> lheEventHandle;
109  iEvent.getByToken(_LHECollection,lheEventHandle);
110  const LHEEventProduct::WGT& wgt = lheEventHandle->weights().at(_LHEweightNumber);
111  tmpGenEvtPtr->weights()[0] = wgt.wgt;
112  }
113  myGenEvent = tmpGenEvtPtr;
114 
115  }
116 
117 
118  //aaply the beams initialization on the first event
119  if (_isFirstEvent){
120  _analysisHandler.init(*myGenEvent);
121  _isFirstEvent = false;
122  }
123 
124  //run the analysis
125  _analysisHandler.analyze(*myGenEvent);
126 
127  //if we have cloned the GenEvent, we delete it
128  if ( _useExternalWeight )
129  delete myGenEvent;
130 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:48
edm::EDGetTokenT< edm::HepMCProduct > _hepmcCollection
Definition: RivetAnalyzer.h:42
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 68 of file RivetAnalyzer.cc.

References AlCaHLTBitMon_QueryRunRegistry::string.

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

Reimplemented from edm::EDAnalyzer.

Definition at line 79 of file RivetAnalyzer.cc.

79  {
80  return;
81 }
void RivetAnalyzer::endJob ( void  )
overridevirtual

List of registered analysis data objects.

Reimplemented from edm::EDAnalyzer.

Definition at line 157 of file RivetAnalyzer.cc.

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

Reimplemented from edm::EDAnalyzer.

Definition at line 133 of file RivetAnalyzer.cc.

References _analysisHandler, _doFinalize, and _outFileName.

133  {
134  if (_doFinalize)
135  _analysisHandler.finalize();
136  else {
137  //if we don't finalize we just want to do the transformation from histograms to DPS
139  //normalizeTree();
140 
141  }
142  _analysisHandler.writeData(_outFileName);
143 
144  return;
145 }
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:48
std::string _outFileName
Definition: RivetAnalyzer.h:50
void RivetAnalyzer::normalizeTree ( )
private

Definition at line 161 of file RivetAnalyzer.cc.

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

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

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 55 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().

DQMStore* RivetAnalyzer::dbe
private

Definition at line 54 of file RivetAnalyzer.h.

Referenced by normalizeTree(), and RivetAnalyzer().