BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Analysis.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef ANALYSIS_H
20#define ANALYSIS_H
21
22#include "TH1D.h"
23#include "TH2D.h"
24#include "TH3D.h"
25#include "BDSBH4DBase.hh"
26
27#include <map>
28#include <string>
29#include <vector>
30
31class HistogramDef;
34class TChain;
35class TFile;
36
44{
45public:
50 Analysis(const std::string& treeNameIn,
51 TChain* chainIn,
52 const std::string& mergedHistogramNameIn,
53 bool perEntryAnalysis = true,
54 bool debugIn = false);
55 virtual ~Analysis();
56
58 virtual void Execute();
59
63 virtual void Process() = 0;
64
66 virtual void UserProcess();
67
69 virtual void SimpleHistograms();
70
73
75 void AccumulatePerEntryHistograms(long int entryNumber);
76
79
82 virtual void Terminate();
83
85 virtual void Write(TFile* outputFile);
86
87protected:
89 void FillHistogram(HistogramDef* definition,
90 std::vector<TH1*>* outputHistograms = nullptr);
91
92 std::string treeName;
93 TChain* chain;
94 std::string mergedHistogramName;
95 std::vector<TH1*> simpleHistograms;
96 std::vector<PerEntryHistogram*> perEntryHistograms;
98 bool debug;
99 long int entries;
100 bool perEntry;
101
102private:
104 Analysis() = delete;
105};
106
107#endif
Base class for any TTree analysis.
Definition: Analysis.hh:44
virtual void Execute()
Method which calls all other methods in order.
Definition: Analysis.cc:64
virtual void Terminate()
Definition: Analysis.cc:123
virtual void SimpleHistograms()
Process histogram definitions from configuration instance.
Definition: Analysis.cc:86
void PreparePerEntryHistograms()
Create structures necessary for per entry histograms.
Definition: Analysis.cc:100
long int entries
Number of entries in the chain.
Definition: Analysis.hh:99
bool debug
Whether debug print out is used or not.
Definition: Analysis.hh:98
HistogramMeanFromFile * histoSum
Merge of per event stored histograms.
Definition: Analysis.hh:97
void AccumulatePerEntryHistograms(long int entryNumber)
Accumulate means and variances for per entry histograms.
Definition: Analysis.cc:111
virtual void Write(TFile *outputFile)
Write rebdsim histograms.
Definition: Analysis.cc:131
Analysis()=delete
No default constructor for this base class.
void TerminatePerEntryHistograms()
Prepare result of per entry histogram accumulation.
Definition: Analysis.cc:117
std::string mergedHistogramName
Name of directory for merged histograms.
Definition: Analysis.hh:94
virtual void UserProcess()
Virtual function for user to overload and use. Does nothing by default.
Definition: Analysis.cc:83
bool perEntry
Whether to analyse each entry in the tree in a for loop or not.
Definition: Analysis.hh:100
void FillHistogram(HistogramDef *definition, std::vector< TH1 * > *outputHistograms=nullptr)
Create an individual histogram based on a definition.
Definition: Analysis.cc:166
virtual void Process()=0
Common specification for a histogram.
Definition: HistogramDef.hh:33
Accumulator to merge pre-made per-entry histograms.
Holder for information to calculate per entry histograms.