#include <FWLiteAnalyzerWrapper.h>
Public Member Functions | |
virtual void | analyze () |
everything which has to be done during the event loop. NOTE: the event will be looped inside this function | |
AnalyzerWrapper (const edm::ParameterSet &cfg, std::string analyzerName, std::string directory="") | |
default constructor | |
virtual void | beginJob () |
everything which has to be done before the event loop | |
virtual void | endJob () |
everything which has to be done after the event loop | |
virtual | ~AnalyzerWrapper () |
default destructor | |
Protected Attributes | |
boost::shared_ptr< T > | analyzer_ |
derived class of type BasicAnalyzer | |
fwlite::TFileService | fileService_ |
TFileService for histogram management. | |
fwlite::InputSource | inputHandler_ |
helper class for input parameter handling | |
int | maxEvents_ |
maximal number of events to be processed (-1 means to loop over all event) | |
fwlite::OutputFiles | outputHandler_ |
helper class for output file handling | |
unsigned int | reportAfter_ |
number of events after which the progress will be reported (0 means no report) |
Definition at line 104 of file FWLiteAnalyzerWrapper.h.
fwlite::AnalyzerWrapper< T >::AnalyzerWrapper | ( | const edm::ParameterSet & | cfg, |
std::string | analyzerName, | ||
std::string | directory = "" |
||
) |
default constructor
default contructor
Definition at line 135 of file FWLiteAnalyzerWrapper.h.
References fwlite::AnalyzerWrapper< T >::analyzer_, dir, fwlite::AnalyzerWrapper< T >::fileService_, edm::ParameterSet::getParameter(), and TFileDirectory::mkdir().
: inputHandler_( cfg ), outputHandler_( cfg ), maxEvents_(inputHandler_.maxEvents()), reportAfter_(inputHandler_.reportAfter()), fileService_( outputHandler_.file() ) { // analysis specific parameters const edm::ParameterSet& ana = cfg.getParameter<edm::ParameterSet>(analyzerName.c_str()); if(directory.empty()){ // create analysis class of type BasicAnalyzer analyzer_ = boost::shared_ptr<T>( new T( ana, fileService_) ); } else{ // create a directory in the file if directory string is non empty TFileDirectory dir = fileService_.mkdir(directory.c_str()); analyzer_ = boost::shared_ptr<T>( new T( ana, dir ) ); } }
virtual fwlite::AnalyzerWrapper< T >::~AnalyzerWrapper | ( | ) | [inline, virtual] |
void fwlite::AnalyzerWrapper< T >::analyze | ( | void | ) | [virtual] |
everything which has to be done during the event loop. NOTE: the event will be looped inside this function
Definition at line 154 of file FWLiteAnalyzerWrapper.h.
References loadConditions::analyzer_, fwlite::ChainEvent::atEnd(), gather_cfg::cout, event(), merge::inputFiles, and fwlite::ChainEvent::toBegin().
Referenced by main().
{ int ievt=0; std::vector<std::string> const & inputFiles = inputHandler_.files(); // loop the vector of input files fwlite::ChainEvent event( inputFiles ); for(event.toBegin(); !event.atEnd(); ++event, ++ievt){ // break loop if maximal number of events is reached if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break; // simple event counter if(reportAfter_!=0 ? (ievt>0 && ievt%reportAfter_==0) : false) std::cout << " processing event: " << ievt << std::endl; // analyze event analyzer_->analyze(event); } }
virtual void fwlite::AnalyzerWrapper< T >::beginJob | ( | void | ) | [inline, virtual] |
everything which has to be done before the event loop
Definition at line 112 of file FWLiteAnalyzerWrapper.h.
References fwlite::AnalyzerWrapper< T >::analyzer_.
Referenced by main().
{ analyzer_->beginJob(); }
virtual void fwlite::AnalyzerWrapper< T >::endJob | ( | void | ) | [inline, virtual] |
everything which has to be done after the event loop
Definition at line 116 of file FWLiteAnalyzerWrapper.h.
References fwlite::AnalyzerWrapper< T >::analyzer_.
Referenced by main().
{ analyzer_->endJob(); }
boost::shared_ptr<T> fwlite::AnalyzerWrapper< T >::analyzer_ [protected] |
derived class of type BasicAnalyzer
Definition at line 130 of file FWLiteAnalyzerWrapper.h.
Referenced by fwlite::AnalyzerWrapper< T >::AnalyzerWrapper(), fwlite::AnalyzerWrapper< T >::beginJob(), and fwlite::AnalyzerWrapper< T >::endJob().
fwlite::TFileService fwlite::AnalyzerWrapper< T >::fileService_ [protected] |
TFileService for histogram management.
Definition at line 128 of file FWLiteAnalyzerWrapper.h.
Referenced by fwlite::AnalyzerWrapper< T >::AnalyzerWrapper().
fwlite::InputSource fwlite::AnalyzerWrapper< T >::inputHandler_ [protected] |
helper class for input parameter handling
Definition at line 120 of file FWLiteAnalyzerWrapper.h.
int fwlite::AnalyzerWrapper< T >::maxEvents_ [protected] |
maximal number of events to be processed (-1 means to loop over all event)
Definition at line 124 of file FWLiteAnalyzerWrapper.h.
fwlite::OutputFiles fwlite::AnalyzerWrapper< T >::outputHandler_ [protected] |
helper class for output file handling
Definition at line 122 of file FWLiteAnalyzerWrapper.h.
unsigned int fwlite::AnalyzerWrapper< T >::reportAfter_ [protected] |
number of events after which the progress will be reported (0 means no report)
Definition at line 126 of file FWLiteAnalyzerWrapper.h.