Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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: Katsuya.Amako@kek.jp)
39// Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
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 <iomanip> // Include from 'system'
48#include <vector> // Include from 'system'
49#include "globals.hh" // Include from 'global'
50#include "Randomize.hh" // Include from 'global'
51
52#include "G4Navigator.hh" // Include from 'geometry'
53#include "G4LogicalVolume.hh" // Include from 'geometry'
54#include "G4VPhysicalVolume.hh" // Include from 'geometry'
55#include "G4ProcessManager.hh" // Include from 'processes'
56#include "G4NoProcess.hh" // Include from 'processes'
57
58#include "G4Track.hh" // Include from 'tracking'
59#include "G4TrackVector.hh" // Include from 'tracking'
60#include "G4TrackStatus.hh" // Include from 'tracking'
61#include "G4StepStatus.hh" // Include from 'tracking'
62#include "G4UserSteppingAction.hh" // Include from 'tracking'
63#include "G4Step.hh" // Include from 'tracking'
64#include "G4StepPoint.hh" // Include from 'tracking'
65#include "G4VSteppingVerbose.hh" // Include from 'tracking'
66#include "G4TouchableHandle.hh" // Include from 'geometry'
67#include "G4TouchableHistoryHandle.hh" // Include from 'geometry'
68
69using G4SelectedAtRestDoItVector = std::vector<G4int>;
70using G4SelectedAlongStepDoItVector = std::vector<G4int>;
71using G4SelectedPostStepDoItVector = std::vector<G4int>;
72
74
76{
77 public:
79
80 public:
81 // Constructor/Destructor
82
84 // SteppingManger should be dynamically allocated, therefore
85 // you need to invoke new() when you call this constructor.
86 // "Secodary track vector" will be dynamically created by this
87 // constructor. G4UserSteppingAction will be also created
88 // in this constructor, and "this" pointer will be passed to
89 // G4UserSteppingAction.
90
92
93 // Get/Set functions
94
95 const G4TrackVector* GetSecondary() const;
96 void SetUserAction(G4UserSteppingAction* apAction);
97 G4Track* GetTrack() const;
98 void SetVerboseLevel(G4int vLevel);
100 G4Step* GetStep() const;
101 void SetNavigator(G4Navigator* value);
102
103 // Other member functions
104
106 // Steers to move the give particle from the TrackingManger by one Step.
107
108 void SetInitialStep(G4Track* valueTrack);
109 // Sets up initial track information (enegry, position, etc) to
110 // the PreStepPoint of the G4Step. This funciton has to be called
111 // just once before the stepping loop in the "TrackingManager".
112
113 void GetProcessNumber();
114
115 // Get methods
116
130 G4Step* GetfStep();
144 std::size_t GetfAtRestDoItProcTriggered();
145 std::size_t GetfAlongStepDoItProcTriggered();
146 std::size_t GetfPostStepDoItProcTriggered();
152 std::size_t GetMAXofAtRestLoops();
153 std::size_t GetMAXofAlongStepLoops();
154 std::size_t GetMAXofPostStepLoops();
165
166 private:
167
168 // Member functions
169
170 void DefinePhysicalStepLength();
171 // Calculate corresponding physical length from the mean free path
172 // left for each discrete physics process. The minimum allowable
173 // step for each continuous process will be also calculated.
174 void InvokeAtRestDoItProcs();
175 void InvokeAlongStepDoItProcs();
176 void InvokePostStepDoItProcs();
177 void InvokePSDIP(size_t); //
178 G4int ProcessSecondariesFromParticleChange();
179 G4double CalculateSafety();
180 // Return the estimated safety value at the PostStepPoint
181
182 // Member data
183
184 static const size_t SizeOfSelectedDoItVector = 100;
185
186 G4bool KillVerbose = false;
187
188 G4UserSteppingAction* fUserSteppingAction = nullptr;
189
190 G4VSteppingVerbose* fVerbose = nullptr;
191
192 G4double PhysicalStep = 0.0;
193 G4double GeometricalStep = 0.0;
194 G4double CorrectedStep = 0.0;
195 G4bool PreStepPointIsGeom = false;
196 G4bool FirstStep = false;
197 G4StepStatus fStepStatus = fUndefined;
198
199 G4double TempInitVelocity = 0.0;
200 G4double TempVelocity = 0.0;
201 G4double Mass = 0.0;
202
203 G4double sumEnergyChange = 0.0;
204
205 G4VParticleChange* fParticleChange = nullptr;
206 G4Track* fTrack = nullptr;
207 G4TrackVector* fSecondary = nullptr;
208 G4Step* fStep = nullptr;
209 G4StepPoint* fPreStepPoint = nullptr;
210 G4StepPoint* fPostStepPoint = nullptr;
211
212 G4VPhysicalVolume* fCurrentVolume = nullptr;
213 G4VSensitiveDetector* fSensitive = nullptr;
214 G4VProcess* fCurrentProcess = nullptr;
215 // The pointer to the process of which DoIt() or
216 // GetPhysicalInteractionLength() has been just executed.
217
218
219 G4ProcessVector* fAtRestDoItVector = nullptr;
220 G4ProcessVector* fAlongStepDoItVector = nullptr;
221 G4ProcessVector* fPostStepDoItVector = nullptr;
222
223 G4ProcessVector* fAtRestGetPhysIntVector = nullptr;
224 G4ProcessVector* fAlongStepGetPhysIntVector = nullptr;
225 G4ProcessVector* fPostStepGetPhysIntVector = nullptr;
226
227 std::size_t MAXofAtRestLoops = 0;
228 std::size_t MAXofAlongStepLoops = 0;
229 std::size_t MAXofPostStepLoops = 0;
230
231 std::size_t fAtRestDoItProcTriggered = 0;
232 std::size_t fAlongStepDoItProcTriggered = 0;
233 std::size_t fPostStepDoItProcTriggered = 0;
234
235 G4int fN2ndariesAtRestDoIt = 0;
236 G4int fN2ndariesAlongStepDoIt = 0;
237 G4int fN2ndariesPostStepDoIt = 0;
238 // These are the numbers of secondaries generated by the process
239 // just executed.
240
241 G4Navigator* fNavigator = nullptr;
242
243 G4int verboseLevel = 0;
244
245 G4SelectedAtRestDoItVector* fSelectedAtRestDoItVector = nullptr;
246 G4SelectedAlongStepDoItVector* fSelectedAlongStepDoItVector = nullptr;
247 G4SelectedPostStepDoItVector* fSelectedPostStepDoItVector = nullptr;
248
249 G4double fPreviousStepSize = 0.0;
250
251 G4TouchableHandle fTouchableHandle;
252
253 G4SteppingControl StepControlFlag = NormalCondition;
254
255 G4double kCarTolerance = 0.0;
256 // Cached geometrical tolerance on surface
257 G4double proposedSafety = 0.0;
258 // This keeps the minimum safety value proposed by AlongStepGPILs.
259 G4ThreeVector endpointSafOrigin;
260 G4double endpointSafety = 0.0;
261 // To get the true safety value at the PostStepPoint, you have
262 // to subtract the distance to 'endpointSafOrigin' from this value.
263 G4double physIntLength = 0.0;
264 G4ForceCondition fCondition = InActivated;
266 // Above three variables are for the method
267 // DefinePhysicalStepLength(). To pass these information to
268 // the method Verbose, they are kept at here. Need a more
269 // elegant mechanism.
270
271 G4NoProcess const* fNoProcess = nullptr;
272 // Used in the InvokeAtRestDoItProcs() method to flag the process
273 // of any stable ion at rest.
274};
275
276//*******************************************************************
277//
278// Inline functions
279//
280//*******************************************************************
281
283 {
284 return PhysicalStep;
285 }
286
288 {
289 return GeometricalStep;
290 }
291
293 {
294 return CorrectedStep;
295 }
296
298 {
299 return PreStepPointIsGeom;
300 }
301
303 {
304 return FirstStep;
305 }
306
308 {
309 return fStepStatus;
310 }
311
313 {
314 return TempInitVelocity;
315 }
316
318 {
319 return TempVelocity;
320 }
321
323 {
324 return Mass;
325 }
326
328 {
329 return sumEnergyChange;
330 }
331
333 {
334 return fParticleChange;
335 }
336
338 {
339 return fTrack;
340 }
341
343 {
344 return fStep->GetfSecondary();
345 }
346
348 {
349 return fStep;
350 }
351
353 {
354 return fPreStepPoint;
355 }
356
358 {
359 return fPostStepPoint;
360 }
361
363 {
364 return fCurrentVolume;
365 }
366
368 {
369 return fSensitive;
370 }
371
373 {
374 return fCurrentProcess;
375 }
376
378 {
379 return fAtRestDoItVector;
380 }
381
383 {
384 return fAlongStepDoItVector;
385 }
386
388 {
389 return fPostStepDoItVector;
390 }
391
393 {
394 return fAtRestGetPhysIntVector;
395 }
396
398 {
399 return fAlongStepGetPhysIntVector;
400 }
401
403 {
404 return fPostStepGetPhysIntVector;
405 }
406
408 {
409 return MAXofAtRestLoops;
410 }
411
413 {
414 return MAXofAlongStepLoops;
415 }
416
418 {
419 return MAXofPostStepLoops;
420 }
421
423 {
424 return fAtRestDoItProcTriggered;
425 }
426
428 {
429 return fAtRestDoItProcTriggered;
430 }
431
433 {
434 return fPostStepDoItProcTriggered;
435 }
436
438 {
439 return fN2ndariesAtRestDoIt;
440 }
441
443 {
444 return fN2ndariesAlongStepDoIt;
445 }
446
448 {
449 return fN2ndariesPostStepDoIt;
450 }
451
453 {
454 return fNavigator;
455 }
456
458 {
459 return verboseLevel;
460 }
461
464 {
465 return fSelectedAtRestDoItVector;
466 }
467
470 {
471 return fSelectedAlongStepDoItVector;
472 }
473
476 {
477 return fSelectedPostStepDoItVector;
478 }
479
481 {
482 return fPreviousStepSize;
483 }
484
486 {
487 return fTouchableHandle;
488 }
489
491 {
492 return StepControlFlag;
493 }
494
496 {
497 return physIntLength;
498 }
499
501 {
502 return fCondition;
503 }
504
506 {
507 return fGPILSelection;
508 }
509
511 {
512 return fStep->GetSecondary();
513 }
514
516 {
517 fNavigator = value;
518 }
519
521 {
522 fUserSteppingAction = apAction;
523 }
524
526 {
527 return fUserSteppingAction;
528 }
529
531 {
532 return fTrack;
533 }
534
536 {
537 verboseLevel = vLevel;
538 }
539
541 {
542 fVerbose = yourVerbose;
543 }
544
546 {
547 return fStep;
548 }
549
550 inline G4double G4SteppingManager::CalculateSafety()
551 {
552 return std::max( endpointSafety -
553 (endpointSafOrigin - fPostStepPoint->GetPosition()).mag(),
554 kCarTolerance );
555 }
556
557#endif
G4ForceCondition
@ InActivated
G4GPILSelection
@ NotCandidateForSelection
class std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
class std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
class std::vector< int, std::allocator< int > > G4SelectedAlongStepDoItVector
G4StepStatus
Definition: G4StepStatus.hh:40
@ fUndefined
Definition: G4StepStatus.hh:55
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
Definition: G4Step.hh:62
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)
G4double GetsumEnergyChange()
G4double GetTempInitVelocity()
G4VSensitiveDetector * GetfSensitive()
G4ProcessVector * GetfAtRestDoItVector()
G4StepStatus GetfStepStatus()
G4Step * GetStep() const
G4SelectedAtRestDoItVector * GetfSelectedAtRestDoItVector()
G4int GetfN2ndariesAlongStepDoIt()
G4double GetnumberOfInteractionLengthLeft()
G4GPILSelection GetfGPILSelection()
G4ProcessVector * GetfAtRestGetPhysIntVector()
G4double GetcurrentMinimumStep()
G4ProcessVector * GetfPostStepGetPhysIntVector()
G4StepPoint * GetfPostStepPoint()
G4UserSteppingAction * GetUserAction()
void SetUserAction(G4UserSteppingAction *apAction)
G4VProcess * GetfCurrentProcess()
const G4TouchableHandle & GetTouchableHandle()
G4double GetGeometricalStep()
const G4TrackVector * GetSecondary() const