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