1 #ifndef _DAnalysisAction_
2 #define _DAnalysisAction_
8 #include "TDirectoryFile.h"
13 #include "JANA/JEventLoop.h"
21 using namespace DAnalysis;
30 DAnalysisAction(
const DReaction* locReaction,
string locActionBaseName,
bool locUseKinFitResultsFlag =
false,
string locActionUniqueString =
"");
42 virtual void Initialize(JEventLoop* locEventLoop) = 0;
49 bool operator()(JEventLoop* locEventLoop);
50 bool operator()(JEventLoop* locEventLoop,
const DParticleCombo* locParticleCombo);
57 virtual bool Perform_Action(JEventLoop* locEventLoop,
const DParticleCombo* locParticleCombo = NULL) = 0;
64 TDirectoryFile* CreateAndChangeTo_ActionDirectory(
void);
65 TDirectoryFile* ChangeTo_BaseDirectory(
void);
66 TDirectoryFile* CreateAndChangeTo_Directory(TDirectoryFile* locBaseDirectory,
string locDirName,
string locDirTitle);
67 TDirectoryFile* CreateAndChangeTo_Directory(
string locDirName,
string locDirTitle);
68 TDirectoryFile* CreateAndChangeTo_Directory(
string locBaseDirectoryPath,
string locDirName,
string locDirTitle);
72 template <
typename DHistType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax)
const;
73 template <
typename DHistType,
typename DBinType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges)
const;
75 template <
typename DHistType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax, Int_t locNumBinsY, Double_t locYRangeMin, Double_t locYRangeMax)
const;
76 template <
typename DHistType,
typename DBinType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges, Int_t locNumBinsY, DBinType* locYBinEdges)
const;
77 template <
typename DHistType,
typename DBinType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges, Int_t locNumBinsY, Double_t locYRangeMin, Double_t locYRangeMax)
const;
78 template <
typename DHistType,
typename DBinType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax, Int_t locNumBinsY, DBinType* locYBinEdges)
const;
80 template <
typename DHistType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax, Int_t locNumBinsY, Double_t locYRangeMin, Double_t locYRangeMax, Int_t locNumBinsZ, Double_t locZRangeMin, Double_t locZRangeMax)
const;
81 template <
typename DHistType,
typename DBinType> DHistType* GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges, Int_t locNumBinsY, DBinType* locYBinEdges, Int_t locNumBinsZ, DBinType* locZBinEdges)
const;
89 void Lock_Action(
void);
90 void Unlock_Action(
void);
106 bool dCalledPriorWithComboFlag =
false;
108 template <
typename DHistType>
bool Check_IsValidTH1(
string locHistName)
const;
109 template <
typename DHistType>
bool Check_IsValidTH2(
string locHistName)
const;
110 template <
typename DHistType>
bool Check_IsValidTH3(
string locHistName)
const;
121 inline bool DAnalysisAction::operator()(JEventLoop* locEventLoop)
124 bool locResult = Perform_Action(locEventLoop,
nullptr);
125 return (dPerformAntiCut ? !locResult : locResult);
128 inline bool DAnalysisAction::operator()(JEventLoop* locEventLoop,
const DParticleCombo* locParticleCombo)
131 bool locResult = Perform_Action(locEventLoop, locParticleCombo);
132 dCalledPriorWithComboFlag =
true;
133 return (dPerformAntiCut ? !locResult : locResult);
136 inline void DAnalysisAction::Lock_Action(
void)
138 pthread_rwlock_wrlock(dActionLock);
141 inline void DAnalysisAction::Unlock_Action(
void)
143 pthread_rwlock_unlock(dActionLock);
146 inline TDirectoryFile* DAnalysisAction::ChangeTo_BaseDirectory(
void)
149 TFile* locFile = (TFile*)gROOT->FindObject(dOutputFileName.c_str());
154 return (TDirectoryFile*)gDirectory;
157 inline TDirectoryFile* DAnalysisAction::CreateAndChangeTo_Directory(TDirectoryFile* locBaseDirectory,
string locDirName,
string locDirTitle)
160 locBaseDirectory->cd();
161 TDirectoryFile* locSubDirectory =
static_cast<TDirectoryFile*
>(locBaseDirectory->Get(locDirName.c_str()));
162 if(locSubDirectory == NULL)
163 locSubDirectory =
new TDirectoryFile(locDirName.c_str(), locDirTitle.c_str());
164 locSubDirectory->cd();
165 return locSubDirectory;
168 inline TDirectoryFile* DAnalysisAction::CreateAndChangeTo_Directory(
string locDirName,
string locDirTitle)
171 return CreateAndChangeTo_Directory((TDirectoryFile*)gDirectory, locDirName, locDirTitle);
174 inline TDirectoryFile* DAnalysisAction::CreateAndChangeTo_Directory(
string locBaseDirectoryPath,
string locDirName,
string locDirTitle)
177 TDirectoryFile* locBaseDirectory =
static_cast<TDirectoryFile*
>(gDirectory->GetDirectory(locBaseDirectoryPath.c_str()));
178 return CreateAndChangeTo_Directory(locBaseDirectory, locDirName, locDirTitle);
181 template <
typename DHistType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax)
const
184 if(!Check_IsValidTH1<DHistType>(locHistName))
187 const char* locHistNameCString = locHistName.c_str();
188 const char* locHistTitleCString = locHistTitle.c_str();
189 TObject*
locHist = gDirectory->Get(locHistNameCString);
191 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXRangeMin, locXRangeMax);
193 return static_cast<DHistType*
>(
locHist);
196 template <
typename DHistType,
typename DBinType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges)
const
199 if(!Check_IsValidTH1<DHistType>(locHistName))
202 const char* locHistNameCString = locHistName.c_str();
203 const char* locHistTitleCString = locHistTitle.c_str();
204 TObject*
locHist = gDirectory->Get(locHistNameCString);
206 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXBinEdges);
208 return static_cast<DHistType*
>(
locHist);
211 template <
typename DHistType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax, Int_t locNumBinsY, Double_t locYRangeMin, Double_t locYRangeMax)
const
214 if(!Check_IsValidTH2<DHistType>(locHistName))
217 const char* locHistNameCString = locHistName.c_str();
218 const char* locHistTitleCString = locHistTitle.c_str();
219 TObject*
locHist = gDirectory->Get(locHistNameCString);
221 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXRangeMin, locXRangeMax, locNumBinsY, locYRangeMin, locYRangeMax);
223 return static_cast<DHistType*
>(
locHist);
226 template <
typename DHistType,
typename DBinType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges, Int_t locNumBinsY, DBinType* locYBinEdges)
const
229 if(!Check_IsValidTH2<DHistType>(locHistName))
232 const char* locHistNameCString = locHistName.c_str();
233 const char* locHistTitleCString = locHistTitle.c_str();
234 TObject*
locHist = gDirectory->Get(locHistNameCString);
236 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXBinEdges, locNumBinsY, locYBinEdges);
238 return static_cast<DHistType*
>(
locHist);
241 template <
typename DHistType,
typename DBinType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges, Int_t locNumBinsY, Double_t locYRangeMin, Double_t locYRangeMax)
const
244 if(!Check_IsValidTH2<DHistType>(locHistName))
247 const char* locHistNameCString = locHistName.c_str();
248 const char* locHistTitleCString = locHistTitle.c_str();
249 TObject*
locHist = gDirectory->Get(locHistNameCString);
251 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXBinEdges, locNumBinsY, locYRangeMin, locYRangeMax);
253 return static_cast<DHistType*
>(
locHist);
256 template <
typename DHistType,
typename DBinType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax, Int_t locNumBinsY, DBinType* locYBinEdges)
const
259 if(!Check_IsValidTH2<DHistType>(locHistName))
262 const char* locHistNameCString = locHistName.c_str();
263 const char* locHistTitleCString = locHistTitle.c_str();
264 TObject*
locHist = gDirectory->Get(locHistNameCString);
266 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXRangeMin, locXRangeMax, locNumBinsY, locYBinEdges);
268 return static_cast<DHistType*
>(
locHist);
271 template <
typename DHistType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, Double_t locXRangeMin, Double_t locXRangeMax, Int_t locNumBinsY, Double_t locYRangeMin, Double_t locYRangeMax, Int_t locNumBinsZ, Double_t locZRangeMin, Double_t locZRangeMax)
const
274 if(!Check_IsValidTH3<DHistType>(locHistName))
277 const char* locHistNameCString = locHistName.c_str();
278 const char* locHistTitleCString = locHistTitle.c_str();
279 TObject*
locHist = gDirectory->Get(locHistNameCString);
281 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXRangeMin, locXRangeMax, locNumBinsY, locYRangeMin, locYRangeMax, locNumBinsZ, locZRangeMin, locZRangeMax);
283 return static_cast<DHistType*
>(
locHist);
286 template <
typename DHistType,
typename DBinType>
inline DHistType* DAnalysisAction::GetOrCreate_Histogram(
string locHistName,
string locHistTitle, Int_t locNumBinsX, DBinType* locXBinEdges, Int_t locNumBinsY, DBinType* locYBinEdges, Int_t locNumBinsZ, DBinType* locZBinEdges)
const
289 if(!Check_IsValidTH3<DHistType>(locHistName))
292 const char* locHistNameCString = locHistName.c_str();
293 const char* locHistTitleCString = locHistTitle.c_str();
294 TObject*
locHist = gDirectory->Get(locHistNameCString);
296 return new DHistType(locHistNameCString, locHistTitleCString, locNumBinsX, locXBinEdges, locNumBinsY, locYBinEdges, locNumBinsZ, locZBinEdges);
298 return static_cast<DHistType*
>(
locHist);
301 template <
typename DHistType>
inline bool DAnalysisAction::Check_IsValidTH3(
string locHistName)
const
303 const char* locTypeName = DHistType::Class()->GetName();
304 const char* locBadTypeName =
"TH3";
305 if((!DHistType::Class()->InheritsFrom(
"TH3")) || (
string(locTypeName) ==
string(locBadTypeName)))
307 cout <<
"ERROR, WRONG CLASS TYPE IN GetOrCreate_Histogram for HISTOGRAM " << locHistName <<
". HISTOGRAM NOT CREATED." << endl;
313 template <
typename DHistType>
inline bool DAnalysisAction::Check_IsValidTH2(
string locHistName)
const
315 const char* locTypeName = DHistType::Class()->GetName();
316 const char* locBadTypeName =
"TH2";
317 if((!DHistType::Class()->InheritsFrom(
"TH2")) || (
string(locTypeName) ==
string(locBadTypeName)) || DHistType::Class()->InheritsFrom(
"TH3"))
319 cout <<
"ERROR, WRONG CLASS TYPE IN GetOrCreate_Histogram for HISTOGRAM " << locHistName <<
". HISTOGRAM NOT CREATED." << endl;
325 template <
typename DHistType>
inline bool DAnalysisAction::Check_IsValidTH1(
string locHistName)
const
327 const char* locTypeName = DHistType::Class()->GetName();
328 const char* locBadTypeName =
"TH1";
329 if((!DHistType::Class()->InheritsFrom(
"TH1")) || (
string(locTypeName) ==
string(locBadTypeName)) || DHistType::Class()->InheritsFrom(
"TH2") || DHistType::Class()->InheritsFrom(
"TH3"))
331 cout <<
"ERROR, WRONG CLASS TYPE IN GetOrCreate_Histogram for HISTOGRAM " << locHistName <<
". HISTOGRAM NOT CREATED." << endl;
339 #endif // _DAnalysisAction_
virtual ~DAnalysisAction(void)
const DReaction * Get_Reaction(void) const
const DReaction * dReaction
bool Get_UseKinFitResultsFlag(void) const
bool dUseKinFitResultsFlag
pthread_rwlock_t * dActionLock
string Get_ActionUniqueString(void) const
bool Get_CalledPriorWithComboFlag(void) const
virtual void Reset_NewEvent(void)
string dActionUniqueString
virtual string Get_ActionName(void) const