BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSSDEnergyDeposition.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 "BDSHitEnergyDeposition.hh"
21#include "BDSSDEnergyDeposition.hh"
22#include "BDSDebug.hh"
23#include "BDSGlobalConstants.hh"
24#include "BDSPhysicalVolumeInfo.hh"
25#include "BDSPhysicalVolumeInfoRegistry.hh"
26#include "BDSStep.hh"
27#include "BDSUtilities.hh"
28
29#include "globals.hh" // geant4 types / globals
30#include "G4AffineTransform.hh"
31#include "G4Event.hh"
32#include "G4EventManager.hh"
33#include "G4LogicalVolume.hh"
34#include "G4ParticleDefinition.hh"
35#include "G4SDManager.hh"
36#include "G4Step.hh"
37#include "G4ThreeVector.hh"
38#include "G4Track.hh"
39#include "G4VPhysicalVolume.hh"
40#include "G4VTouchable.hh"
41#include "Randomize.hh"
42
43BDSSDEnergyDeposition::BDSSDEnergyDeposition(const G4String& name,
44 G4bool storeExtrasIn,
45 G4bool killedParticleMassAddedToElossIn):
46 BDSSensitiveDetector("energy_counter/"+name),
47 storeExtras(storeExtrasIn),
48 killedParticleMassAddedToEloss(killedParticleMassAddedToElossIn),
49 colName(name),
50 hits(nullptr),
51 HCIDe(-1),
52 auxNavigator(new BDSAuxiliaryNavigator())
53{
54 collectionName.insert(colName);
55}
56
57BDSSDEnergyDeposition::~BDSSDEnergyDeposition()
58{
59 delete auxNavigator;
60}
61
62void BDSSDEnergyDeposition::Initialize(G4HCofThisEvent* HCE)
63{
65 if (HCIDe < 0)
66 {HCIDe = G4SDManager::GetSDMpointer()->GetCollectionID(hits);}
67 HCE->AddHitsCollection(HCIDe,hits);
68
69#ifdef BDSDEBUG
70 G4cout << __METHOD_NAME__ << "Hits Collection ID: " << HCIDe << G4endl;
71#endif
72}
73
75 G4TouchableHistory* /*th*/)
76{
77 // Get the energy deposited along the step
78 G4double energy = aStep->GetTotalEnergyDeposit();
79
80 //if the energy is 0, don't do anything
81 if (!BDS::IsFinite(energy))
82 {return false;}
83
84 G4Track* track = aStep->GetTrack();
85 G4int parentID = track->GetParentID(); // needed later on too
86 G4int ptype = track->GetDefinition()->GetPDGEncoding();
87
88 // step points - used many times
89 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
90 G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
91
92 G4double preStepKineticEnergy = preStepPoint->GetKineticEnergy();
93
94 // attribute the energy deposition to a uniformly random position along the step - correct!
95 // random distance - store to use twice to ensure global and local represent the same point
96 G4double randDist = G4UniformRand();
97
98 // global coordinate positions of the step
99 const G4ThreeVector& posbefore = preStepPoint->GetPosition();
100 const G4ThreeVector& posafter = postStepPoint->GetPosition();
101 G4ThreeVector eDepPos = posbefore + randDist*(posafter - posbefore);
102
103 // calculate local coordinates
104 BDSStep stepLocal = auxNavigator->ConvertToLocal(aStep);
105 const G4ThreeVector& posbeforelocal = stepLocal.PreStepPoint();
106 const G4ThreeVector& posafterlocal = stepLocal.PostStepPoint();
107 G4ThreeVector eDepPosLocal = posbeforelocal + randDist*(posafterlocal - posbeforelocal);
108 G4double stepLength = (posafterlocal - posbeforelocal).mag();
109
110 // global
111 G4double X = eDepPos.x();
112 G4double Y = eDepPos.y();
113 G4double Z = eDepPos.z();
114 // local
115 G4double x = eDepPosLocal.x();
116 G4double y = eDepPosLocal.y();
117 G4double z = eDepPosLocal.z();
118
119 // Just as the energy deposition is attributed to a uniformly random
120 // point between the preStep and the postStep positions, attribute the
121 // deposition to random time between preStep and postStep times,
122 // using the same random number as for the position.
123 G4double preGlobalTime = preStepPoint->GetGlobalTime();
124 G4double postGlobalTime = postStepPoint->GetGlobalTime();
125 G4double globalTime = preGlobalTime + randDist * (postGlobalTime - preGlobalTime);
126
127 // get the s coordinate (central s + local z)
128 // volume is from curvilinear coordinate parallel geometry
129 BDSPhysicalVolumeInfo* theInfo = BDSPhysicalVolumeInfoRegistry::Instance()->GetInfo(stepLocal.VolumeForTransform());
130 G4int beamlineIndex = -1;
131
132 // declare lambda for updating parameters if info found (avoid duplication of code)
133 G4double sBefore = -1000;
134 G4double sAfter = -1000;
135 auto UpdateParams = [&](BDSPhysicalVolumeInfo* info)
136 {
137 G4double sCentre = info->GetSPos();
138 sAfter = sCentre + posafterlocal.z();
139 sBefore = sCentre + posbeforelocal.z();
140 beamlineIndex = info->GetBeamlineIndex();
141 };
142
143 if (theInfo)
144 {UpdateParams(theInfo);}
145 else
146 {
147 // Try again but with the pre step point only
148 G4ThreeVector unitDirection = (posafter - posbefore).unit();
149 BDSStep stepLocal2 = auxNavigator->ConvertToLocal(posbefore, unitDirection);
151 if (theInfo)
152 {UpdateParams(theInfo);}
153 else
154 {
155 // Try yet again with just a slight shift (100um is bigger than any padding space).
156 G4ThreeVector shiftedPos = posbefore + 0.1*CLHEP::mm*unitDirection;
157 stepLocal2 = auxNavigator->ConvertToLocal(shiftedPos, unitDirection);
159 if (theInfo)
160 {UpdateParams(theInfo);}
161 else
162 {
163#ifdef BDSDEBUG
164 G4cerr << "No volume info for ";
165 auto vol = stepLocal.VolumeForTransform();
166 if (vol)
167 {G4cerr << vol->GetName() << G4endl;}
168 else
169 {G4cerr << "Unknown" << G4endl;}
170#endif
171 // unphysical default value to allow easy identification in output
172 sAfter = -1000;
173 sBefore = -1000;
174 beamlineIndex = -2;
175 }
176 }
177 }
178
179 G4double sHit = sBefore + randDist*(sAfter - sBefore);
180
181 G4double weight = track->GetWeight();
182 G4int trackID = track->GetTrackID();
183 G4int turnsTaken = BDSGlobalConstants::Instance()->TurnsTaken();
184
185 G4int postStepProcessType = -1;
186 G4int postStepProcessSubType = -1;
187 if (storeExtras)
188 {// physics processes
189 const G4StepPoint* postPoint = aStep->GetPostStepPoint();
190 const G4VProcess* postProcess = postPoint->GetProcessDefinedStep();
191 if (postProcess)
192 {
193 postStepProcessType = postProcess->GetProcessType();
194 postStepProcessSubType = postProcess->GetProcessSubType();
195 }
196 }
197
198 //create hits and put in hits collection of the event
200 sHit,
201 weight,
203 preStepKineticEnergy,
204 X, Y, Z,
205 x, y, z,
206 globalTime,
207 ptype,
208 trackID,
209 parentID,
210 turnsTaken,
211 stepLength,
212 beamlineIndex,
213 postStepProcessType,
214 postStepProcessSubType);
215
216 // don't worry, won't add 0 energy tracks as filtered at top by if statement
217 hits->insert(hit);
218
219 return true;
220}
221
222G4bool BDSSDEnergyDeposition::ProcessHitsTrack(const G4Track* track,
223 G4TouchableHistory* /*th*/)
224{
225 G4int parentID = track->GetParentID(); // needed later on too
226 G4int ptype = track->GetDefinition()->GetPDGEncoding();
227
228 // depending on our option, include the rest mass - for backwards compatibility
229 G4double energy = killedParticleMassAddedToEloss ? track->GetTotalEnergy() : track->GetKineticEnergy();
230
231 G4double globalTime = track->GetGlobalTime();
232 G4double weight = track->GetWeight();
233 G4int trackID = track->GetTrackID();
234 G4double preStepKineticEnergy = track->GetKineticEnergy();
235
236 // if the energy is 0, don't do anything
237 if (!BDS::IsFinite(energy))
238 {return false;}
239
240 G4double stepLength = 0;
241 const G4ThreeVector& posGlobal = track->GetPosition();
242 G4double X = posGlobal.x();
243 G4double Y = posGlobal.y();
244 G4double Z = posGlobal.z();
245
246 // calculate local coordinates
247 const G4ThreeVector& momGlobalUnit = track->GetMomentumDirection();
248 BDSStep stepLocal = auxNavigator->ConvertToLocal(posGlobal, momGlobalUnit, 1*CLHEP::mm, true, 1*CLHEP::mm);
249 G4ThreeVector posLocal = stepLocal.PreStepPoint();
250 // local
251 G4double x = posLocal.x();
252 G4double y = posLocal.y();
253 G4double z = posLocal.z();
254
255 // get the s coordinate (central s + local z)
256 // volume is from curvilinear coordinate parallel geometry
258 G4int beamlineIndex = -1;
259
260 // declare lambda for updating parameters if info found (avoid duplication of code)
261 G4double sBefore = -1000;
262 G4double sAfter = -1000;
263 auto UpdateParams = [&](BDSPhysicalVolumeInfo* info)
264 {
265 G4double sCentre = info->GetSPos();
266 sAfter = sCentre + posLocal.z();
267 sBefore = sCentre + posLocal.z();
268 beamlineIndex = info->GetBeamlineIndex();
269 };
270
271 if (theInfo)
272 {UpdateParams(theInfo);}
273 else
274 {
275 // Try yet again with just a slight shift (100um is bigger than any padding space).
276 G4ThreeVector shiftedPos = posGlobal + 0.1*CLHEP::mm * momGlobalUnit;
277 BDSStep stepLocal2 = auxNavigator->ConvertToLocal(shiftedPos, momGlobalUnit);
279 if (theInfo)
280 {UpdateParams(theInfo);}
281 else
282 {
283#ifdef BDSDEBUG
284 G4cerr << "No volume info for ";
285 auto vol = stepLocal.VolumeForTransform();
286 if (vol)
287 {G4cerr << vol->GetName() << G4endl;}
288 else
289 {G4cerr << "Unknown" << G4endl;}
290#endif
291 // unphysical default value to allow easy identification in output
292 sAfter = -1000;
293 sBefore = -1000;
294 beamlineIndex = -2;
295 }
296 }
297 G4double sHit = sBefore; // duplicate
298
299 G4int turnsTaken = BDSGlobalConstants::Instance()->TurnsTaken();
300
301 G4int postStepProcessType = -1;
302 G4int postStepProcessSubType = -1;
303 if (storeExtras)
304 {// physics processes
305 const G4Step* aStep = track->GetStep();
306 if (aStep)
307 {
308 const G4StepPoint* postPoint = aStep->GetPostStepPoint();
309 const G4VProcess* postProcess = postPoint->GetProcessDefinedStep();
310 if (postProcess)
311 {
312 postStepProcessType = postProcess->GetProcessType();
313 postStepProcessSubType = postProcess->GetProcessSubType();
314 }
315 }
316 }
317
318 //create hits and put in hits collection of the event
320 sHit,
321 weight,
323 preStepKineticEnergy,
324 X, Y, Z,
325 x, y, z,
326 globalTime,
327 ptype,
328 trackID,
329 parentID,
330 turnsTaken,
331 stepLength,
332 beamlineIndex,
333 postStepProcessType,
334 postStepProcessSubType);
335
336 // don't worry, won't add 0 energy tracks as filtered at top by if statement
337 hits->insert(hit);
338
339 return true;
340}
341
343{
344 BDSHitEnergyDeposition* lastHit = hits->GetVector()->back();
345 return dynamic_cast<G4VHit*>(lastHit);
346}
Extra G4Navigator to get coordinate transforms.
BDSStep ConvertToLocal(G4Step const *const step, G4bool useCurvilinear=true) const
static BDSGlobalConstants * Instance()
Access method.
Information recorded for a single piece of energy deposition.
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.
G4double GetSPos() const
Get the s position coordinate of the logical volume.
G4bool killedParticleMassAddedToEloss
In the case of a G4Track being deposited.
G4String colName
Collection name.
G4bool storeExtras
Whether to store extra information.
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *th)
BDSAuxiliaryNavigator * auxNavigator
Navigator for checking points in read out geometry.
virtual G4VHit * last() const
Provide access to last hit.
virtual G4bool ProcessHitsTrack(const G4Track *track, G4TouchableHistory *th)
Virtual class to define interface for ordered multi-sensitive detector.
A simple class to represent the positions of a step.
Definition: BDSStep.hh:33
G4VPhysicalVolume * VolumeForTransform() const
Accessor.
Definition: BDSStep.hh:44
G4ThreeVector PreStepPoint() const
Accessor.
Definition: BDSStep.hh:42
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())
const char * GetName(int)
Name of element.
Definition: python.cc:38