Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SteppingManager Class Reference

#include <G4SteppingManager.hh>

Public Types

using ProfilerConfig = G4Step::ProfilerConfig
 

Public Member Functions

 G4SteppingManager ()
 
 ~G4SteppingManager ()
 
const G4TrackVectorGetSecondary () const
 
void SetUserAction (G4UserSteppingAction *apAction)
 
G4TrackGetTrack () const
 
void SetVerboseLevel (G4int vLevel)
 
void SetVerbose (G4VSteppingVerbose *)
 
G4StepGetStep () const
 
void SetNavigator (G4Navigator *value)
 
G4StepStatus Stepping ()
 
void SetInitialStep (G4Track *valueTrack)
 
void GetProcessNumber ()
 
G4double GetPhysicalStep ()
 
G4double GetGeometricalStep ()
 
G4double GetCorrectedStep ()
 
G4bool GetPreStepPointIsGeom ()
 
G4bool GetFirstStep ()
 
G4StepStatus GetfStepStatus ()
 
G4double GetTempInitVelocity ()
 
G4double GetTempVelocity ()
 
G4double GetMass ()
 
G4double GetsumEnergyChange ()
 
G4VParticleChangeGetfParticleChange ()
 
G4TrackGetfTrack ()
 
G4TrackVectorGetfSecondary ()
 
G4StepGetfStep ()
 
G4StepPointGetfPreStepPoint ()
 
G4StepPointGetfPostStepPoint ()
 
G4VPhysicalVolumeGetfCurrentVolume ()
 
G4VSensitiveDetectorGetfSensitive ()
 
G4VProcessGetfCurrentProcess ()
 
G4ProcessVectorGetfAtRestDoItVector ()
 
G4ProcessVectorGetfAlongStepDoItVector ()
 
G4ProcessVectorGetfPostStepDoItVector ()
 
G4ProcessVectorGetfAlongStepGetPhysIntVector ()
 
G4ProcessVectorGetfPostStepGetPhysIntVector ()
 
G4ProcessVectorGetfAtRestGetPhysIntVector ()
 
G4double GetcurrentMinimumStep ()
 
G4double GetnumberOfInteractionLengthLeft ()
 
std::size_t GetfAtRestDoItProcTriggered ()
 
std::size_t GetfAlongStepDoItProcTriggered ()
 
std::size_t GetfPostStepDoItProcTriggered ()
 
G4int GetfN2ndariesAtRestDoIt ()
 
G4int GetfN2ndariesAlongStepDoIt ()
 
G4int GetfN2ndariesPostStepDoIt ()
 
G4NavigatorGetfNavigator ()
 
G4int GetverboseLevel ()
 
std::size_t GetMAXofAtRestLoops ()
 
std::size_t GetMAXofAlongStepLoops ()
 
std::size_t GetMAXofPostStepLoops ()
 
G4SelectedAtRestDoItVectorGetfSelectedAtRestDoItVector ()
 
G4SelectedAlongStepDoItVectorGetfSelectedAlongStepDoItVector ()
 
G4SelectedPostStepDoItVectorGetfSelectedPostStepDoItVector ()
 
G4double GetfPreviousStepSize ()
 
const G4TouchableHandleGetTouchableHandle ()
 
G4SteppingControl GetStepControlFlag ()
 
G4UserSteppingActionGetUserAction ()
 
G4double GetphysIntLength ()
 
G4ForceCondition GetfCondition ()
 
G4GPILSelection GetfGPILSelection ()
 

Detailed Description

Definition at line 74 of file G4SteppingManager.hh.

Member Typedef Documentation

◆ ProfilerConfig

Constructor & Destructor Documentation

◆ G4SteppingManager()

G4SteppingManager::G4SteppingManager ( )

Definition at line 48 of file G4SteppingManager.cc.

