Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SteppingManager.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4SteppingManager
27//
28// Class description:
29//
30// This is the class which plays an essential role in tracking particles.
31// It takes cares of all message passing between objects in the different
32// categories (for example, geometry - including transportation, interactions
33// in matter, etc). Its public method 'stepping' steers to step the particle.
34// Only used within the Geant4 kernel.
35
36// Contact:
37// Questions and comments to this code should be sent to
38// Katsuya Amako (e-mail: [email protected])
39// Takashi Sasaki (e-mail: [email protected])
40//
41// History:
42// 12.03.1998, H.Kurashige - modified for use of G4ParticleChange
43// --------------------------------------------------------------------
44#ifndef G4SteppingManager_hh
45#define G4SteppingManager_hh 1
46
47#include "G4LogicalVolume.hh" // Include from 'geometry'
48#include "G4Navigator.hh" // Include from 'geometry'
49#include "G4NoProcess.hh" // Include from 'processes'
50#include "G4ProcessManager.hh" // Include from 'processes'
51#include "G4Step.hh" // Include from 'tracking'
52#include "G4StepPoint.hh" // Include from 'tracking'
53#include "G4StepStatus.hh" // Include from 'tracking'
54#include "G4TouchableHandle.hh" // Include from 'geometry'
55#include "G4Track.hh" // Include from 'tracking'
56#include "G4TrackStatus.hh" // Include from 'tracking'
57#include "G4TrackVector.hh" // Include from 'tracking'
58#include "G4UserSteppingAction.hh" // Include from 'tracking'
59#include "G4VPhysicalVolume.hh" // Include from 'geometry'
60#include "G4VSteppingVerbose.hh" // Include from 'tracking'
61#include "Randomize.hh" // Include from 'global'
62#include "globals.hh" // Include from 'global'
63
64#include <iomanip> // Include from 'system'
65#include <vector> // Include from 'system'
66
67using G4SelectedAtRestDoItVector = std::vector<G4int>;
68using G4SelectedAlongStepDoItVector = std::vector<G4int>;
69using G4SelectedPostStepDoItVector = std::vector<G4int>;
70
72
74{
75 public:
77
78 public:
79 // Constructor/Destructor
80
81 // SteppingManger should be dynamically allocated, therefore
82 // you need to invoke new() when you call this constructor.
83 // "Secodary track vector" will be dynamically created by this
84 // constructor. G4UserSteppingAction will be also created
85 // in this constructor, and "this" pointer will be passed to
86 // G4UserSteppingAction.
89
90 // Get/Set functions
91
92 const G4TrackVector* GetSecondary() const;
93 void SetUserAction(G4UserSteppingAction* apAction);
94 G4Track* GetTrack() const;
95 void SetVerboseLevel(G4int vLevel);
97 G4Step* GetStep() const;
98 void SetNavigator(G4Navigator* value);
99
100 // Other member functions
101
102 // Steers to move the give particle from the TrackingManger by one Step.
104
105 // Sets up initial track information (enegry, position, etc) to
106 // the PreStepPoint of the G4Step. This funciton has to be called
107 // just once before the stepping loop in the "TrackingManager".
108 void SetInitialStep(G4Track* valueTrack);
109
110 // Get methods
111 void GetProcessNumber();
125 G4Step* GetfStep();
139 std::size_t GetfAtRestDoItProcTriggered();
140 std::size_t GetfAlongStepDoItProcTriggered();
141 std::size_t GetfPostStepDoItProcTriggered();
147 std::size_t GetMAXofAtRestLoops();
148 std::size_t GetMAXofAlongStepLoops();
149 std::size_t GetMAXofPostStepLoops();
160
161 private:
162 // Member functions
163
164 // Calculate corresponding physical length from the mean free path
165 // left for each discrete physics process. The minimum allowable
166 // step for each continuous process will be also calculated.
167 void DefinePhysicalStepLength();
168
169 void InvokeAtRestDoItProcs();
170 void InvokeAlongStepDoItProcs();
171 void InvokePostStepDoItProcs();
172 void InvokePSDIP(size_t); //
173 G4int ProcessSecondariesFromParticleChange();
174
175 // Return the estimated safety value at the PostStepPoint
176 G4double CalculateSafety();
177
178 // Member data
179
180 static const size_t SizeOfSelectedDoItVector = 100;
181
182 G4bool KillVerbose = false;
183
184 G4UserSteppingAction* fUserSteppingAction = nullptr;
185
186 G4VSteppingVerbose* fVerbose = nullptr;
187
188 G4double PhysicalStep = 0.0;
189 G4double GeometricalStep = 0.0;
190 G4double CorrectedStep = 0.0;
191 G4bool PreStepPointIsGeom = false;
192 G4bool FirstStep = false;
193 G4StepStatus fStepStatus = fUndefined;
194
195 G4double TempInitVelocity = 0.0;
196 G4double TempVelocity = 0.0;
197 G4double Mass = 0.0;
198
199 G4double sumEnergyChange = 0.0;
200
201 G4VParticleChange* fParticleChange = nullptr;
202 G4Track* fTrack = nullptr;
203 G4TrackVector* fSecondary = nullptr;
204 G4Step* fStep = nullptr;
205 G4StepPoint* fPreStepPoint = nullptr;
206 G4StepPoint* fPostStepPoint = nullptr;
207
208 G4VPhysicalVolume* fCurrentVolume = nullptr;
209 G4VSensitiveDetector* fSensitive = nullptr;
210 G4VProcess* fCurrentProcess = nullptr; // Pointer to process of which DoIt() or
211 // GetPhysicalInteractionLength() has been just executed.
212
213 G4ProcessVector* fAtRestDoItVector = nullptr;
214 G4ProcessVector* fAlongStepDoItVector = nullptr;
215 G4ProcessVector* fPostStepDoItVector = nullptr;
216
217 G4ProcessVector* fAtRestGetPhysIntVector = nullptr;
218 G4ProcessVector* fAlongStepGetPhysIntVector = nullptr;
219 G4ProcessVector* fPostStepGetPhysIntVector = nullptr;
220
221 std::size_t MAXofAtRestLoops = 0;
222 std::size_t MAXofAlongStepLoops = 0;
223 std::size_t MAXofPostStepLoops = 0;
224
225 std::size_t fAtRestDoItProcTriggered = 0;
226 std::size_t fAlongStepDoItProcTriggered = 0;
227 std::size_t fPostStepDoItProcTriggered = 0;
228
229 G4int fN2ndariesAtRestDoIt = 0;
230 G4int fN2ndariesAlongStepDoIt = 0;
231 G4int fN2ndariesPostStepDoIt = 0;
232 // These are the numbers of secondaries generated by the process
233 // just executed.
234
235 G4Navigator* fNavigator = nullptr;
236
237 G4int verboseLevel = 0;
238
239 G4SelectedAtRestDoItVector* fSelectedAtRestDoItVector = nullptr;
240 G4SelectedAlongStepDoItVector* fSelectedAlongStepDoItVector = nullptr;
241 G4SelectedPostStepDoItVector* fSelectedPostStepDoItVector = nullptr;
242
243 G4double fPreviousStepSize = 0.0;
244
245 G4TouchableHandle fTouchableHandle;
246
247 G4SteppingControl StepControlFlag = NormalCondition;
248
249 G4double kCarTolerance = 0.0; // Cached geometrical tolerance on surface
250 G4double proposedSafety = 0.0; // This keeps the minimum safety value proposed by AlongStepGPILs.
251 G4ThreeVector endpointSafOrigin;
252 G4double endpointSafety = 0.0; // To get the true safety value at the PostStepPoint, you have to
253 // subtract the distance to 'endpointSafOrigin' from this value.
254 G4double physIntLength = 0.0;
255 G4ForceCondition fCondition = InActivated;
257 // Above three variables are for the method
258 // DefinePhysicalStepLength(). To pass these information to
259 // the method Verbose, they are kept at here. Need a more
260 // elegant mechanism.
261
262 G4NoProcess const* fNoProcess = nullptr; // Used in InvokeAtRestDoItProcs() method to flag the
263 // process of any stable ion at rest.
264};
265
266//*******************************************************************
267//
268// Inline functions
269//
270//*******************************************************************
271
272inline G4double G4SteppingManager::GetPhysicalStep() { return PhysicalStep; }
273
274inline G4double G4SteppingManager::GetGeometricalStep() { return GeometricalStep; }
275
276inline G4double G4SteppingManager::GetCorrectedStep() { return CorrectedStep; }
277
278inline G4bool G4SteppingManager::GetPreStepPointIsGeom() { return PreStepPointIsGeom; }
279
280inline G4bool G4SteppingManager::GetFirstStep() { return FirstStep; }
281
282inline G4StepStatus G4SteppingManager::GetfStepStatus() { return fStepStatus; }
283
284inline G4double G4SteppingManager::GetTempInitVelocity() { return TempInitVelocity; }
285
286inline G4double G4SteppingManager::GetTempVelocity() { return TempVelocity; }
287
288inline G4double G4SteppingManager::GetMass() { return Mass; }
289
290inline G4double G4SteppingManager::GetsumEnergyChange() { return sumEnergyChange; }
291
292inline G4VParticleChange* G4SteppingManager::GetfParticleChange() { return fParticleChange; }
293
294inline G4Track* G4SteppingManager::GetfTrack() { return fTrack; }
295
297
298inline G4Step* G4SteppingManager::GetfStep() { return fStep; }
299
300inline G4StepPoint* G4SteppingManager::GetfPreStepPoint() { return fPreStepPoint; }
301
302inline G4StepPoint* G4SteppingManager::GetfPostStepPoint() { return fPostStepPoint; }
303
304inline G4VPhysicalVolume* G4SteppingManager::GetfCurrentVolume() { return fCurrentVolume; }
305
307
308inline G4VProcess* G4SteppingManager::GetfCurrentProcess() { return fCurrentProcess; }
309
310inline G4ProcessVector* G4SteppingManager::GetfAtRestDoItVector() { return fAtRestDoItVector; }
311
313{
314 return fAlongStepDoItVector;
315}
316
317inline G4ProcessVector* G4SteppingManager::GetfPostStepDoItVector() { return fPostStepDoItVector; }
318
320{
321 return fAtRestGetPhysIntVector;
322}
323
325{
326 return fAlongStepGetPhysIntVector;
327}
328
330{
331 return fPostStepGetPhysIntVector;
332}
333
334inline size_t G4SteppingManager::GetMAXofAtRestLoops() { return MAXofAtRestLoops; }
335
336inline size_t G4SteppingManager::GetMAXofAlongStepLoops() { return MAXofAlongStepLoops; }
337
338inline size_t G4SteppingManager::GetMAXofPostStepLoops() { return MAXofPostStepLoops; }
339
340inline size_t G4SteppingManager::GetfAtRestDoItProcTriggered() { return fAtRestDoItProcTriggered; }
341
343{
344 return fAtRestDoItProcTriggered;
345}
346
348{
349 return fPostStepDoItProcTriggered;
350}
351
352inline G4int G4SteppingManager::GetfN2ndariesAtRestDoIt() { return fN2ndariesAtRestDoIt; }
353
354inline G4int G4SteppingManager::GetfN2ndariesAlongStepDoIt() { return fN2ndariesAlongStepDoIt; }
355
356inline G4int G4SteppingManager::GetfN2ndariesPostStepDoIt() { return fN2ndariesPostStepDoIt; }
357
358inline G4Navigator* G4SteppingManager::GetfNavigator() { return fNavigator; }
359
360inline G4int G4SteppingManager::GetverboseLevel() { return verboseLevel; }
361
363{
364 return fSelectedAtRestDoItVector;
365}
366
368{
369 return fSelectedAlongStepDoItVector;
370}
371
373{
374 return fSelectedPostStepDoItVector;
375}
376
377inline G4double G4SteppingManager::GetfPreviousStepSize() { return fPreviousStepSize; }
378
379inline const G4TouchableHandle& G4SteppingManager::GetTouchableHandle() { return fTouchableHandle; }
380
382
383inline G4double G4SteppingManager::GetphysIntLength() { return physIntLength; }
384
386
387inline G4GPILSelection G4SteppingManager::GetfGPILSelection() { return fGPILSelection; }
388
390{
391 return fStep->GetSecondary();
392}
393
394inline void G4SteppingManager::SetNavigator(G4Navigator* value) { fNavigator = value; }
395
397{
398 fUserSteppingAction = apAction;
399}
400
401inline G4UserSteppingAction* G4SteppingManager::GetUserAction() { return fUserSteppingAction; }
402
403inline G4Track* G4SteppingManager::GetTrack() const { return fTrack; }
404
405inline void G4SteppingManager::SetVerboseLevel(G4int vLevel) { verboseLevel = vLevel; }
406
408{
409 fVerbose = yourVerbose;
410}
411
412inline G4Step* G4SteppingManager::GetStep() const { return fStep; }
413
414inline G4double G4SteppingManager::CalculateSafety()
415{
416 return std::max(
417 endpointSafety - (endpointSafOrigin - fPostStepPoint->GetPosition()).mag(), kCarTolerance);
418}
419
420#endif
G4ForceCondition
@ InActivated
G4GPILSelection
@ NotCandidateForSelection
std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
std::vector< int, std::allocator< int > > G4SelectedAlongStepDoItVector
G4StepStatus
@ fUndefined
G4SteppingControl
@ NormalCondition
std::vector< G4int > G4SelectedPostStepDoItVector
std::vector< G4int > G4SelectedAtRestDoItVector
std::vector< G4int > G4SelectedAlongStepDoItVector
std::vector< G4Track * > G4TrackVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4ThreeVector & GetPosition() const
G4TrackVector * GetfSecondary()
G4ProfilerConfig< G4ProfileType::Step > ProfilerConfig
Definition G4Step.hh:65
const G4TrackVector * GetSecondary() const
G4ProcessVector * GetfPostStepDoItVector()
std::size_t GetfAtRestDoItProcTriggered()
std::size_t GetMAXofAtRestLoops()
G4TrackVector * GetfSecondary()
void SetVerbose(G4VSteppingVerbose *)
std::size_t GetMAXofAlongStepLoops()
G4SteppingControl GetStepControlFlag()
std::size_t GetMAXofPostStepLoops()
G4SelectedPostStepDoItVector * GetfSelectedPostStepDoItVector()
G4StepStatus Stepping()
G4Navigator * GetfNavigator()
G4Track * GetTrack() const
G4VPhysicalVolume * GetfCurrentVolume()
G4StepPoint * GetfPreStepPoint()
G4ProcessVector * GetfAlongStepDoItVector()
G4VParticleChange * GetfParticleChange()
std::size_t GetfPostStepDoItProcTriggered()
void SetNavigator(G4Navigator *value)
std::size_t GetfAlongStepDoItProcTriggered()
void SetVerboseLevel(G4int vLevel)
G4double GetfPreviousStepSize()
G4SelectedAlongStepDoItVector * GetfSelectedAlongStepDoItVector()
G4ForceCondition GetfCondition()
G4ProcessVector * GetfAlongStepGetPhysIntVector()
void SetInitialStep(G4Track *valueTrack)
G4VSensitiveDetector * GetfSensitive()
G4ProcessVector * GetfAtRestDoItVector()
G4StepStatus GetfStepStatus()
G4Step * GetStep() const
G4SelectedAtRestDoItVector * GetfSelectedAtRestDoItVector()
G4double GetnumberOfInteractionLengthLeft()
G4GPILSelection GetfGPILSelection()
G4ProcessVector * GetfAtRestGetPhysIntVector()
G4double GetcurrentMinimumStep()
G4ProcessVector * GetfPostStepGetPhysIntVector()
G4StepPoint * GetfPostStepPoint()
G4UserSteppingAction * GetUserAction()
void SetUserAction(G4UserSteppingAction *apAction)
G4VProcess * GetfCurrentProcess()
const G4TouchableHandle & GetTouchableHandle()
const G4TrackVector * GetSecondary() const