Hall A ROOT/C++ Analyzer (podd)
Loading...
Searching...
No Matches
THaAnalyzer.h
Go to the documentation of this file.
1#ifndef Podd_THaAnalyzer_h_
2#define Podd_THaAnalyzer_h_
3
5//
6// THaAnalyzer
7//
9
10#include "TObject.h"
11#include "TString.h"
12#include <vector>
13
14class THaEvent;
15class THaRunBase;
16class THaOutput;
17class TList;
18class TIter;
19class TFile;
20class TDatime;
21class THaCut;
22class THaBenchmark;
23class THaEvData;
24class THaPostProcess;
25class THaCrateMap;
28class THaApparatus;
29class THaSpectrometer;
32namespace Podd {
33 class InterStageModule;
34}
35
36class THaAnalyzer : public TObject {
37
38public:
40 virtual ~THaAnalyzer();
41
43 virtual Int_t AddPostProcess( THaPostProcess* module );
44 virtual void Close();
45 virtual Int_t Init( THaRunBase* run );
46 Int_t Init( THaRunBase& run ) { return Init( &run ); }
47 virtual Int_t Process( THaRunBase* run=nullptr );
48 Int_t Process( THaRunBase& run ) { return Process(&run); }
49 virtual void Print( Option_t* opt="" ) const;
50
51 void EnableBenchmarks( Bool_t b = true );
52 void EnableHelicity( Bool_t b = true );
53 void EnableOtherEvents( Bool_t b = true );
54 void EnableOverwrite( Bool_t b = true );
55 void EnablePhysicsEvents( Bool_t b = true );
56 void EnableRunUpdate( Bool_t b = true );
57 void EnableSlowControl( Bool_t b = true );
58 const char* GetOutFileName() const { return fOutFileName.Data(); }
59 const char* GetCutFileName() const { return fCutFileName.Data(); }
60 const char* GetOdefFileName() const { return fOdefFileName.Data(); }
61 const char* GetSummaryFileName() const { return fSummaryFileName.Data(); }
62 TFile* GetOutFile() const { return fFile; }
64 THaEvent* GetEvent() const { return fEvent; }
65 THaEvData* GetDecoder() const;
66 const std::vector<THaApparatus*>&
67 GetApps() const { return fApps; }
68 const std::vector<THaPhysicsModule*>&
69 GetPhysics() const { return fPhysics; }
72 const std::vector<THaEvtTypeHandler*>&
73 GetEvtHandlers() const { return fEvtHandlers; }
74 const std::vector<THaPostProcess*>&
75 GetPostProcess() const { return fPostProcess; }
78 Bool_t PhysicsEnabled() const { return fDoPhysics; }
81 virtual Int_t SetCountMode( Int_t mode );
82 void SetCrateMapFileName( const char* name );
84 void SetOutFile( const char* name ) { fOutFileName = name; }
85 void SetCutFile( const char* name ) { fCutFileName = name; }
86 void SetOdefFile( const char* name ) { fOdefFileName = name; }
87 void SetSummaryFile( const char* name ) { fSummaryFileName = name; }
88 void SetCompressionLevel( Int_t level ) { fCompress = level; }
89 void SetMarkInterval( UInt_t interval ) { fMarkInterval = interval; }
90 void SetVerbosity( Int_t level ) { fVerbose = level; }
91 void SetCodaVersion(Int_t vers);
92
93 // Set the EPICS event type
94 void SetEpicsEvtType(Int_t itype);
95 void AddEpicsEvtType(Int_t itype);
96
97 static THaAnalyzer* GetInstance() { return fgAnalyzer; }
98
99 // Return codes for analysis routines inside event loop
100 // These should be ordered by severity
102
103 enum {
106 };
107
108 // For SetCountMode
110
111protected:
112 // Test and histogram blocks
113 class Stage_t {
114 public:
115 Stage_t( Int_t _key, Int_t _countkey, const char* _name )
116 : key(_key), countkey(_countkey), name(_name), cut_list(nullptr),
117 hist_list(nullptr), master_cut(nullptr) {}
120 const char* name;
124 };
125 // Statistics counters and message texts
126 enum {
131 };
132 class Counter_t {
133 public:
134 Counter_t( Int_t _key, const char* _description )
135 : key(_key), count(0), description(_description) {}
138 const char* description;
139 };
140
141 TFile* fFile; //The ROOT output file.
142 THaOutput* fOutput; //Flexible ROOT output (tree, histograms)
143 THaEpicsEvtHandler* fEpicsHandler; // EPICS event handler used by THaOutput
144 TString fOutFileName; //Name of output ROOT file.
145 TString fCutFileName; //Name of cut definition file to load
146 TString fLoadedCutFileName;//Name of last loaded cut definition file
147 TString fOdefFileName; //Name of output definition file
148 TString fSummaryFileName; //Name of test/cut statistics output file
149 THaEvent* fEvent; //The event structure to be written to file.
150 Int_t fWantCodaVers; //Version of CODA assumed for file
151 std::vector<Stage_t> fStages; //Parameters for analysis stages
152 std::vector<Counter_t> fCounters;//Statistics counters
153 UInt_t fNev; //Number of events read during most recent replay //FIXME: make ULong64_t
154 UInt_t fMarkInterval; //Interval for printing event numbers
155 Int_t fCompress; //Compression level for ROOT output file
156 Int_t fVerbose; //Verbosity level
157 Int_t fCountMode; //Event counting mode (see ECountMode)
158 THaBenchmark* fBench; //Counters for timing statistics
159 THaEvent* fPrevEvent; //Event structure from last Init()
160 THaRunBase* fRun; //Pointer to current run
161 THaEvData* fEvData; //Instance of decoder used by us
162
163 // Lists of processing modules defined for current analysis
164 std::vector<THaApparatus*> fApps; // Apparatuses
165 std::vector<THaSpectrometer*> fSpectrometers; // Spectrometers
166 std::vector<THaPhysicsModule*> fPhysics; // Physics modules
167 // Modules in the following lists are owned by this THaAnalyzer
168 std::vector<Podd::InterStageModule*> fInterStage; // Inter-stage modules
169 std::vector<THaEvtTypeHandler*> fEvtHandlers; // Event type handlers
170 std::vector<THaPostProcess*> fPostProcess; // Post-processing mods
171 // Combined list of fApps, fInterStage and fPhysics for PhysicsAnalysis.
172 // Does not include fPostProcess and fEvtHandlers.
173 std::vector<THaAnalysisObject*> fAnalysisModules; // Analysis modules
174
175 // Status and control flags
176 Bool_t fIsInit; // Init() called successfully
177 Bool_t fAnalysisStarted; // Process() run and output file open
178 Bool_t fLocalEvent; // fEvent allocated by this object
179 Bool_t fUpdateRun; // Update run parameters during replay
180 Bool_t fOverwrite; // Overwrite existing output files
181 Bool_t fDoBench; // Collect detailed timing statistics
182 Bool_t fDoHelicity; // Enable helicity decoding
183 Bool_t fDoPhysics; // Enable physics event processing
184 Bool_t fDoOtherEvents; // Enable other event processing
185 Bool_t fDoSlowControl; // Enable slow control processing
186
187 // Variables used by analysis functions
188 Bool_t fFirstPhysics; // Status flag for physics analysis
189
190 // Main analysis functions
191 virtual Int_t BeginAnalysis();
192 virtual Int_t DoInit( THaRunBase* run );
193 virtual Int_t EndAnalysis();
194 virtual Int_t MainAnalysis();
195 virtual Int_t PhysicsAnalysis( Int_t code );
196 virtual Int_t SlowControlAnalysis( Int_t code );
197 virtual Int_t OtherAnalysis( Int_t code );
198 virtual Int_t PostProcess( Int_t code );
199 virtual Int_t ReadOneEvent();
200
201 // Support methods & data
202 void ClearCounters();
203 UInt_t GetCount( Int_t which ) const;
204 UInt_t Incr( Int_t which );
205 virtual bool EvalStage( int n );
206 virtual void InitCounters();
207 virtual void InitCuts();
208 virtual void InitStages();
209 virtual Int_t InitModules( const std::vector<THaAnalysisObject*>& module_list,
210 TDatime& run_time );
211 virtual Int_t InitOutput( const std::vector<THaAnalysisObject*>& module_list );
212
214 virtual void PrepareModuleList();
215 virtual void PrintCounters() const;
216 virtual void PrintExitStatus( EExitStatus status ) const;
217 virtual void PrintRunSummary() const;
218 virtual void PrintCutSummary() const;
219 virtual void PrintTimingSummary() const;
220 virtual void PrintSummary( EExitStatus exit_status ) const;
221
222 static THaAnalyzer* fgAnalyzer; //Pointer to instance of this class
223
224 TObject* fExtra; // Additional member data (for binary compat.)
225
226 // In-class constants
227 static const char* const kMasterCutName;
228 static const char* const kDefaultOdefFile;
229
230private:
233
234 ClassDef(THaAnalyzer,0) //Hall A Analyzer Standard Event Loop
235
236};
237
238//---------------- inlines ----------------------------------------------------
239inline UInt_t THaAnalyzer::GetCount( Int_t which ) const
240{
241 return fCounters[which].count;
242}
243
244//_____________________________________________________________________________
246{
247 return ++(fCounters[which].count);
248}
249
250#endif
int Int_t
unsigned int UInt_t
bool Bool_t
const char Option_t
#define ClassDef(name, id)
char name[80]
static const vector< ModuleDef > module_list
Definition THaDetMap.cxx:32
const char * description
Counter_t(Int_t _key, const char *_description)
Stage_t(Int_t _key, Int_t _countkey, const char *_name)
THaEvent * GetEvent() const
Definition THaAnalyzer.h:64
virtual void PrintRunSummary() const
void SetCutFile(const char *name)
Definition THaAnalyzer.h:85
THaBenchmark * fBench
virtual bool EvalStage(int n)
UInt_t GetCount(Int_t which) const
virtual ~THaAnalyzer()
const std::vector< THaApparatus * > & GetApps() const
Definition THaAnalyzer.h:67
void SetOutFile(const char *name)
Definition THaAnalyzer.h:84
virtual void Print(Option_t *opt="") const
static const char *const kDefaultOdefFile
std::vector< THaApparatus * > fApps
virtual void PrintCounters() const
virtual Int_t PhysicsAnalysis(Int_t code)
Bool_t fFirstPhysics
void SetOdefFile(const char *name)
Definition THaAnalyzer.h:86
virtual Int_t SetCountMode(Int_t mode)
const std::vector< THaPhysicsModule * > & GetPhysics() const
Definition THaAnalyzer.h:69
const char * GetCutFileName() const
Definition THaAnalyzer.h:59
void SetEpicsEvtType(Int_t itype)
TString fOutFileName
std::vector< THaAnalysisObject * > fAnalysisModules
static THaAnalyzer * fgAnalyzer
virtual Int_t MainAnalysis()
Int_t GetCompressionLevel() const
Definition THaAnalyzer.h:63
void SetMarkInterval(UInt_t interval)
Definition THaAnalyzer.h:89
virtual Int_t Process(THaRunBase *run=nullptr)
THaEvData * fEvData
Int_t fWantCodaVers
THaEvent * fEvent
TObject * fExtra
Bool_t HasStarted() const
Definition THaAnalyzer.h:76
TString fCutFileName
THaAnalyzer & operator=(const THaAnalyzer &)
UInt_t Incr(Int_t which)
void EnableRunUpdate(Bool_t b=true)
std::vector< THaPhysicsModule * > fPhysics
Bool_t fDoSlowControl
virtual Int_t BeginAnalysis()
virtual Int_t ReadOneEvent()
Bool_t fDoOtherEvents
static THaAnalyzer * GetInstance()
Definition THaAnalyzer.h:97
THaEvent * fPrevEvent
std::vector< THaEvtTypeHandler * > fEvtHandlers
std::vector< Stage_t > fStages
TString fSummaryFileName
void SetEvent(THaEvent *event)
Definition THaAnalyzer.h:83
virtual Int_t InitModules(const std::vector< THaAnalysisObject * > &module_list, TDatime &run_time)
void EnablePhysicsEvents(Bool_t b=true)
Bool_t fIsInit
void SetSummaryFile(const char *name)
Definition THaAnalyzer.h:87
virtual void InitCounters()
TFile * GetOutFile() const
Definition THaAnalyzer.h:62
UInt_t fMarkInterval
void ClearCounters()
virtual void PrepareModuleList()
Bool_t fDoBench
void EnableOverwrite(Bool_t b=true)
const char * GetOdefFileName() const
Definition THaAnalyzer.h:60
TString fOdefFileName
Int_t fCompress
virtual void InitCuts()
std::vector< THaPostProcess * > fPostProcess
void SetCrateMapFileName(const char *name)
virtual void Close()
THaEvData * GetDecoder() const
const std::vector< THaPostProcess * > & GetPostProcess() const
Definition THaAnalyzer.h:75
virtual void PrintExitStatus(EExitStatus status) const
virtual Int_t Init(THaRunBase *run)
virtual void PrintSummary(EExitStatus exit_status) const
THaAnalyzer(const THaAnalyzer &)
const std::vector< THaEvtTypeHandler * > & GetEvtHandlers() const
Definition THaAnalyzer.h:73
const char * GetSummaryFileName() const
Definition THaAnalyzer.h:61
Bool_t SlowControlEnabled() const
Definition THaAnalyzer.h:80
THaOutput * fOutput
void SetCodaVersion(Int_t vers)
const char * GetOutFileName() const
Definition THaAnalyzer.h:58
THaRunBase * fRun
void EnableBenchmarks(Bool_t b=true)
virtual void PrintCutSummary() const
void EnableOtherEvents(Bool_t b=true)
virtual Int_t DoInit(THaRunBase *run)
TFile * fFile
std::vector< Podd::InterStageModule * > fInterStage
static const char *const kMasterCutName
Bool_t fOverwrite
virtual Int_t SlowControlAnalysis(Int_t code)
THaEpicsEvtHandler * fEpicsHandler
Bool_t OtherEventsEnabled() const
Definition THaAnalyzer.h:79
virtual void PrintTimingSummary() const
Bool_t fLocalEvent
std::vector< THaSpectrometer * > fSpectrometers
Bool_t fDoPhysics
Bool_t PhysicsEnabled() const
Definition THaAnalyzer.h:78
Bool_t fUpdateRun
virtual Int_t PostProcess(Int_t code)
Int_t fCountMode
virtual Int_t InitOutput(const std::vector< THaAnalysisObject * > &module_list)
void EnableHelicity(Bool_t b=true)
std::vector< Counter_t > fCounters
Int_t Process(THaRunBase &run)
Definition THaAnalyzer.h:48
virtual Int_t AddInterStage(Podd::InterStageModule *module)
Int_t Init(THaRunBase &run)
Definition THaAnalyzer.h:46
void SetVerbosity(Int_t level)
Definition THaAnalyzer.h:90
virtual Int_t AddPostProcess(THaPostProcess *module)
virtual void InitStages()
virtual Int_t EndAnalysis()
TString fLoadedCutFileName
void EnableSlowControl(Bool_t b=true)
virtual Int_t OtherAnalysis(Int_t code)
void SetCompressionLevel(Int_t level)
Definition THaAnalyzer.h:88
Bool_t HelicityEnabled() const
Definition THaAnalyzer.h:77
void AddEpicsEvtType(Int_t itype)
Bool_t fAnalysisStarted
THaEpicsEvtHandler * GetEpicsEvtHandler() const
Definition THaAnalyzer.h:71
Bool_t fDoHelicity
const char * Data() const