#include <RivetAnalyzer.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
virtual void | beginRun (const edm::Run &, const edm::EventSetup &) |
virtual void | endJob () |
virtual void | endRun (const edm::Run &, const edm::EventSetup &) |
RivetAnalyzer (const edm::ParameterSet &) | |
virtual | ~RivetAnalyzer () |
Private Attributes | |
Rivet::AnalysisHandler | _analysisHandler |
edm::InputTag | _hepmcCollection |
bool | _isFirstEvent |
std::string | _outFileName |
Definition at line 16 of file RivetAnalyzer.h.
RivetAnalyzer::RivetAnalyzer | ( | const edm::ParameterSet & | pset | ) |
Definition at line 22 of file RivetAnalyzer.cc.
References _analysisHandler, _hepmcCollection, and edm::ParameterSet::getParameter().
: _analysisHandler(), _isFirstEvent(true), _outFileName(pset.getParameter<std::string>("OutputFile")) { //retrive the analysis name from paarmeter set std::vector<std::string> analysisNames = pset.getParameter<std::vector<std::string> >("AnalysisNames"); _hepmcCollection = pset.getParameter<edm::InputTag>("HepMCCollection"); //get the analyses _analysisHandler.addAnalyses(analysisNames); //go through the analyses and check those that need the cross section const std::set< AnaHandle, AnaHandleLess > & analyses = _analysisHandler.analyses(); std::set< AnaHandle, AnaHandleLess >::const_iterator ibeg = analyses.begin(); std::set< AnaHandle, AnaHandleLess >::const_iterator iend = analyses.end(); std::set< AnaHandle, AnaHandleLess >::const_iterator iana; double xsection = -1.; xsection = pset.getParameter<double>("CrossSection"); for (iana = ibeg; iana != iend; ++iana){ if ((*iana)->needsCrossSection()) (*iana)->setCrossSection(xsection); } }
RivetAnalyzer::~RivetAnalyzer | ( | ) | [virtual] |
Definition at line 49 of file RivetAnalyzer.cc.
{ }
void RivetAnalyzer::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 67 of file RivetAnalyzer.cc.
References _analysisHandler, _hepmcCollection, _isFirstEvent, and edm::Event::getByLabel().
{ //get the hepmc product from the event edm::Handle<HepMCProduct> evt; iEvent.getByLabel(_hepmcCollection, evt); // get HepMC GenEvent const HepMC::GenEvent *myGenEvent = evt->GetEvent(); //aaply the beams initialization on the first event if (_isFirstEvent){ _analysisHandler.init(*myGenEvent); _isFirstEvent = false; } //run the analysis _analysisHandler.analyze(*myGenEvent); }
void RivetAnalyzer::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 52 of file RivetAnalyzer.cc.
{ //set the environment, very ugly but rivet is monolithic when it comes to paths char * cmsswbase = getenv("CMSSW_BASE"); char * cmsswrelease = getenv("CMSSW_RELEASE_BASE"); std::string rivetref, rivetinfo; rivetref = "RIVET_REF_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data"; rivetinfo = "RIVET_INFO_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data"; putenv(strdup(rivetref.c_str())); putenv(strdup(rivetinfo.c_str())); }
void RivetAnalyzer::beginRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
void RivetAnalyzer::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 91 of file RivetAnalyzer.cc.
References _analysisHandler, and _outFileName.
{ _analysisHandler.finalize(); _analysisHandler.writeData(_outFileName); }
void RivetAnalyzer::endRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Rivet::AnalysisHandler RivetAnalyzer::_analysisHandler [private] |
Definition at line 36 of file RivetAnalyzer.h.
Referenced by analyze(), endJob(), and RivetAnalyzer().
edm::InputTag RivetAnalyzer::_hepmcCollection [private] |
Definition at line 35 of file RivetAnalyzer.h.
Referenced by analyze(), and RivetAnalyzer().
bool RivetAnalyzer::_isFirstEvent [private] |
Definition at line 37 of file RivetAnalyzer.h.
Referenced by analyze().
std::string RivetAnalyzer::_outFileName [private] |
Definition at line 38 of file RivetAnalyzer.h.
Referenced by endJob().