BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSOutputROOTEventSampler.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 BDSOUTPUTROOTEVENTSAMPLER_H
20#define BDSOUTPUTROOTEVENTSAMPLER_H
21
22#include "Rtypes.h"
23#include "TObject.h"
24
25#include <string>
26#include <vector>
27
29class TTree;
30
31#ifndef __ROOTBUILD__
32#include "globals.hh"
34class BDSHitSampler;
36#endif
37
44template<class U> class BDSOutputROOTEventSampler: public TObject
45{
46public:
47 std::string samplerName; //|| Don't split the header
48
49 int n;
50 std::vector<U> energy;
51 std::vector<U> x;
52 std::vector<U> y;
53 U z;
54 std::vector<U> xp;
55 std::vector<U> yp;
56 std::vector<U> zp;
57 std::vector<U> p;
58 std::vector<U> T;
59
60 std::vector<U> weight;
61 std::vector<int> partID;
62 std::vector<int> parentID;
63 std::vector<int> trackID;
64 int modelID;
65 std::vector<int> turnNumber;
66
67 U S;
68
69 std::vector<U> r;
70 std::vector<U> rp;
71 std::vector<U> phi;
72 std::vector<U> phip;
73 std::vector<U> theta;
74
76 std::vector<int> charge;
77 std::vector<U> kineticEnergy;
78 std::vector<U> mass;
79 std::vector<U> rigidity;
80 std::vector<bool> isIon;
81 std::vector<int> ionA;
82 std::vector<int> ionZ;
83 std::vector<int> nElectrons;
85
87 std::vector<U> getKineticEnergy();
88 std::vector<U> getMass();
89 std::vector<U> getRigidity();
90 std::vector<bool> getIsIon();
91 std::vector<int> getIonA();
92 std::vector<int> getIonZ();
94
96 explicit BDSOutputROOTEventSampler(std::string samplerNameIn);
98#ifndef __ROOTBUILD__
99 void Fill(const BDSHitSampler* hit,
100 G4bool storeMass = false,
101 G4bool storeCharge = false,
102 G4bool storePolarCoords = false,
103 G4bool storeElectrons = false,
104 G4bool storeRigidity = false,
105 G4bool storeKineticEnergy = false);
109 void Fill(const BDSParticleCoordsFull& coords,
110 G4double momentumIn,
111 G4double charge,
112 G4int pdgID,
113 G4int turnsTaken,
114 G4int beamlineIndex,
115 G4int nElectronsIn,
116 G4double massIn,
117 G4double rigidityIn,
118 G4bool fillIon = true,
119 G4bool* isIon = nullptr,
120 G4int* ionA = nullptr,
121 G4int* ionZ = nullptr);
122 void FillPolarCoords(const BDSParticleCoordsFull& coords);
123 void Fill(const BDSPrimaryVertexInformationV* vertexInfos,
124 const G4int turnsTaken);
125#endif
126 void Fill(const BDSOutputROOTEventSampler<U>* other);
127
129 inline void FillIon() {isIon = getIsIon(); ionA = getIonA(); ionZ = getIonZ();}
131
132 void SetBranchAddress(TTree *);
133 virtual void Flush();
134
135 static BDSOutputROOTParticleData* particleTable;
136
137 ClassDef(BDSOutputROOTEventSampler,5);
138};
139
140// This unusually has to be in the header because it's a templated static member, so we need
141// the compiler to generate the intialisation for the ones required.
143
144#endif
The information recorded from a particle impacting a sampler.
Information stored per sampler per event.
std::vector< int > getIonA()
Function to calculate on the fly the parameters.
std::vector< bool > isIon
These are not filled by default.
std::vector< int > getIonZ()
Function to calculate on the fly the parameters.
std::vector< int > nElectrons
These are not filled by default.
std::vector< int > ionA
These are not filled by default.
std::vector< U > getRigidity()
Function to calculate on the fly the parameters.
std::vector< U > getMass()
Function to calculate on the fly the parameters.
std::vector< U > rigidity
These are not filled by default.
std::vector< bool > getIsIon()
Function to calculate on the fly the parameters.
std::vector< int > charge
These are not filled by default.
void FillIon()
Calculate and fill calculated variables.
virtual void Flush()
Clean Sampler.
std::vector< U > kineticEnergy
These are not filled by default.
void FillPolarCoords(const BDSParticleCoordsFull &coords)
Calculate polar coords and fill.
std::vector< U > mass
These are not filled by default.
std::vector< int > ionZ
These are not filled by default.
std::vector< U > getKineticEnergy()
Function to calculate on the fly the parameters.
Geant4 particle data for particles used in simulation.
A set of particle coordinates including energy and weight.
Full set of coordinates for association with primary vertex. Vector version.