50{
51 // Construct simple 'has-a' related objects
52
53 fStep = new G4Step();
54 fSecondary = fStep->NewSecondaryVector();
55 fPreStepPoint = fStep->GetPreStepPoint();
56 fPostStepPoint = fStep->GetPostStepPoint();
57
58#ifdef G4VERBOSE
60 {
61 fVerbose = new G4SteppingVerbose();
63 fVerbose -> SetManager(this);
64 KillVerbose = true;
65 }
66 else
67 {
69 fVerbose -> SetManager(this);
70 KillVerbose = false;
71 }
72#endif
73
75 ->GetNavigatorForTracking());
76
77 fSelectedAtRestDoItVector
78 = new G4SelectedAtRestDoItVector(SizeOfSelectedDoItVector,0);
79 fSelectedAlongStepDoItVector
80 = new G4SelectedAlongStepDoItVector(SizeOfSelectedDoItVector,0);
81 fSelectedPostStepDoItVector
82 = new G4SelectedPostStepDoItVector(SizeOfSelectedDoItVector,0);
83
85 ->GetNavigatorForTracking());
86
87 physIntLength = DBL_MAX;
88 kCarTolerance = 0.5*G4GeometryTolerance::GetInstance()
90}
std::vector< G4int > G4SelectedPostStepDoItVector
std::vector< G4int > G4SelectedAtRestDoItVector
std::vector< G4int > G4SelectedAlongStepDoItVector
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4TrackVector * NewSecondaryVector()
G4StepPoint * GetPostStepPoint() const
void SetNavigator(G4Navigator *value)
static G4TransportationManager * GetTransportationManager()
static G4VSteppingVerbose * GetInstance()
static void SetInstance(G4VSteppingVerbose *Instance)
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4SteppingManager()

G4SteppingManager::~G4SteppingManager ( )

Definition at line 93 of file G4SteppingManager.cc.

95{
96 fTouchableHandle = 0;
97
98 // Destruct simple 'has-a' objects
99 //
100 fStep->DeleteSecondaryVector();
101
102 // delete fSecondary;
103 delete fStep;
104 delete fSelectedAtRestDoItVector;
105 delete fSelectedAlongStepDoItVector;
106 delete fSelectedPostStepDoItVector;
107 delete fUserSteppingAction;
108 #ifdef G4VERBOSE
109 if(KillVerbose) delete fVerbose;
110 #endif
111}
void DeleteSecondaryVector()

Member Function Documentation

◆ GetCorrectedStep()

G4double G4SteppingManager::GetCorrectedStep ( )
inline

Definition at line 287 of file G4SteppingManager.hh.

