BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSTrajectoryPoint.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 BDSTRAJECTORYPOINT_H
20#define BDSTRAJECTORYPOINT_H
21#include "BDSTrajectoryPointIon.hh"
22#include "BDSTrajectoryPointLocal.hh"
23#include "BDSTrajectoryPointLink.hh"
24
25#include "globals.hh" // geant4 types / globals
26#include "G4Allocator.hh"
27#include "G4TrajectoryPoint.hh"
28
29#include <ostream>
30
31class G4Material;
32class G4Step;
33class G4Track;
34
36class BDSBeamline;
37
44class BDSTrajectoryPoint: public G4TrajectoryPoint
45{
46public:
49
52 BDSTrajectoryPoint(const G4Step* step,
53 G4bool storeExtrasLocalIn,
54 G4bool storeExtrasLinkIn,
55 G4bool storeExtrasIonIn);
56
59 BDSTrajectoryPoint(const G4Track* track,
60 G4bool storeExtrasLocalIn,
61 G4bool storeExtrasLinkIn,
62 G4bool storeExtrasIonIn);
63
66
67 virtual ~BDSTrajectoryPoint();
68
69 inline void *operator new(size_t);
70 inline void operator delete(void *aTrajectoryPoint);
71 inline int operator==(const BDSTrajectoryPoint& right) const {return (this==&right);};
72
73 inline void DeleteExtraLocal() {delete extraLocal; extraLocal = nullptr;}
74 inline void DeleteExtraLinks() {delete extraLink; extraLink = nullptr;}
75 inline void DeleteExtraIon() {delete extraIon; extraIon = nullptr;}
76
82 G4bool IsScatteringPoint() const;
83
85 G4bool NotTransportationLimitedStep() const;
86
88 static G4bool IsScatteringPoint(const G4Step* step);
89 static G4bool IsScatteringPoint(G4int postProcessType,
91 G4double totalEnergyDeposit);
93
95 inline G4int GetPreProcessType() const {return preProcessType;}
96 inline G4int GetPreProcessSubType() const {return preProcessSubType;}
97 inline G4int GetPostProcessType() const {return postProcessType;}
98 inline G4int GetPostProcessSubType() const {return postProcessSubType;}
99 inline G4double GetPreWeight() const {return preWeight;}
100 inline G4double GetPostWeight() const {return postWeight;}
101 inline G4double GetPreEnergy() const {return preEnergy;}
102 inline G4double GetPostEnergy() const {return postEnergy;}
103 inline G4double GetEnergyDeposit() const {return energyDeposit;}
104 inline G4ThreeVector GetPreMomentum() const {return preMomentum;}
105 inline G4ThreeVector GetPostMomentum() const {return postMomentum;}
106 inline G4double GetPreS() const {return preS;}
107 inline G4double GetPostS() const {return postS;}
108 inline G4double GetPreGlobalTime() const {return preGlobalTime;}
109 inline G4double GetPostGlobalTime() const {return postGlobalTime;}
110 inline G4int GetBeamLineIndex() const {return beamlineIndex;}
111 inline BDSBeamline* GetBeamLine() const {return beamline;}
112 inline G4ThreeVector GetPrePosLocal() const {return prePosLocal;}
113 inline G4ThreeVector GetPostPosLocal() const {return postPosLocal;}
114 inline G4Material* GetMaterial() const {return material;}
115
118 inline void SetMaterial(G4Material* materialIn) {material = materialIn;}
119
121 inline G4ThreeVector GetPositionLocal() const {return extraLocal ? extraLocal->positionLocal : G4ThreeVector();}
122 inline G4ThreeVector GetMomentumLocal() const {return extraLocal ? extraLocal->momentumLocal : G4ThreeVector();}
124
126 inline G4int GetCharge() const {return extraLink ? extraLink->charge : 0;}
127 inline G4double GetKineticEnergy() const {return preEnergy;}
128 inline G4int GetTurnsTaken() const {return extraLink ? extraLink->turnsTaken : 0;}
129 inline G4double GetMass() const {return extraLink ? extraLink->mass : 0;}
130 inline G4double GetRigidity() const {return extraLink ? extraLink->rigidity : 0;}
132
134 inline G4bool GetIsIon() const {return extraIon ? extraIon->isIon : false;}
135 inline G4int GetIonA() const {return extraIon ? extraIon->ionA : 0;}
136 inline G4int GetIonZ() const {return extraIon ? extraIon->ionZ : 0;}
137 inline G4int GetNElectrons() const {return extraIon ? extraIon->nElectrons : 0;}
139
141 G4double PrePosR() const;
142 G4double PostPosR() const;
144
146 friend std::ostream& operator<< (std::ostream &out, BDSTrajectoryPoint const &p);
147
148 // logic taken from BDSExtent - implement only one operator; rest come for free
150 inline G4bool operator< (const BDSTrajectoryPoint& other) const;
151 inline G4bool operator> (const BDSTrajectoryPoint& other) const {return other < (*this);}
152 inline G4bool operator<=(const BDSTrajectoryPoint& other) const {return !((*this) > other);}
153 inline G4bool operator>=(const BDSTrajectoryPoint& other) const {return !((*this) < other);}
155
156 BDSTrajectoryPointLocal* extraLocal;
157 BDSTrajectoryPointLink* extraLink;
158 BDSTrajectoryPointIon* extraIon;
159
165
166private:
169 void InitialiseVariables();
170
172 void StoreExtrasLink(const G4Track* track);
173
175 void StoreExtrasIon(const G4Track* track);
176
181
182 G4double preWeight;
183 G4double postWeight;
184 G4double preEnergy;
185 G4double postEnergy;
186 G4ThreeVector preMomentum;
187 G4ThreeVector postMomentum;
188 G4double energyDeposit;
189 G4double preS;
190 G4double postS;
191 G4double preGlobalTime;
192 G4double postGlobalTime;
195 G4ThreeVector prePosLocal;
196 G4ThreeVector postPosLocal;
197 G4Material* material;
198
202};
203
204extern G4Allocator<BDSTrajectoryPoint> bdsTrajectoryPointAllocator;
205
206inline void* BDSTrajectoryPoint::operator new(size_t)
207{
208 void *aTrajectoryPoint;
209 aTrajectoryPoint = (void *) bdsTrajectoryPointAllocator.MallocSingle();
210 return aTrajectoryPoint;
211}
212
213inline void BDSTrajectoryPoint::operator delete(void *aTrajectoryPoint)
214{
215 bdsTrajectoryPointAllocator.FreeSingle((BDSTrajectoryPoint *) aTrajectoryPoint);
216}
217
218inline G4bool BDSTrajectoryPoint::operator< (const BDSTrajectoryPoint& other) const
219{
220 // can't test position without knowledge of beam line direction etc - too difficult / inaccurate
221 // for now, this is simplistic
222 // TODO deal with multiple beam lines and s coordinate change at join point
223 return (preS < other.preS) && (preGlobalTime < other.preGlobalTime);
224}
225
226#endif
Extra G4Navigator to get coordinate transforms.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
Extra information recorded for a single piece of energy deposition.
Extra information recorded for a single piece of energy deposition.
A Point in a trajectory with extra information.
G4double preS
Global curvilinear S coordinate of pre-step point.
G4double GetKineticEnergy() const
Accessor for the extra information links.
G4double GetPreEnergy() const
Accessor.
G4double GetPostS() const
Accessor.
G4int GetCharge() const
Accessor for the extra information links.
G4ThreeVector postMomentum
Momentum of post-step point.
BDSBeamline * beamline
Beam line (if any) point belongs to (always mass world).
BDSBeamline * GetBeamLine() const
Accessor.
void StoreExtrasIon(const G4Track *track)
Utility function to prepare and fill extra ion variables.
G4int beamlineIndex
Index to beam line element in the mass world beam line.
G4double GetPreGlobalTime() const
Accessor.
G4bool operator>=(const BDSTrajectoryPoint &other) const
Comparison operator.
friend std::ostream & operator<<(std::ostream &out, BDSTrajectoryPoint const &p)
Output stream.
G4int GetPostProcessType() const
Accessor.
G4bool operator<(const BDSTrajectoryPoint &other) const
Comparison operator.
G4ThreeVector GetMomentumLocal() const
Accessor.
void SetMaterial(G4Material *materialIn)
G4int GetPreProcessSubType() const
Accessor.
G4bool GetIsIon() const
Accessor for the extra information ions.
G4double PostPosR() const
Return the transverse local radius in x,y.
G4ThreeVector GetPostPosLocal() const
Accessor.
G4double GetPreWeight() const
Accessor.
G4double postGlobalTime
Time since event started of post-step point.
G4double GetRigidity() const
Accessor for the extra information links.
G4double postEnergy
Kinetic energy of post step point.
G4int GetTurnsTaken() const
Accessor for the extra information links.
G4double PrePosR() const
Return the transverse local radius in x,y.
G4double GetEnergyDeposit() const
Accessor.
G4bool operator<=(const BDSTrajectoryPoint &other) const
Comparison operator.
void StoreExtrasLink(const G4Track *track)
Utility function to prepare and fill extra link variables.
G4double GetMass() const
Accessor for the extra information links.
static BDSAuxiliaryNavigator * auxNavigator
G4double GetPreS() const
Accessor.
G4Material * material
Material point for pre-step point.
G4bool operator>(const BDSTrajectoryPoint &other) const
Comparison operator.
BDSTrajectoryPoint()
Default constructor.
G4double postWeight
Weight associated with post step point.
G4int GetPreProcessType() const
Accessor.
G4ThreeVector prePosLocal
Local coordinates of pre-step point.
G4int GetIonZ() const
Accessor for the extra information ions.
G4ThreeVector GetPositionLocal() const
Accessor for the extra information local.
G4int postProcessSubType
Process sub type of post step point.
G4ThreeVector postPosLocal
Local coordinates of post-step point.
G4double GetPostEnergy() const
Accessor.
G4double GetPostWeight() const
Accessor.
G4int GetPostProcessSubType() const
Accessor.
G4bool NotTransportationLimitedStep() const
Return true if step isn't defined by transportation processes.
G4ThreeVector GetPostMomentum() const
Accessor.
G4int preProcessType
Process type of pre-step point.
G4double energyDeposit
Total energy deposited during step.
G4bool IsScatteringPoint() const
G4int GetBeamLineIndex() const
Accessor.
G4double preWeight
Weight associated with pre-step point.
G4int GetNElectrons() const
Accessor for the extra information ions.
G4int GetIonA() const
Accessor for the extra information ions.
G4double GetPostGlobalTime() const
Accessor.
G4ThreeVector preMomentum
Momentum of pre-step point.
G4Material * GetMaterial() const
Accessor.
G4int preProcessSubType
Process sub type of pre-step point.
G4double preGlobalTime
Time since event started of pre-step point.
G4ThreeVector GetPreMomentum() const
Accessor.
G4ThreeVector GetPrePosLocal() const
Accessor.
G4double preEnergy
Kinetic energy of pre-step point.
G4double postS
Global curvilinear S coordinate of post step point.
G4int postProcessType
Process type of post step point.
static G4double dEThresholdForScattering