BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSOutputROOTEventLoss.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 BDSOUTPUTROOTEVENTLOSS_H
20#define BDSOUTPUTROOTEVENTLOSS_H
21
22#ifndef __ROOTBUILD__
26#endif
27
28#include "TObject.h"
29
30#include <vector>
31
38class BDSOutputROOTEventLoss: public TObject
39{
40public:
43
44 int n;
45 std::vector<float> energy;
46 std::vector<float> S;
47 std::vector<float> weight;
48
49 std::vector<int> partID;
50 std::vector<int> trackID;
51 std::vector<int> parentID;
52
53 std::vector<int> modelID;
54 std::vector<int> turn;
55
57 std::vector<float> x;
58 std::vector<float> y;
59 std::vector<float> z;
61
63 std::vector<float> X;
64 std::vector<float> Y;
65 std::vector<float> Z;
67
68 std::vector<float> T;
69
70 std::vector<float> stepLength;
71 std::vector<float> preStepKineticEnergy;
72
73 std::vector<int> postStepProcessType;
74 std::vector<int> postStepProcessSubType;
75
77 void Fill(const BDSOutputROOTEventLoss* other);
78 virtual void Flush();
79
80#ifndef __ROOTBUILD__
81 BDSOutputROOTEventLoss(bool storeTurnIn,
82 bool storeLinksIn,
83 bool storeModleIDIn,
84 bool storeLocalIn,
85 bool storeGobalIn,
86 bool storeTimeIn,
87 bool storeStepLengthIn,
88 bool storePreStepKineticEnergyIn,
89 bool storePhysicsProcessesIn);
90 void Fill(const BDSTrajectoryPointHit* hit);
91 void Fill(const BDSHitEnergyDeposition* hit);
92
93 bool storeTurn = false;
94 bool storeLinks = false;
95 bool storeModelID = false;
96 bool storeLocal = false;
97 bool storeGlobal = false;
98 bool storeTime = false;
99 bool storeStepLength = false;
101 bool storePhysicsProcesses = false;
102#endif
103
104 ClassDef(BDSOutputROOTEventLoss,5);
105};
106
107#endif
Information recorded for a single piece of energy deposition.
Data stored for energy deposition hits per event.
bool storePreStepKineticEnergy
Whether to store pre step kinetic energy.
std::vector< float > y
Local coordinate.
bool storeLinks
Whether to store links between Eloss and model and trajectors.
std::vector< float > preStepKineticEnergy
Kinetic energy in GeV at pre step point.
bool storeLocal
Whether to store local coordinates.
bool storeGlobal
Whether to store global coordinates.
bool storeTurn
Store turn number.
std::vector< int > modelID
Geometry model index.
std::vector< float > weight
Weight associated with loss.
bool storeTime
Whether to store global time.
std::vector< int > trackID
TrackID that created the deposit.
std::vector< int > parentID
ParentID that created the deposit.
std::vector< float > z
Local coordinate.
std::vector< float > x
Local coordinate.
std::vector< int > turn
Turn number.
std::vector< float > stepLength
Step length taken for hit.
bool storeModelID
Whether to store the beam line index.
void Fill(const BDSOutputROOTEventLoss *other)
Fill from another instance.
std::vector< float > energy
Energy deposited in step.
bool storeStepLength
Whether to store step length.
std::vector< int > partID
ParticleID that create the deposit.
std::vector< float > S
Global curvilinear S coordinate.
std::vector< float > T
Global time (time since beginning of event).
A summary trajectory object of a loss point.
A Point in a trajectory with extra information.