288 {
289 return CorrectedStep;
290 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetcurrentMinimumStep()

G4double G4SteppingManager::GetcurrentMinimumStep ( )

◆ GetfAlongStepDoItProcTriggered()

size_t G4SteppingManager::GetfAlongStepDoItProcTriggered ( )
inline

Definition at line 422 of file G4SteppingManager.hh.

423 {
424 return fAtRestDoItProcTriggered;
425 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAlongStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepDoItVector ( )
inline

Definition at line 377 of file G4SteppingManager.hh.

378 {
379 return fAlongStepDoItVector;
380 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAlongStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepGetPhysIntVector ( )
inline

Definition at line 392 of file G4SteppingManager.hh.

393 {
394 return fAlongStepGetPhysIntVector;
395 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAtRestDoItProcTriggered()

size_t G4SteppingManager::GetfAtRestDoItProcTriggered ( )
inline

Definition at line 417 of file G4SteppingManager.hh.

418 {
419 return fAtRestDoItProcTriggered;
420 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAtRestDoItVector()

G4ProcessVector * G4SteppingManager::GetfAtRestDoItVector ( )
inline

Definition at line 372 of file G4SteppingManager.hh.

373 {
374 return fAtRestDoItVector;
375 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAtRestGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAtRestGetPhysIntVector ( )
inline

Definition at line 387 of file G4SteppingManager.hh.

388 {
389 return fAtRestGetPhysIntVector;
390 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfCondition()

G4ForceCondition G4SteppingManager::GetfCondition ( )
inline

Definition at line 495 of file G4SteppingManager.hh.

496 {
497 return fCondition;
498 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfCurrentProcess()

G4VProcess * G4SteppingManager::GetfCurrentProcess ( )
inline

Definition at line 367 of file G4SteppingManager.hh.

368 {
369 return fCurrentProcess;
370 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfCurrentVolume()

G4VPhysicalVolume * G4SteppingManager::GetfCurrentVolume ( )
inline

Definition at line 357 of file G4SteppingManager.hh.

358 {
359 return fCurrentVolume;
360 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfGPILSelection()

G4GPILSelection G4SteppingManager::GetfGPILSelection ( )
inline

Definition at line 500 of file G4SteppingManager.hh.

501 {
502 return fGPILSelection;
503 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetFirstStep()

G4bool G4SteppingManager::GetFirstStep ( )
inline

Definition at line 297 of file G4SteppingManager.hh.

298 {
299 return FirstStep;
300 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfN2ndariesAlongStepDoIt()

G4int G4SteppingManager::GetfN2ndariesAlongStepDoIt ( )
inline

Definition at line 437 of file G4SteppingManager.hh.

438 {
439 return fN2ndariesAlongStepDoIt;
440 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfN2ndariesAtRestDoIt()

G4int G4SteppingManager::GetfN2ndariesAtRestDoIt ( )
inline

Definition at line 432 of file G4SteppingManager.hh.

433 {
434 return fN2ndariesAtRestDoIt;
435 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfN2ndariesPostStepDoIt()

G4int G4SteppingManager::GetfN2ndariesPostStepDoIt ( )
inline

Definition at line 442 of file G4SteppingManager.hh.

443 {
444 return fN2ndariesPostStepDoIt;
445 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfNavigator()

G4Navigator * G4SteppingManager::GetfNavigator ( )
inline

Definition at line 447 of file G4SteppingManager.hh.

448 {
449 return fNavigator;
450 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfParticleChange()

G4VParticleChange * G4SteppingManager::GetfParticleChange ( )
inline

Definition at line 327 of file G4SteppingManager.hh.

328 {
329 return fParticleChange;
330 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepDoItProcTriggered()

size_t G4SteppingManager::GetfPostStepDoItProcTriggered ( )
inline

Definition at line 427 of file G4SteppingManager.hh.

428 {
429 return fPostStepDoItProcTriggered;
430 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfPostStepDoItVector ( )
inline

Definition at line 382 of file G4SteppingManager.hh.

383 {
384 return fPostStepDoItVector;
385 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfPostStepGetPhysIntVector ( )
inline

Definition at line 397 of file G4SteppingManager.hh.

398 {
399 return fPostStepGetPhysIntVector;
400 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepPoint()

G4StepPoint * G4SteppingManager::GetfPostStepPoint ( )
inline

Definition at line 352 of file G4SteppingManager.hh.

353 {
354 return fPostStepPoint;
355 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPreStepPoint()

G4StepPoint * G4SteppingManager::GetfPreStepPoint ( )
inline

Definition at line 347 of file G4SteppingManager.hh.

348 {
349 return fPreStepPoint;
350 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPreviousStepSize()

G4double G4SteppingManager::GetfPreviousStepSize ( )
inline

Definition at line 475 of file G4SteppingManager.hh.

476 {
477 return fPreviousStepSize;
478 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSecondary()

G4TrackVector * G4SteppingManager::GetfSecondary ( )
inline

Definition at line 337 of file G4SteppingManager.hh.

338 {
339 return fStep->GetfSecondary();
340 }
G4TrackVector * GetfSecondary()

Referenced by G4VSteppingVerbose::CopyState(), and G4TrackingManager::GimmeSecondaries().

◆ GetfSelectedAlongStepDoItVector()

G4SelectedAlongStepDoItVector * G4SteppingManager::GetfSelectedAlongStepDoItVector ( )
inline

Definition at line 464 of file G4SteppingManager.hh.

465 {
466 return fSelectedAlongStepDoItVector;
467 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSelectedAtRestDoItVector()

G4SelectedAtRestDoItVector * G4SteppingManager::GetfSelectedAtRestDoItVector ( )
inline

Definition at line 458 of file G4SteppingManager.hh.

459 {
460 return fSelectedAtRestDoItVector;
461 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSelectedPostStepDoItVector()

G4SelectedPostStepDoItVector * G4SteppingManager::GetfSelectedPostStepDoItVector ( )
inline

Definition at line 470 of file G4SteppingManager.hh.

471 {
472 return fSelectedPostStepDoItVector;
473 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSensitive()

G4VSensitiveDetector * G4SteppingManager::GetfSensitive ( )
inline

Definition at line 362 of file G4SteppingManager.hh.

363 {
364 return fSensitive;
365 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfStep()

G4Step * G4SteppingManager::GetfStep ( )
inline

Definition at line 342 of file G4SteppingManager.hh.

343 {
344 return fStep;
345 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfStepStatus()

G4StepStatus G4SteppingManager::GetfStepStatus ( )
inline

Definition at line 302 of file G4SteppingManager.hh.

303 {
304 return fStepStatus;
305 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfTrack()

G4Track * G4SteppingManager::GetfTrack ( )
inline

Definition at line 332 of file G4SteppingManager.hh.

333 {
334 return fTrack;
335 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetGeometricalStep()

G4double G4SteppingManager::GetGeometricalStep ( )
inline

Definition at line 282 of file G4SteppingManager.hh.

283 {
284 return GeometricalStep;
285 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMass()

G4double G4SteppingManager::GetMass ( )
inline

Definition at line 317 of file G4SteppingManager.hh.

318 {
319 return Mass;
320 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMAXofAlongStepLoops()

size_t G4SteppingManager::GetMAXofAlongStepLoops ( )
inline

Definition at line 407 of file G4SteppingManager.hh.

408 {
409 return MAXofAlongStepLoops;
410 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMAXofAtRestLoops()

size_t G4SteppingManager::GetMAXofAtRestLoops ( )
inline

Definition at line 402 of file G4SteppingManager.hh.

403 {
404 return MAXofAtRestLoops;
405 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMAXofPostStepLoops()

size_t G4SteppingManager::GetMAXofPostStepLoops ( )
inline

Definition at line 412 of file G4SteppingManager.hh.

413 {
414 return MAXofPostStepLoops;
415 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetnumberOfInteractionLengthLeft()

G4double G4SteppingManager::GetnumberOfInteractionLengthLeft ( )

◆ GetPhysicalStep()

G4double G4SteppingManager::GetPhysicalStep ( )
inline

Definition at line 277 of file G4SteppingManager.hh.

278 {
279 return PhysicalStep;
280 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetphysIntLength()

G4double G4SteppingManager::GetphysIntLength ( )
inline

Definition at line 490 of file G4SteppingManager.hh.

491 {
492 return physIntLength;
493 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetPreStepPointIsGeom()

G4bool G4SteppingManager::GetPreStepPointIsGeom ( )
inline

Definition at line 292 of file G4SteppingManager.hh.

293 {
294 return PreStepPointIsGeom;
295 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetProcessNumber()

void G4SteppingManager::GetProcessNumber ( )

Definition at line 50 of file G4SteppingManager2.cc.

52{
53 #ifdef debug
54 G4cout << "G4SteppingManager::GetProcessNumber: is called track="
55 << fTrack << G4endl;
56 #endif
57
59 if(pm == nullptr)
60 {
61 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
62 << " ProcessManager is NULL for particle = "
63 << fTrack->GetDefinition()->GetParticleName() << ", PDG_code = "
64 << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
65 G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011",
66 FatalException, "Process Manager is not found.");
67 return;
68 }
69
70 // AtRestDoits
71 //
72 MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
73 fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
74 fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
75
76 #ifdef debug
77 G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest="
78 << MAXofAtRestLoops << G4endl;
79 #endif
80
81 // AlongStepDoits
82 //
83 MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
84 fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
85 fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
86
87 #ifdef debug
88 G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp="
89 << MAXofAlongStepLoops << G4endl;
90 #endif
91
92 // PostStepDoits
93 //
94 MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
95 fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
96 fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
97
98 #ifdef debug
99 G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep="
100 << MAXofPostStepLoops << G4endl;
101 #endif
102
103 if (SizeOfSelectedDoItVector<MAXofAtRestLoops ||
104 SizeOfSelectedDoItVector<MAXofAlongStepLoops ||
105 SizeOfSelectedDoItVector<MAXofPostStepLoops )
106 {
107 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
108 << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
109 << " ; is smaller then one of MAXofAtRestLoops= "
110 << MAXofAtRestLoops << G4endl
111 << " or MAXofAlongStepLoops= " << MAXofAlongStepLoops
112 << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
113 G4Exception("G4SteppingManager::GetProcessNumber()",
114 "Tracking0012", FatalException,
115 "The array size is smaller than the actual No of processes.");
116 }
117}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ typeGPIL
@ typeDoIt
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t entries() const
G4ParticleDefinition * GetDefinition() const

Referenced by G4ErrorPropagator::InitG4Track(), and G4TrackingManager::ProcessOneTrack().

◆ GetSecondary()

const G4TrackVector * G4SteppingManager::GetSecondary ( ) const
inline

Definition at line 505 of file G4SteppingManager.hh.

506 {
507 return fStep->GetSecondary();
508 }
const G4TrackVector * GetSecondary() const

◆ GetStep()

G4Step * G4SteppingManager::GetStep ( ) const
inline

Definition at line 540 of file G4SteppingManager.hh.

541 {
542 return fStep;
543 }

Referenced by G4ErrorPropagator::InitG4Track(), and G4TrackingManager::ProcessOneTrack().

◆ GetStepControlFlag()

G4SteppingControl G4SteppingManager::GetStepControlFlag ( )
inline

Definition at line 485 of file G4SteppingManager.hh.

486 {
487 return StepControlFlag;
488 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetsumEnergyChange()

G4double G4SteppingManager::GetsumEnergyChange ( )
inline

Definition at line 322 of file G4SteppingManager.hh.

323 {
324 return sumEnergyChange;
325 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTempInitVelocity()

G4double G4SteppingManager::GetTempInitVelocity ( )
inline

Definition at line 307 of file G4SteppingManager.hh.

308 {
309 return TempInitVelocity;
310 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTempVelocity()

G4double G4SteppingManager::GetTempVelocity ( )
inline

Definition at line 312 of file G4SteppingManager.hh.

313 {
314 return TempVelocity;
315 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTouchableHandle()

const G4TouchableHandle & G4SteppingManager::GetTouchableHandle ( )
inline

Definition at line 480 of file G4SteppingManager.hh.

481 {
482 return fTouchableHandle;
483 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTrack()

G4Track * G4SteppingManager::GetTrack ( ) const
inline

Definition at line 525 of file G4SteppingManager.hh.

526 {
527 return fTrack;
528 }

Referenced by G4TrackingMessenger::SetNewValue().

◆ GetUserAction()

G4UserSteppingAction * G4SteppingManager::GetUserAction ( )
inline

Definition at line 520 of file G4SteppingManager.hh.

521 {
522 return fUserSteppingAction;
523 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetverboseLevel()

G4int G4SteppingManager::GetverboseLevel ( )
inline

Definition at line 452 of file G4SteppingManager.hh.

453 {
454 return verboseLevel;
455 }

Referenced by G4VSteppingVerbose::CopyState().

◆ SetInitialStep()

void G4SteppingManager::SetInitialStep ( G4Track valueTrack)

Definition at line 272 of file G4SteppingManager.cc.

274{
275 // Set up several local variables
276 //
277 PreStepPointIsGeom = false;
278 FirstStep = true;
279 fParticleChange = nullptr;
280 fPreviousStepSize = 0.;
281 fStepStatus = fUndefined;
282
283 fTrack = valueTrack;
284 Mass = fTrack->GetDynamicParticle()->GetMass();
285
286 PhysicalStep = 0.;
287 GeometricalStep = 0.;
288 CorrectedStep = 0.;
289 PreStepPointIsGeom = false;
290 FirstStep = false;
291
292 TempInitVelocity = 0.;
293 TempVelocity = 0.;
294 sumEnergyChange = 0.;
295
296 // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
297 // set the track state to 'Alive'
298 //
299 if ( (fTrack->GetTrackStatus()==fSuspend)
300 || (fTrack->GetTrackStatus()==fPostponeToNextEvent) )
301 {
302 fTrack->SetTrackStatus(fAlive);
303 }
304
305 // If the primary track has 'zero' kinetic energy, set the track
306 // state to 'StopButAlive'
307 //
308 if(fTrack->GetKineticEnergy() <= 0.0)
309 {
310 fTrack->SetTrackStatus( fStopButAlive );
311 }
312
313 // Set Touchable to track and a private attribute of G4SteppingManager
314
315 if ( ! fTrack->GetTouchableHandle() )
316 {
317 G4ThreeVector direction= fTrack->GetMomentumDirection();
318 fNavigator->LocateGlobalPointAndSetup( fTrack->GetPosition(),
319 &direction, false, false );
320 fTouchableHandle = fNavigator->CreateTouchableHistory();
321 fTrack->SetTouchableHandle( fTouchableHandle );
322 fTrack->SetNextTouchableHandle( fTouchableHandle );
323 }
324 else
325 {
326 fTrack->SetNextTouchableHandle( fTouchableHandle = fTrack->GetTouchableHandle() );
327 G4VPhysicalVolume* oldTopVolume = fTrack->GetTouchableHandle()->GetVolume();
328 G4VPhysicalVolume* newTopVolume =
329 fNavigator->ResetHierarchyAndLocate( fTrack->GetPosition(),
330 fTrack->GetMomentumDirection(),
331 *((G4TouchableHistory*)fTrack->GetTouchableHandle()()) );
332 if ( newTopVolume != oldTopVolume
333 || oldTopVolume->GetRegularStructureId() == 1 )
334 {
335 fTouchableHandle = fNavigator->CreateTouchableHistory();
336 fTrack->SetTouchableHandle( fTouchableHandle );
337 fTrack->SetNextTouchableHandle( fTouchableHandle );
338 }
339 }
340
341 // Set OriginTouchableHandle for primary track
342 //
343 if(fTrack->GetParentID()==0)
344 {
346 }
347
348 // Set vertex information of G4Track at here
349 //
350 if ( fTrack->GetCurrentStepNumber() == 0 )
351 {
352 fTrack->SetVertexPosition( fTrack->GetPosition() );
354 fTrack->SetVertexKineticEnergy( fTrack->GetKineticEnergy() );
356 }
357
358 // Initial set up for attributes of 'G4SteppingManager'
359 fCurrentVolume = fTouchableHandle->GetVolume();
360
361 // If track is already outside the world boundary, kill it
362 //
363 if( fCurrentVolume==nullptr )
364 {
365 // If the track is a primary, stop processing
366 if(fTrack->GetParentID()==0)
367 {
368 G4cerr << "ERROR - G4SteppingManager::SetInitialStep()" << G4endl
369 << " Primary particle starting at - "
370 << fTrack->GetPosition()
371 << " - is outside of the world volume." << G4endl;
372 G4Exception("G4SteppingManager::SetInitialStep()", "Tracking0010",
373 FatalException, "Primary vertex outside of the world!");
374 }
375
376 fTrack->SetTrackStatus( fStopAndKill );
377 G4cout << "WARNING - G4SteppingManager::SetInitialStep()" << G4endl
378 << " Initial track position is outside world! - "
379 << fTrack->GetPosition() << G4endl;
380 }
381 else
382 {
383 // Initial set up for attributes of 'Step'
384 fStep->InitializeStep( fTrack );
385 }
386
387 #ifdef G4VERBOSE
388 if(verboseLevel>0) fVerbose->TrackingStarted();
389 #endif
390}
@ fUndefined
Definition: G4StepStatus.hh:55
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
G4double GetMass() const
G4TouchableHistory * CreateTouchableHistory() const
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=nullptr, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:126
virtual G4VPhysicalVolume * ResetHierarchyAndLocate(const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
Definition: G4Navigator.cc:96
void InitializeStep(G4Track *aValue)
G4TrackStatus GetTrackStatus() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4VPhysicalVolume * GetVolume() const
void SetVertexPosition(const G4ThreeVector &aValue)
void SetVertexMomentumDirection(const G4ThreeVector &aValue)
void SetNextTouchableHandle(const G4TouchableHandle &apValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4int GetCurrentStepNumber() const
void SetOriginTouchableHandle(const G4TouchableHandle &apValue)
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetVertexKineticEnergy(const G4double aValue)
G4int GetParentID() const
void SetLogicalVolumeAtVertex(const G4LogicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetRegularStructureId() const =0
virtual void TrackingStarted()=0
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:41

Referenced by G4ErrorPropagator::InitG4Track(), and G4TrackingManager::ProcessOneTrack().

◆ SetNavigator()

void G4SteppingManager::SetNavigator ( G4Navigator value)
inline

Definition at line 510 of file G4SteppingManager.hh.

511 {
512 fNavigator = value;
513 }

Referenced by G4SteppingManager().

◆ SetUserAction()

void G4SteppingManager::SetUserAction ( G4UserSteppingAction apAction)
inline

Definition at line 515 of file G4SteppingManager.hh.

516 {
517 fUserSteppingAction = apAction;
518 }

Referenced by G4TrackingManager::SetUserAction().

◆ SetVerbose()

void G4SteppingManager::SetVerbose ( G4VSteppingVerbose yourVerbose)
inline

Definition at line 535 of file G4SteppingManager.hh.

536 {
537 fVerbose = yourVerbose;
538 }

◆ SetVerboseLevel()

void G4SteppingManager::SetVerboseLevel ( G4int  vLevel)
inline

Definition at line 530 of file G4SteppingManager.hh.

531 {
532 verboseLevel = vLevel;
533 }

Referenced by G4ErrorPropagatorManager::SetSteppingManagerVerboseLevel().

◆ Stepping()

G4StepStatus G4SteppingManager::Stepping ( )

Definition at line 114 of file G4SteppingManager.cc.

116{
117#ifdef GEANT4_USE_TIMEMORY
118 ProfilerConfig profiler{ fStep };
119#endif
120
121 //--------
122 // Prelude
123 //--------
124 #ifdef G4VERBOSE
125 if(verboseLevel>0)
126 {
127 fVerbose->NewStep();
128 }
129 else if (verboseLevel==-1)
130 {
132 }
133 else
134 {
136 }
137 #endif
138
139 // Store last PostStepPoint to PreStepPoint, and swap current and nex
140 // volume information of G4Track. Reset total energy deposit in one Step.
141 //
142 fStep->CopyPostToPreStepPoint();
144
145 // Switch next touchable in track to current one
146 //
148
149 // Reset the secondary particles
150 //
151 fN2ndariesAtRestDoIt = 0;
152 fN2ndariesAlongStepDoIt = 0;
153 fN2ndariesPostStepDoIt = 0;
154
155 // Set the volume before it is used (in DefineStepLength() for User Limit)
156 //
157 fCurrentVolume = fStep->GetPreStepPoint()->GetPhysicalVolume();
158
159 // Reset the step's auxiliary points vector pointer
160 //
162
163 //-----------------
164 // AtRest Processes
165 //-----------------
166
167 if( fTrack->GetTrackStatus() == fStopButAlive )
168 {
169 if( MAXofAtRestLoops>0 )
170 {
171 InvokeAtRestDoItProcs();
172 fStepStatus = fAtRestDoItProc;
173 fStep->GetPostStepPoint()->SetStepStatus( fStepStatus );
174
175 #ifdef G4VERBOSE
176 if(verboseLevel>0) fVerbose->AtRestDoItInvoked();
177 #endif
178
179 }
180 // Make sure the track is killed
181 //
182 fTrack->SetTrackStatus( fStopAndKill );
183 }
184
185 //---------------------------------
186 // AlongStep and PostStep Processes
187 //---------------------------------
188
189 else
190 {
191 // Find minimum Step length demanded by active disc./cont. processes
192 DefinePhysicalStepLength();
193
194 // Store the Step length (geometrical length) to G4Step and G4Track
195 fStep->SetStepLength( PhysicalStep );
196 fTrack->SetStepLength( PhysicalStep );
197 G4double GeomStepLength = PhysicalStep;
198
199 // Store StepStatus to PostStepPoint
200 fStep->GetPostStepPoint()->SetStepStatus( fStepStatus );
201
202 // Invoke AlongStepDoIt
203 InvokeAlongStepDoItProcs();
204
205 // Update track by taking into account all changes by AlongStepDoIt
206 fStep->UpdateTrack();
207
208 // Update safety after invocation of all AlongStepDoIts
209 endpointSafOrigin= fPostStepPoint->GetPosition();
210 // endpointSafety= std::max( proposedSafety - GeomStepLength, 0.);
211 endpointSafety= std::max( proposedSafety - GeomStepLength, kCarTolerance);
212
213 fStep->GetPostStepPoint()->SetSafety( endpointSafety );
214
215 #ifdef G4VERBOSE
216 if(verboseLevel>0) fVerbose->AlongStepDoItAllDone();
217 #endif
218
219 // Invoke PostStepDoIt
220 InvokePostStepDoItProcs();
221
222 #ifdef G4VERBOSE
223 if(verboseLevel>0) fVerbose->PostStepDoItAllDone();
224 #endif
225 }
226
227 //-------
228 // Finale
229 //-------
230
231 // Update 'TrackLength' and remeber the Step length of the current Step
232 //
233 fTrack->AddTrackLength(fStep->GetStepLength());
234 fPreviousStepSize = fStep->GetStepLength();
235 fStep->SetTrack(fTrack);
236
237 #ifdef G4VERBOSE
238 if(verboseLevel>0) fVerbose->StepInfo();
239 #endif
240
241 // Send G4Step information to Hit/Dig if the volume is sensitive
242 //
243 fCurrentVolume = fStep->GetPreStepPoint()->GetPhysicalVolume();
244 StepControlFlag = fStep->GetControlFlag();
245 if( fCurrentVolume != nullptr && StepControlFlag != AvoidHitInvocation )
246 {
247 fSensitive = fStep->GetPreStepPoint()->GetSensitiveDetector();
248 if( fSensitive != 0 )
249 {
250 fSensitive->Hit(fStep);
251 }
252 }
253
254 // User intervention process
255 //
256 if( fUserSteppingAction != nullptr )
257 {
258 fUserSteppingAction->UserSteppingAction(fStep);
259 }
260 G4UserSteppingAction* regionalAction = fStep->GetPreStepPoint()
263 if(regionalAction)
264 regionalAction->UserSteppingAction(fStep);
265
266 // Stepping process finish. Return the value of the StepStatus
267 //
268 return fStepStatus;
269}
@ fAtRestDoItProc
Definition: G4StepStatus.hh:45
@ AvoidHitInvocation
double G4double
Definition: G4Types.hh:83
G4Region * GetRegion() const
G4UserSteppingAction * GetRegionalSteppingAction() const
Definition: G4Region.cc:148
void SetSafety(const G4double aValue)
void SetStepStatus(const G4StepStatus aValue)
const G4ThreeVector & GetPosition() const
G4VSensitiveDetector * GetSensitiveDetector() const
G4VPhysicalVolume * GetPhysicalVolume() const
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *vec)
G4SteppingControl GetControlFlag() const
void UpdateTrack()
void ResetTotalEnergyDeposit()
void SetStepLength(G4double value)
void CopyPostToPreStepPoint()
G4double GetStepLength() const
void SetTrack(G4Track *value)
G4Step::ProfilerConfig ProfilerConfig
void SetStepLength(G4double value)
const G4TouchableHandle & GetNextTouchableHandle() const
void AddTrackLength(const G4double aValue)
virtual void UserSteppingAction(const G4Step *)
G4bool Hit(G4Step *aStep)
virtual void AlongStepDoItAllDone()=0
virtual void PostStepDoItAllDone()=0
virtual void AtRestDoItInvoked()=0
virtual void StepInfo()=0
static void SetSilent(G4int fSilent)
virtual void NewStep()=0

Referenced by G4ErrorPropagator::MakeOneStep(), and G4TrackingManager::ProcessOneTrack().


The documentation for this class was generated from the following files: