BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Base class for any TTree analysis. More...
#include <Analysis.hh>
Public Member Functions | |
Analysis (const std::string &treeNameIn, TChain *chainIn, const std::string &mergedHistogramNameIn, bool perEntryAnalysis=true, bool debugIn=false) | |
virtual void | Execute () |
Method which calls all other methods in order. More... | |
virtual void | Process ()=0 |
virtual void | UserProcess () |
Virtual function for user to overload and use. Does nothing by default. More... | |
virtual void | SimpleHistograms () |
Process histogram definitions from configuration instance. More... | |
void | PreparePerEntryHistograms () |
Create structures necessary for per entry histograms. More... | |
void | AccumulatePerEntryHistograms (long int entryNumber) |
Accumulate means and variances for per entry histograms. More... | |
void | TerminatePerEntryHistograms () |
Prepare result of per entry histogram accumulation. More... | |
virtual void | Terminate () |
virtual void | Write (TFile *outputFile) |
Write rebdsim histograms. More... | |
Protected Member Functions | |
void | FillHistogram (HistogramDef *definition, std::vector< TH1 * > *outputHistograms=nullptr) |
Create an individual histogram based on a definition. More... | |
Protected Attributes | |
std::string | treeName |
TChain * | chain |
std::string | mergedHistogramName |
Name of directory for merged histograms. More... | |
std::vector< TH1 * > | simpleHistograms |
std::vector< PerEntryHistogram * > | perEntryHistograms |
HistogramMeanFromFile * | histoSum |
Merge of per event stored histograms. More... | |
bool | debug |
Whether debug print out is used or not. More... | |
long int | entries |
Number of entries in the chain. More... | |
bool | perEntry |
Whether to analyse each entry in the tree in a for loop or not. More... | |
Private Member Functions | |
Analysis ()=delete | |
No default constructor for this base class. | |
Base class for any TTree analysis.
Definition at line 43 of file Analysis.hh.
Analysis::Analysis | ( | const std::string & | treeNameIn, |
TChain * | chainIn, | ||
const std::string & | mergedHistogramNameIn, | ||
bool | perEntryAnalysis = true , |
||
bool | debugIn = false |
||
) |
Requires the name of the tree to analyse, the chain of files to operate on, the name of the directory to create in the output root file for the combined existing histograms from that tree, whether to operate on each entry in the tree and whether or not we're in debug mode.
Definition at line 43 of file Analysis.cc.
|
virtual |
Definition at line 57 of file Analysis.cc.
void Analysis::AccumulatePerEntryHistograms | ( | long int | entryNumber | ) |
Accumulate means and variances for per entry histograms.
Definition at line 111 of file Analysis.cc.
Referenced by BeamAnalysis::Process(), EventAnalysis::Process(), ModelAnalysis::Process(), OptionsAnalysis::Process(), and RunAnalysis::Process().
|
virtual |
Method which calls all other methods in order.
Reimplemented in EventAnalysis.
Definition at line 64 of file Analysis.cc.
References perEntry, PreparePerEntryHistograms(), Process(), SimpleHistograms(), and Terminate().
|
protected |
Create an individual histogram based on a definition.
Definition at line 166 of file Analysis.cc.
References HistogramFactory::CreateHistogram().
Referenced by EventAnalysis::FillHistogram(), and SimpleHistograms().
void Analysis::PreparePerEntryHistograms | ( | ) |
Create structures necessary for per entry histograms.
Definition at line 100 of file Analysis.cc.
References Config::Instance().
Referenced by Execute(), and EventAnalysis::Execute().
|
pure virtual |
Operate on each entry in the tree. Pure virtual as it is not known what analysis will be formed in any derived class. This is only called if perEntry is true (by default it is).
Implemented in BeamAnalysis, EventAnalysis, ModelAnalysis, OptionsAnalysis, and RunAnalysis.
Referenced by Execute().
|
virtual |
Process histogram definitions from configuration instance.
Reimplemented in EventAnalysis.
Definition at line 86 of file Analysis.cc.
References FillHistogram(), Config::HistogramDefinitionsSimple(), and Config::Instance().
Referenced by Execute(), and EventAnalysis::SimpleHistograms().
|
virtual |
Optional final action after Process() and SimpleHistograms(). The version in this base class terminates the histogram merges if there are any in histoSum.
Reimplemented in EventAnalysis.
Definition at line 123 of file Analysis.cc.
References histoSum, perEntry, HistogramMeanFromFile::Terminate(), and TerminatePerEntryHistograms().
Referenced by Execute(), and EventAnalysis::Terminate().
void Analysis::TerminatePerEntryHistograms | ( | ) |
Prepare result of per entry histogram accumulation.
Definition at line 117 of file Analysis.cc.
Referenced by Terminate().
|
virtual |
Virtual function for user to overload and use. Does nothing by default.
Definition at line 83 of file Analysis.cc.
Referenced by EventAnalysis::Process(), and RunAnalysis::Process().
|
virtual |
Write rebdsim histograms.
Reimplemented in EventAnalysis.
Definition at line 131 of file Analysis.cc.
References histoSum, and HistogramMeanFromFile::Write().
Referenced by EventAnalysis::Write().
|
protected |
Definition at line 93 of file Analysis.hh.
|
protected |
Whether debug print out is used or not.
Definition at line 98 of file Analysis.hh.
Referenced by EventAnalysis::EventAnalysis(), EventAnalysis::Process(), and RunAnalysis::Process().
|
protected |
Number of entries in the chain.
Definition at line 99 of file Analysis.hh.
Referenced by EventAnalysisOrbit::ExtractOrbit(), BeamAnalysis::Process(), EventAnalysis::Process(), ModelAnalysis::Process(), OptionsAnalysis::Process(), and EventAnalysis::SetPrintModuloFraction().
|
protected |
Merge of per event stored histograms.
Definition at line 97 of file Analysis.hh.
Referenced by EventAnalysis::Process(), RunAnalysis::Process(), Terminate(), and Write().
|
protected |
Name of directory for merged histograms.
Definition at line 94 of file Analysis.hh.
|
protected |
Whether to analyse each entry in the tree in a for loop or not.
Definition at line 100 of file Analysis.hh.
Referenced by Execute(), EventAnalysis::Execute(), and Terminate().
|
protected |
Definition at line 96 of file Analysis.hh.
|
protected |
Definition at line 95 of file Analysis.hh.
|
protected |
Definition at line 92 of file Analysis.hh.