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 Attributes
RivetAnalyzer Class Reference

#include <RivetAnalyzer.h>

Inheritance diagram for RivetAnalyzer:
edm::EDAnalyzer

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 ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

Rivet::AnalysisHandler _analysisHandler
 
edm::InputTag _hepmcCollection
 
bool _isFirstEvent
 
std::string _outFileName
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 16 of file RivetAnalyzer.h.

Constructor & Destructor Documentation

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

Definition at line 22 of file RivetAnalyzer.cc.

References _analysisHandler, _hepmcCollection, and edm::ParameterSet::getParameter().

22  :
24 _isFirstEvent(true),
25 _outFileName(pset.getParameter<std::string>("OutputFile"))
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 = pset.getParameter<edm::InputTag>("HepMCCollection");
31 
32  //get the analyses
33  _analysisHandler.addAnalyses(analysisNames);
34 
35  //go through the analyses and check those that need the cross section
36  const std::set< AnaHandle, AnaHandleLess > & analyses = _analysisHandler.analyses();
37 
38  std::set< AnaHandle, AnaHandleLess >::const_iterator ibeg = analyses.begin();
39  std::set< AnaHandle, AnaHandleLess >::const_iterator iend = analyses.end();
40  std::set< AnaHandle, AnaHandleLess >::const_iterator iana;
41  double xsection = -1.;
42  xsection = pset.getParameter<double>("CrossSection");
43  for (iana = ibeg; iana != iend; ++iana){
44  if ((*iana)->needsCrossSection())
45  (*iana)->setCrossSection(xsection);
46  }
47 }
T getParameter(std::string const &) const
edm::InputTag _hepmcCollection
Definition: RivetAnalyzer.h:35
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:36
std::string _outFileName
Definition: RivetAnalyzer.h:38
RivetAnalyzer::~RivetAnalyzer ( )
virtual

Definition at line 49 of file RivetAnalyzer.cc.

49  {
50 }

Member Function Documentation

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

67  {
68 
69  //get the hepmc product from the event
71  iEvent.getByLabel(_hepmcCollection, evt);
72 
73  // get HepMC GenEvent
74  const HepMC::GenEvent *myGenEvent = evt->GetEvent();
75 
76  //aaply the beams initialization on the first event
77  if (_isFirstEvent){
78  _analysisHandler.init(*myGenEvent);
79  _isFirstEvent = false;
80  }
81 
82  //run the analysis
83  _analysisHandler.analyze(*myGenEvent);
84 }
edm::InputTag _hepmcCollection
Definition: RivetAnalyzer.h:35
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:36
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
void RivetAnalyzer::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 52 of file RivetAnalyzer.cc.

52  {
53  //set the environment, very ugly but rivet is monolithic when it comes to paths
54  char * cmsswbase = getenv("CMSSW_BASE");
55  char * cmsswrelease = getenv("CMSSW_RELEASE_BASE");
56  std::string rivetref, rivetinfo;
57  rivetref = "RIVET_REF_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
58  rivetinfo = "RIVET_INFO_PATH=" + string(cmsswbase) + "/src/GeneratorInterface/RivetInterface/data:" + string(cmsswrelease) + "/src/GeneratorInterface/RivetInterface/data";
59  putenv(strdup(rivetref.c_str()));
60  putenv(strdup(rivetinfo.c_str()));
61 }
void RivetAnalyzer::beginRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 63 of file RivetAnalyzer.cc.

63  {
64  return;
65 }
void RivetAnalyzer::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 91 of file RivetAnalyzer.cc.

References _analysisHandler, and _outFileName.

91  {
92  _analysisHandler.finalize();
93  _analysisHandler.writeData(_outFileName);
94 }
Rivet::AnalysisHandler _analysisHandler
Definition: RivetAnalyzer.h:36
std::string _outFileName
Definition: RivetAnalyzer.h:38
void RivetAnalyzer::endRun ( const edm::Run iRun,
const edm::EventSetup iSetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 87 of file RivetAnalyzer.cc.

87  {
88  return;
89 }

Member Data Documentation

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