BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSOutputROOTEventHistograms.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 BDSOUTPUTROOTEVENTHISTOGRAMS_H
20#define BDSOUTPUTROOTEVENTHISTOGRAMS_H
21
22#include <string>
23#include <vector>
24
25#ifdef USE_BOOST
26#include <boost/histogram.hpp>
27#endif
28
29#ifndef __ROOTBUILD__
30#include "globals.hh"
31#endif
32
33#include "Rtypes.h"
34#include "TObject.h"
35#include "BDSBH4DBase.hh"
36
37class TH1D;
38class TH2D;
39class TH3D;
40
47class BDSOutputROOTEventHistograms: public TObject
48{
49public:
53 BDSOutputROOTEventHistograms(std::vector<TH1D*>& histogram1DIn,
54 std::vector<TH2D*>& histogram2DIn,
55 std::vector<TH3D*>& histogram3DIn,
56 std::vector<BDSBH4DBase*>& histograms4DIn);
58
60 int Create1DHistogramSTD(std::string name, std::string title,
61 int nbins, double xmin, double xmax);
62
63#ifndef __ROOTBUILD__
64 G4int Create1DHistogram(G4String name, G4String title,
65 G4int nbins, G4double xmin, G4double xmax);
66 G4int Create1DHistogram(G4String name, G4String title,
67 std::vector<double>& edges);
68 G4int Create2DHistogram(G4String name, G4String title,
69 G4int nxbins, G4double xmin, G4double xmax,
70 G4int nybins, G4double ymin, G4double ymax);
71 G4int Create2DHistogram(G4String name, G4String title,
72 std::vector<double>& xedges,
73 std::vector<double>& yedges);
74 G4int Create3DHistogram(G4String name, G4String title,
75 G4int nxbins, G4double xmin, G4double xmax,
76 G4int nybins, G4double ymin, G4double ymax,
77 G4int nzbins, G4double zmin, G4double zmax);
78 G4int Create3DHistogram(G4String name, G4String title,
79 std::vector<double>& xedges,
80 std::vector<double>& yedges,
81 std::vector<double>& zedges);
82 G4int Create4DHistogram(const G4String& name,
83 const G4String& title,
84 const G4String& eScale,
85 const std::vector<double>& eBinsEdges,
86 unsigned int nxbins, G4double xmin, G4double xmax,
87 unsigned int nybins, G4double ymin, G4double ymax,
88 unsigned int nzbins, G4double zmin, G4double zmax,
89 unsigned int nebins, G4double emin, G4double emax);
90
91 void Fill1DHistogram(G4int histoId, G4double value, G4double weight = 1.0);
92 void Fill2DHistogram(G4int histoId, G4double xValue, G4double yValue, G4double weight = 1.0);
93 void Fill3DHistogram(G4int histoId, G4double xValue, G4double yValue, G4double zValue, G4double weight = 1.0);
94 void Fill4DHistogram(G4int histoId, G4double xValue, G4double yValue, G4double zvalue, G4double eValue);
95
98 void Set3DHistogramBinContent(G4int histoId,
99 G4int globalBinID,
100 G4double value);
101
102 void Set4DHistogramBinContent(G4int histoId,
103 G4int x,
104 G4int y,
105 G4int z,
106 G4int e,
107 G4double value);
108
110 void AccumulateHistogram3D(G4int histoId,
111 TH3D* otherHistogram);
112 void AccumulateHistogram4D(G4int histoId,
113 BDSBH4DBase* otherHistogram);
114#endif
116 virtual void Flush();
117
119 void Fill(const BDSOutputROOTEventHistograms* rhs);
120
123
125 std::vector<TH1D*>& Get1DHistograms() {return histograms1D;}
126 std::vector<TH2D*>& Get2DHistograms() {return histograms2D;}
127 std::vector<TH3D*>& Get3DHistograms() {return histograms3D;}
128 std::vector<BDSBH4DBase*>& Get4DHistograms() {return histograms4D;}
129 TH1D* Get1DHistogram(int iHisto) const {return histograms1D[iHisto];}
130 TH2D* Get2DHistogram(int iHisto) const {return histograms2D[iHisto];}
131 TH3D* Get3DHistogram(int iHisto) const {return histograms3D[iHisto];}
132 BDSBH4DBase* Get4DHistogram(int iHisto) const {return histograms4D[iHisto];}
134
135private:
136 std::vector<TH1D*> histograms1D;
137 std::vector<TH2D*> histograms2D;
138 std::vector<TH3D*> histograms3D;
139 std::vector<BDSBH4DBase*> histograms4D;
140
142};
143
144#endif
Base class for the 4D histogram classes.
Definition: BDSBH4DBase.hh:37
Holder for a set of histograms to be stored.
std::vector< BDSBH4DBase * > & Get4DHistograms()
Accessors.
std::vector< TH2D * > & Get2DHistograms()
Accessors.
std::vector< TH3D * > & Get3DHistograms()
Accessors.
TH2D * Get2DHistogram(int iHisto) const
Accessors.
virtual void Flush()
Flush the contents.
std::vector< TH1D * > & Get1DHistograms()
Accessors.
void AccumulateHistogram3D(G4int histoId, TH3D *otherHistogram)
Add the values from one supplied 3D histogram to another. Uses TH3-Add().
BDSBH4DBase * Get4DHistogram(int iHisto) const
Accessors.
TH1D * Get1DHistogram(int iHisto) const
Accessors.
int Create1DHistogramSTD(std::string name, std::string title, int nbins, double xmin, double xmax)
Interface function to create a 1D histogram using only standard types.
void Fill(const BDSOutputROOTEventHistograms *rhs)
Copy (using the TH->Clone) method from another instance.
void FillSimple(const BDSOutputROOTEventHistograms *rhs)
Copy (without using the TH->Clone) method from another instance. (Quicker).
void Set3DHistogramBinContent(G4int histoId, G4int globalBinID, G4double value)
TH3D * Get3DHistogram(int iHisto) const
Accessors.