CMS 3D CMS Logo

FWLiteAnalyzerWrapper.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_FWLite_interface_FWLiteAnalyzerWrapper_h
2 #define PhysicsTools_FWLite_interface_FWLiteAnalyzerWrapper_h
3 
4 #include <string>
5 #include <vector>
6 #include <iostream>
7 
8 #include <TFile.h>
9 #include <TSystem.h>
10 
17 
98 namespace fwlite {
99 
100  template <class T>
102  public:
106  virtual ~AnalyzerWrapper(){};
108  virtual void beginJob() { analyzer_->beginJob(); }
110  virtual void analyze();
112  virtual void endJob() { analyzer_->endJob(); }
113 
114  protected:
122  unsigned int reportAfter_;
126  std::shared_ptr<T> analyzer_;
127  };
128 
130  template <class T>
132  : inputHandler_(cfg),
133  outputHandler_(cfg),
134  maxEvents_(inputHandler_.maxEvents()),
135  reportAfter_(inputHandler_.reportAfter()),
136  fileService_(outputHandler_.file()) {
137  // analysis specific parameters
138  const edm::ParameterSet& ana = cfg.getParameter<edm::ParameterSet>(analyzerName.c_str());
139  if (directory.empty()) {
140  // create analysis class of type BasicAnalyzer
141  analyzer_ = std::shared_ptr<T>(new T(ana, fileService_));
142  } else {
143  // create a directory in the file if directory string is non empty
145  analyzer_ = std::shared_ptr<T>(new T(ana, dir));
146  }
147  }
148 
150  template <class T>
152  int ievt = 0;
153  std::vector<std::string> const& inputFiles = inputHandler_.files();
154  // loop the vector of input files
156  for (event.toBegin(); !event.atEnd(); ++event, ++ievt) {
157  // break loop if maximal number of events is reached
158  if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false)
159  break;
160  // simple event counter
161  if (reportAfter_ != 0 ? (ievt > 0 && ievt % reportAfter_ == 0) : false)
162  std::cout << " processing event: " << ievt << std::endl;
163  // analyze event
164  analyzer_->analyze(event);
165  }
166  }
167 } // namespace fwlite
168 
169 #endif
fwlite::OutputFiles outputHandler_
helper class for output file handling
fwlite::TFileService fileService_
TFileService for histogram management.
fwlite::InputSource inputHandler_
helper class for input parameter handling
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
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
virtual void analyze()
everything which has to be done during the event loop. NOTE: the event will be looped inside this fun...
int maxEvents_
maximal number of events to be processed (-1 means to loop over all event)
virtual ~AnalyzerWrapper()
default destructor
std::shared_ptr< T > analyzer_
derived class of type BasicAnalyzer
AnalyzerWrapper(const edm::ParameterSet &cfg, std::string analyzerName, std::string directory="")
default constructor
long double T
Definition: event.py:1
unsigned int reportAfter_
number of events after which the progress will be reported (0 means no report)