BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSTrajectoryPoint.cc
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#include "BDSAuxiliaryNavigator.hh"
20#include "BDSDebug.hh"
21#include "BDSGlobalConstants.hh"
22#include "BDSPhysicalVolumeInfoRegistry.hh"
23#include "BDSPhysicalVolumeInfo.hh"
24#include "BDSStep.hh"
25#include "BDSTrajectoryPoint.hh"
26#include "BDSTrajectoryPointIon.hh"
27#include "BDSTrajectoryPointLocal.hh"
28#include "BDSTrajectoryPointLink.hh"
29#include "BDSPhysicalConstants.hh"
30#include "BDSPhysicsUtilities.hh"
31#include "BDSUtilities.hh"
32#ifdef BDSDEBUG_H
33#include "BDSProcessMap.hh"
34#endif
35
36#include "globals.hh"
37#include "G4Allocator.hh"
38#include "G4ProcessType.hh"
39#include "G4Step.hh"
40#include "G4ThreeVector.hh"
41#include "G4Track.hh"
42#include "G4TransportationProcessType.hh"
43#include "G4VProcess.hh"
44
45#include <ostream>
46
47class G4Material;
48
49G4Allocator<BDSTrajectoryPoint> bdsTrajectoryPointAllocator;
50
52
53// Don't use transform caching in the aux navigator as it's used for all over the geometry here.
55
57 G4TrajectoryPoint(G4ThreeVector())
58{
60}
61
63 G4bool storeExtrasLocal,
64 G4bool storeExtrasLink,
65 G4bool storeExtrasIon):
66 G4TrajectoryPoint(track->GetPosition())
67{
69
70 // Need to store the creator process
71 const G4VProcess* creatorProcess = track->GetCreatorProcess();
72 if (creatorProcess)
73 {
74 preProcessType = track->GetCreatorProcess()->GetProcessType();
75 preProcessSubType = track->GetCreatorProcess()->GetProcessSubType();
78 }
79
80 preWeight = track->GetWeight();
82 energyDeposit = 0.0; // does not lose any energy
83 preEnergy = track->GetKineticEnergy();
85 preMomentum = track->GetMomentum();
87 preGlobalTime = track->GetGlobalTime();
89 // when we construct a trajectory from a track it hasn't taken a step yet,
90 // so we don't know the material
91 material = nullptr;
92
93#ifdef BDSDEBUG
94 G4cout << __METHOD_NAME__ << "Process (main|sub) (" << BDSProcessMap::Instance()->GetProcessName(postProcessType, postProcessSubType) << ")" << G4endl;
95#endif
96
97 // s position for pre and post step point
98 // with a track, we're at the start and have no step - use 1nm for step to aid geometrical lookup
99 BDSStep localPosition = auxNavigator->ConvertToLocal(track->GetPosition(),
100 track->GetMomentumDirection(),
101 1*CLHEP::nm,
102 true);
103 prePosLocal = localPosition.PreStepPoint();
106 if (info)
107 {
108 G4double sCentre = info->GetSPos();
109 preS = sCentre + localPosition.PreStepPoint().z();
110 postS = sCentre + localPosition.PostStepPoint().z();
113 }
114
115 if (storeExtrasLocal)
116 {extraLocal = new BDSTrajectoryPointLocal(prePosLocal, localPosition.PostStepPoint());}
117
118 if (storeExtrasLink)
119 {StoreExtrasLink(track);}
120
121 if (storeExtrasIon)
122 {StoreExtrasIon(track);}
123}
124
126 G4bool storeExtrasLocal,
127 G4bool storeExtrasLink,
128 G4bool storeExtrasIon):
129 G4TrajectoryPoint(step->GetPostStepPoint()->GetPosition())
130{
132
133 const G4StepPoint* prePoint = step->GetPreStepPoint();
134 const G4StepPoint* postPoint = step->GetPostStepPoint();
135 const G4VProcess* preProcess = prePoint->GetProcessDefinedStep();
136 const G4VProcess* postProcess = postPoint->GetProcessDefinedStep();
137
138 if (preProcess)
139 {
140 preProcessType = preProcess->GetProcessType();
141 preProcessSubType = preProcess->GetProcessSubType();
142 }
143 if (postProcess)
144 {
145 postProcessType = postProcess->GetProcessType();
146 postProcessSubType = postProcess->GetProcessSubType();
147 }
148
149 preWeight = prePoint->GetWeight();
150 postWeight = postPoint->GetWeight();
151 energyDeposit = step->GetTotalEnergyDeposit();
152 preEnergy = prePoint->GetKineticEnergy();
153 postEnergy = postPoint->GetKineticEnergy();
154 preMomentum = prePoint->GetMomentum();
155 postMomentum = postPoint->GetMomentum();
156 preGlobalTime = prePoint->GetGlobalTime();
157 postGlobalTime = postPoint->GetGlobalTime();
158 material = prePoint->GetMaterial();
159
160#ifdef BDSDEBUG
161 G4cout << __METHOD_NAME__ << BDSProcessMap::Instance()->GetProcessName(postProcessType, postProcessSubType) << G4endl;
162#endif
163
164 // get local coordinates and volume for transform
165 BDSStep localPosition = auxNavigator->ConvertToLocal(step);
166 prePosLocal = localPosition.PreStepPoint();
167 postPosLocal = localPosition.PostStepPoint();
169 if (info)
170 {
171 G4double sCentre = info->GetSPos();
172 preS = sCentre + localPosition.PreStepPoint().z();
173 postS = sCentre + localPosition.PostStepPoint().z();
176 }
177
178 if (storeExtrasLocal)
179 {extraLocal = new BDSTrajectoryPointLocal(prePosLocal, localPosition.PostStepPoint());}
180
181 G4Track* track = step->GetTrack();
182 if (storeExtrasLink)
183 {StoreExtrasLink(track);}
184
185 if (storeExtrasIon)
186 {StoreExtrasIon(track);}
187}
188
190 G4TrajectoryPoint(static_cast<const G4TrajectoryPoint&>(other))
191{
192 extraLocal = other.extraLocal ? new BDSTrajectoryPointLocal(*other.extraLocal) : nullptr;
193 extraLink = other.extraLink ? new BDSTrajectoryPointLink(*other.extraLink) : nullptr;
194 extraIon = other.extraIon ? new BDSTrajectoryPointIon(*other.extraIon) : nullptr;
199
200 preWeight = other.preWeight;
201 postWeight = other.postWeight;
202 preEnergy = other.preEnergy;
203 postEnergy = other.postEnergy;
204 preMomentum = other.preMomentum;
207 preS = other.preS;
208 postS = other.postS;
212 beamline = other.beamline;
213 prePosLocal = other.prePosLocal;
215 material = other.material;
216}
217
218BDSTrajectoryPoint::~BDSTrajectoryPoint()
219{
220 delete extraLocal;
221 delete extraLink;
222 delete extraIon;
223}
224
226{
227 preProcessType = -1;
229 postProcessType = -1;
231 preWeight = -1.;
232 postWeight = -1.;
233 preEnergy = -1.;
234 postEnergy = -1.;
235 preMomentum = G4ThreeVector();
236 postMomentum = G4ThreeVector();
237 energyDeposit = 0.0;
238 preS = -1000;
239 postS = -1000;
240 preGlobalTime = 0;
241 postGlobalTime = 0;
242 beamlineIndex = -1;
243 beamline = nullptr;
244 prePosLocal = G4ThreeVector();
245 postPosLocal = G4ThreeVector();
246 material = nullptr;
247 extraLocal = nullptr;
248 extraLink = nullptr;
249 extraIon = nullptr;
250}
251
252void BDSTrajectoryPoint::StoreExtrasLink(const G4Track* track)
253{
254 const G4DynamicParticle* dynamicParticleDef = track->GetDynamicParticle();
255 G4double charge = dynamicParticleDef->GetCharge();
256 G4double rigidity = 0;
257 if (BDS::IsFinite(charge))
258 {rigidity = BDS::Rigidity(track->GetMomentum().mag(), charge);}
259 extraLink = new BDSTrajectoryPointLink((G4int)charge,
260 BDSGlobalConstants::Instance()->TurnsTaken(),
261 dynamicParticleDef->GetMass(),
262 rigidity);
263}
264
265void BDSTrajectoryPoint::StoreExtrasIon(const G4Track* track)
266{
267 const G4ParticleDefinition* particleDef = track->GetParticleDefinition();
268 const G4DynamicParticle* dynamicParticleDef = track->GetDynamicParticle();
269 G4bool isIon = BDS::IsIon(dynamicParticleDef);
270 G4int nElectrons = dynamicParticleDef->GetTotalOccupancy();
271 extraIon = new BDSTrajectoryPointIon(isIon,
272 particleDef->GetAtomicMass(),
273 particleDef->GetAtomicNumber(),
274 nElectrons);
275}
276
277
279{
280 // use general static function
284
285#ifdef BDSDEBUG
286 if (isScatteringPoint)
287 {
288 G4cout << "Interaction point found at " << GetPreS()/CLHEP::m << " m - "
290 }
291#endif
292 return isScatteringPoint;
293}
294
296{
297 // if prestep process doesn't exist it won't be set and will
298 // default to value in InitialiseVariables which is -1
299 G4bool preStep = (preProcessType != 1 /* transportation */ &&
300 preProcessType != 10 /* parallel world */);
301 G4bool posStep = (postProcessType != 1 /* transportation */ &&
302 postProcessType != 10 /* parallel world */);
303 return preStep || posStep;
304}
305
306std::ostream& operator<< (std::ostream& out, BDSTrajectoryPoint const &p)
307{
308 out << p.GetPosition();
309 return out;
310}
311
313{
314 return std::hypot(prePosLocal.x(), prePosLocal.y());
315}
316
318{
319 return std::hypot(postPosLocal.x(), postPosLocal.y());
320}
321
323{
324 const G4StepPoint* postPoint = step->GetPostStepPoint();
325 const G4VProcess* postProcess = postPoint->GetProcessDefinedStep();
326
327 G4double totalEnergyDeposit = step->GetTotalEnergyDeposit();
328
329 G4Track* t = step->GetTrack();
330 if (t->GetCurrentStepNumber() == 1 && t->GetStepLength() < 1e-5 && totalEnergyDeposit < 1e-5)
331 {return false;} // ignore the first really small step
332
333 G4int postProcessType = -1;
334 G4int postProcessSubType = -1;
335 if (postProcess)
336 {
337 postProcessType = postProcess->GetProcessType();
338 postProcessSubType = postProcess->GetProcessSubType();
339 }
340
342}
343
344G4bool BDSTrajectoryPoint::IsScatteringPoint(G4int postProcessType,
345 G4int postProcessSubType,
346 G4double totalEnergyDeposit)
347{
348 // test against things we want to exclude like tracking - these are not
349 // points of scattering
350 G4bool initialised = postProcessType != -1;
351 G4bool notTransportation = postProcessType != G4ProcessType::fTransportation;
352 G4bool notGeneral = (postProcessType != G4ProcessType::fGeneral) && (postProcessSubType != STEP_LIMITER);
353 G4bool notParallel = postProcessType != G4ProcessType::fParallel;
354 G4bool notUndefined = postProcessType != G4ProcessType::fNotDefined; // for crystal channelling
355
356 // energy can change in transportation step (EM)
357 if (totalEnergyDeposit > dEThresholdForScattering)
358 {return true;}
359
360 G4bool result = initialised && notTransportation && notGeneral && notParallel && notUndefined;
361 return result;
362}
Extra G4Navigator to get coordinate transforms.
BDSStep ConvertToLocal(G4Step const *const step, G4bool useCurvilinear=true) const
static BDSGlobalConstants * Instance()
Access method.
BDSPhysicalVolumeInfo * GetInfo(G4VPhysicalVolume *logicalVolume, G4bool isTunnel=false)
static BDSPhysicalVolumeInfoRegistry * Instance()
Singleton accessor.
A class holding any information pertaining to a particular physical volume in a BDSIM geant4 model.
BDSBeamline * GetBeamlineMassWorld() const
Accessor.
G4double GetSPos() const
Get the s position coordinate of the logical volume.
G4int GetBeamlineMassWorldIndex() const
Accessor.
static BDSProcessMap * Instance()
Singleton accessor.
G4String GetProcessName(const G4int &type, const G4int &subType=-1) const
Despatched function to operator() for getting the name of processes.
A simple class to represent the positions of a step.
Definition: BDSStep.hh:33
G4ThreeVector PostStepPoint() const
Accessor.
Definition: BDSStep.hh:43
G4VPhysicalVolume * VolumeForTransform() const
Accessor.
Definition: BDSStep.hh:44
G4ThreeVector PreStepPoint() const
Accessor.
Definition: BDSStep.hh:42
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.
G4ThreeVector postMomentum
Momentum of post-step point.
BDSBeamline * beamline
Beam line (if any) point belongs to (always mass world).
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.
G4int GetPostProcessType() const
Accessor.
G4double PostPosR() const
Return the transverse local radius in x,y.
G4double postGlobalTime
Time since event started of post-step point.
G4double postEnergy
Kinetic energy of post step point.
G4double PrePosR() const
Return the transverse local radius in x,y.
void StoreExtrasLink(const G4Track *track)
Utility function to prepare and fill extra link variables.
static BDSAuxiliaryNavigator * auxNavigator
G4double GetPreS() const
Accessor.
G4Material * material
Material point for pre-step point.
BDSTrajectoryPoint()
Default constructor.
G4double postWeight
Weight associated with post step point.
G4ThreeVector prePosLocal
Local coordinates of pre-step point.
G4int postProcessSubType
Process sub type of post step point.
G4ThreeVector postPosLocal
Local coordinates of post-step point.
G4int GetPostProcessSubType() const
Accessor.
G4bool NotTransportationLimitedStep() const
Return true if step isn't defined by transportation processes.
G4int preProcessType
Process type of pre-step point.
G4double energyDeposit
Total energy deposited during step.
G4bool IsScatteringPoint() const
G4double preWeight
Weight associated with pre-step point.
G4ThreeVector preMomentum
Momentum of pre-step point.
G4int preProcessSubType
Process sub type of pre-step point.
G4double preGlobalTime
Time since event started of pre-step point.
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
G4double Rigidity(G4double momentumMagnitude, G4double charge)
Calculate the rigidity for a total momentum and charge.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())
G4bool IsIon(const G4ParticleDefinition *particle)
Whether a particle is an ion. A proton is counted NOT as an ion.