BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
HistogramAccumulator.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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 HISTOGRAMACCUMULATOR_H
20#define HISTOGRAMACCUMULATOR_H
21
22#include <string>
23
24#include "Rtypes.h" // for classdef
25
26class TH1;
27
58{
59public:
62
66 HistogramAccumulator(TH1* baseHistogram,
67 const std::string& resultHistNameIn,
68 const std::string& resultHistTitleIn);
69
73 HistogramAccumulator(TH1* baseHistogram,
74 int nDimensionsIn,
75 const std::string& resultHistName,
76 const std::string& resultHistTitle);
77
80 virtual ~HistogramAccumulator();
81
85 virtual void Accumulate(TH1* newValue);
86
89 virtual TH1* Terminate();
90
92 inline TH1* Result() const {return result;}
93
96 inline void AddNEmptyEntries(unsigned long i){n += i;}
97
99 inline unsigned long N() const {return n;}
100
101protected:
107 virtual void AccumulateSingleValue(double oldMean,
108 double oldVari,
109 double x,
110 double xVari,
111 unsigned long nEntriesAccumulated,
112 unsigned long nEntriesToAccumulate,
113 double& newMean,
114 double& newVari) const;
115
117 unsigned long n;
119 const std::string resultHistName;
120 const std::string resultHistTitle;
121 TH1* mean;
122 TH1* variance;
123 TH1* result;
124
125 ClassDef(HistogramAccumulator,1);
126};
127
128#endif
Class to accumulate and merge histograms in different ways.
const std::string resultHistName
Name for resultant histogram.
unsigned long N() const
Access currently accumulated number of entries.
int nDimensions
Number of dimensions.
virtual void AccumulateSingleValue(double oldMean, double oldVari, double x, double xVari, unsigned long nEntriesAccumulated, unsigned long nEntriesToAccumulate, double &newMean, double &newVari) const
bool terminated
Whether this instance has been finished.
const std::string resultHistTitle
Title for resultant histogram.
unsigned long n
Counter.
void AddNEmptyEntries(unsigned long i)
HistogramAccumulator()
Default constructor only for ROOT reflexivity - not intended for use.
virtual void Accumulate(TH1 *newValue)
TH1 * Result() const
Accessor.