|
||||
File indexing completed on 2024-09-07 04:37:17
0001 #ifndef PhysicsTools_FWLite_interface_FWLiteAnalyzerWrapper_h 0002 #define PhysicsTools_FWLite_interface_FWLiteAnalyzerWrapper_h 0003 0004 #include <string> 0005 #include <vector> 0006 #include <iostream> 0007 0008 #include <TFile.h> 0009 #include <TSystem.h> 0010 0011 #include "DataFormats/FWLite/interface/ChainEvent.h" 0012 #include "DataFormats/FWLite/interface/InputSource.h" 0013 #include "DataFormats/FWLite/interface/OutputFiles.h" 0014 #include "FWCore/FWLite/interface/FWLiteEnabler.h" 0015 #include "FWCore/ParameterSet/interface/ProcessDesc.h" 0016 #include "PhysicsTools/FWLite/interface/TFileService.h" 0017 0018 /** 0019 \class FWLiteAnalyzerWrapper FWLiteAnalyzerWrapper.h "PhysicsTools/FWLite/interface/FWLiteAnalyzerWrapper.h" 0020 \brief Wrapper class for classes of type BasicAnalyzer to "convert" them into a full a basic FWLiteAnalyzer 0021 0022 This template class is a wrapper round classes of type BasicAnalyzer as defined in in the 0023 BasicAnalyzer.h file of this package. From this class the wrapper expects the following 0024 member functions: 0025 0026 + a contructor with a const edm::ParameterSet& and a TFileDirectory& as input. 0027 + a beginJob function 0028 + a endJob function 0029 + a analyze function with an const edm::EventBase& as input 0030 0031 these functions are called within the wrapper. The wrapper translates the common class into 0032 a basic FWLiteAnalyzer, which provides basic functionality of reading configuration files 0033 and event looping. An example of an application given in the PatExamples package is shown 0034 below: 0035 0036 #include "PhysicsTools/PatExamples/interface/BasicMuonAnalyzer.h" 0037 #include "PhysicsTools/FWLite/interface/FWLiteAnalyzerWrapper.h" 0038 0039 typedef fwlite::AnalyzerWrapper<BasicMuonAnalyzer> WrappedFWLiteAnalyzer; 0040 0041 int main(int argc, char* argv[]) 0042 { 0043 // load framework libraries 0044 gSystem->Load( "libFWCoreFWLite" ); 0045 FWLiteEnabler::enable(); 0046 0047 // only allow one argument for this simple example which should be the 0048 // the python cfg file 0049 if ( argc < 2 ) { 0050 std::cout << "Usage : " << argv[0] << " [parameters.py]" << std::endl; 0051 return 0; 0052 } 0053 0054 // get the python configuration 0055 PythonProcessDesc builder(argv[1]); 0056 WrappedFWLiteAnalyzer ana(*(builder.processDesc()->getProcessPSet()), std::string("MuonAnalyzer"), std::string("analyzeBasicPat")); 0057 ana.beginJob(); 0058 ana.analyze(); 0059 ana.endJob(); 0060 return 0; 0061 } 0062 0063 The configuration file for this FWLiteAnalyzer is expected to have the following structure: 0064 0065 import FWCore.ParameterSet.Config as cms 0066 0067 process = cms.Process("FWLitePlots") 0068 0069 process.fwliteInput = cms.PSet( 0070 fileNames = cms.untracked.vstring('file:patTuple.root'), ## mandatory 0071 maxEvents = cms.int32(-1), ## optional 0072 outputEvery = cms.uint32(10), ## optional 0073 ) 0074 0075 process.fwliteOutput = cms.PSet( 0076 fileName = cms.untracked.string('outputHistos.root') ## mandatory 0077 ) 0078 0079 process.muonAnalyzer = cms.PSet( 0080 muons = cms.InputTag('cleanPatMuons') ## input for the simple example above 0081 ) 0082 0083 0084 where the parameters maxEvents and 0085 reportAfter are optional. If omitted all events in the file(s) will be looped and no progress 0086 report will be given. More input files can be given as a vector of strings. Potential histograms 0087 per default will be written directely into the file without any furhter directory structure. If 0088 the class is instantiated with an additional directory string a new directory with the 0089 corresponding name will be created and the histograms will be added to this directory. 0090 With these wrapper classes we have the use case in mind that you keep classes, which easily can 0091 be used both within the full framework and within FWLite. 0092 0093 NOTE: in the current implementation this wrapper class only supports basic event looping. For 0094 anytasks of more complexity we recommend you to start from a FWLiteAnalyzer class from the very 0095 beginning and just to stay within FWLite. 0096 */ 0097 0098 namespace fwlite { 0099 0100 template <class T> 0101 class AnalyzerWrapper { 0102 public: 0103 /// default constructor 0104 AnalyzerWrapper(const edm::ParameterSet& cfg, std::string analyzerName, std::string directory = ""); 0105 /// default destructor 0106 virtual ~AnalyzerWrapper() {} 0107 /// everything which has to be done before the event loop 0108 virtual void beginJob() { analyzer_->beginJob(); } 0109 /// everything which has to be done during the event loop. NOTE: the event will be looped inside this function 0110 virtual void analyze(); 0111 /// everything which has to be done after the event loop 0112 virtual void endJob() { analyzer_->endJob(); } 0113 0114 protected: 0115 /// helper class for input parameter handling 0116 fwlite::InputSource inputHandler_; 0117 /// helper class for output file handling 0118 fwlite::OutputFiles outputHandler_; 0119 /// maximal number of events to be processed (-1 means to loop over all event) 0120 int maxEvents_; 0121 /// number of events after which the progress will be reported (0 means no report) 0122 unsigned int reportAfter_; 0123 /// TFileService for histogram management 0124 fwlite::TFileService fileService_; 0125 /// derived class of type BasicAnalyzer 0126 std::shared_ptr<T> analyzer_; 0127 }; 0128 0129 /// default contructor 0130 template <class T> 0131 AnalyzerWrapper<T>::AnalyzerWrapper(const edm::ParameterSet& cfg, std::string analyzerName, std::string directory) 0132 : inputHandler_(cfg), 0133 outputHandler_(cfg), 0134 maxEvents_(inputHandler_.maxEvents()), 0135 reportAfter_(inputHandler_.reportAfter()), 0136 fileService_(outputHandler_.file()) { 0137 // analysis specific parameters 0138 const edm::ParameterSet& ana = cfg.getParameter<edm::ParameterSet>(analyzerName.c_str()); 0139 if (directory.empty()) { 0140 // create analysis class of type BasicAnalyzer 0141 analyzer_ = std::shared_ptr<T>(new T(ana, fileService_)); 0142 } else { 0143 // create a directory in the file if directory string is non empty 0144 TFileDirectory dir = fileService_.mkdir(directory); 0145 analyzer_ = std::shared_ptr<T>(new T(ana, dir)); 0146 } 0147 } 0148 0149 /// everything which has to be done during the event loop. NOTE: the event will be looped inside this function 0150 template <class T> 0151 void AnalyzerWrapper<T>::analyze() { 0152 int ievt = 0; 0153 std::vector<std::string> const& inputFiles = inputHandler_.files(); 0154 // loop the vector of input files 0155 fwlite::ChainEvent event(inputFiles); 0156 for (event.toBegin(); !event.atEnd(); ++event, ++ievt) { 0157 // break loop if maximal number of events is reached 0158 if (maxEvents_ > 0 ? ievt + 1 > maxEvents_ : false) 0159 break; 0160 // simple event counter 0161 if (reportAfter_ != 0 ? (ievt > 0 && ievt % reportAfter_ == 0) : false) 0162 std::cout << " processing event: " << ievt << std::endl; 0163 // analyze event 0164 analyzer_->analyze(event); 0165 } 0166 } 0167 } // namespace fwlite 0168 0169 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.2.1 LXR engine. The LXR team |