Geant4 11.1.1
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 75 of file G4SteppingManager.hh.

Member Typedef Documentation

◆ ProfilerConfig

Constructor & Destructor Documentation

◆ G4SteppingManager()

G4SteppingManager::G4SteppingManager ( )

Definition at line 52 of file G4SteppingManager.cc.

54{
55 // Construct simple 'has-a' related objects
56
57 fStep = new G4Step();
58 fSecondary = fStep->NewSecondaryVector();
59 fPreStepPoint = fStep->GetPreStepPoint();
60 fPostStepPoint = fStep->GetPostStepPoint();
61
62#ifdef G4VERBOSE
64 if(fVerbose==nullptr)
65 {
67 {
69 if(prec > 0)
70 { fVerbose = new G4SteppingVerboseWithUnits(prec); }
71 else
72 { fVerbose = new G4SteppingVerbose(); }
73 }
74 else
76 KillVerbose = true;
77 }
78 else
79 { KillVerbose = false; }
80 fVerbose -> SetManager(this);
81#endif
82
84 ->GetNavigatorForTracking());
85
86 fSelectedAtRestDoItVector
87 = new G4SelectedAtRestDoItVector(SizeOfSelectedDoItVector,0);
88 fSelectedAlongStepDoItVector
89 = new G4SelectedAlongStepDoItVector(SizeOfSelectedDoItVector,0);
90 fSelectedPostStepDoItVector
91 = new G4SelectedPostStepDoItVector(SizeOfSelectedDoItVector,0);
92
94 ->GetNavigatorForTracking());
95
96 physIntLength = DBL_MAX;
97 kCarTolerance = 0.5*G4GeometryTolerance::GetInstance()
99
100 fNoProcess = new G4NoProcess;
101}
std::vector< G4int > G4SelectedPostStepDoItVector
std::vector< G4int > G4SelectedAtRestDoItVector
std::vector< G4int > G4SelectedAlongStepDoItVector
int G4int
Definition: G4Types.hh:85
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4TrackVector * NewSecondaryVector()
G4StepPoint * GetPostStepPoint() const
void SetNavigator(G4Navigator *value)
static G4int BestUnitPrecision()
static G4TransportationManager * GetTransportationManager()
static G4VSteppingVerbose * GetInstance()
virtual G4VSteppingVerbose * Clone()
static G4VSteppingVerbose * GetMasterInstance()
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4SteppingManager()

G4SteppingManager::~G4SteppingManager ( )

Definition at line 104 of file G4SteppingManager.cc.

106{
107 fTouchableHandle = 0;
108
109 // Destruct simple 'has-a' objects
110 //
111 fStep->DeleteSecondaryVector();
112
113 // delete fSecondary;
114 delete fStep;
115 delete fSelectedAtRestDoItVector;
116 delete fSelectedAlongStepDoItVector;
117 delete fSelectedPostStepDoItVector;
118 delete fUserSteppingAction;
119 #ifdef G4VERBOSE
120 if(KillVerbose) delete fVerbose;
121 #endif
122}
void DeleteSecondaryVector()

Member Function Documentation

◆ GetCorrectedStep()

G4double G4SteppingManager::GetCorrectedStep ( )
inline

Definition at line 292 of file G4SteppingManager.hh.

293 {
294 return CorrectedStep;
295 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetcurrentMinimumStep()

G4double G4SteppingManager::GetcurrentMinimumStep ( )

◆ GetfAlongStepDoItProcTriggered()

size_t G4SteppingManager::GetfAlongStepDoItProcTriggered ( )
inline

Definition at line 427 of file G4SteppingManager.hh.

428 {
429 return fAtRestDoItProcTriggered;
430 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAlongStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepDoItVector ( )
inline

Definition at line 382 of file G4SteppingManager.hh.

383 {
384 return fAlongStepDoItVector;
385 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAlongStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepGetPhysIntVector ( )
inline

Definition at line 397 of file G4SteppingManager.hh.

398 {
399 return fAlongStepGetPhysIntVector;
400 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAtRestDoItProcTriggered()

size_t G4SteppingManager::GetfAtRestDoItProcTriggered ( )
inline

Definition at line 422 of file G4SteppingManager.hh.

423 {
424 return fAtRestDoItProcTriggered;
425 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAtRestDoItVector()

G4ProcessVector * G4SteppingManager::GetfAtRestDoItVector ( )
inline

Definition at line 377 of file G4SteppingManager.hh.

378 {
379 return fAtRestDoItVector;
380 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfAtRestGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAtRestGetPhysIntVector ( )
inline

Definition at line 392 of file G4SteppingManager.hh.

393 {
394 return fAtRestGetPhysIntVector;
395 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfCondition()

G4ForceCondition G4SteppingManager::GetfCondition ( )
inline

Definition at line 500 of file G4SteppingManager.hh.

501 {
502 return fCondition;
503 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfCurrentProcess()

G4VProcess * G4SteppingManager::GetfCurrentProcess ( )
inline

Definition at line 372 of file G4SteppingManager.hh.

373 {
374 return fCurrentProcess;
375 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfCurrentVolume()

G4VPhysicalVolume * G4SteppingManager::GetfCurrentVolume ( )
inline

Definition at line 362 of file G4SteppingManager.hh.

363 {
364 return fCurrentVolume;
365 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfGPILSelection()

G4GPILSelection G4SteppingManager::GetfGPILSelection ( )
inline

Definition at line 505 of file G4SteppingManager.hh.

506 {
507 return fGPILSelection;
508 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetFirstStep()

G4bool G4SteppingManager::GetFirstStep ( )
inline

Definition at line 302 of file G4SteppingManager.hh.

303 {
304 return FirstStep;
305 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfN2ndariesAlongStepDoIt()

G4int G4SteppingManager::GetfN2ndariesAlongStepDoIt ( )
inline

Definition at line 442 of file G4SteppingManager.hh.

443 {
444 return fN2ndariesAlongStepDoIt;
445 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfN2ndariesAtRestDoIt()

G4int G4SteppingManager::GetfN2ndariesAtRestDoIt ( )
inline

Definition at line 437 of file G4SteppingManager.hh.

438 {
439 return fN2ndariesAtRestDoIt;
440 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfN2ndariesPostStepDoIt()

G4int G4SteppingManager::GetfN2ndariesPostStepDoIt ( )
inline

Definition at line 447 of file G4SteppingManager.hh.

448 {
449 return fN2ndariesPostStepDoIt;
450 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfNavigator()

G4Navigator * G4SteppingManager::GetfNavigator ( )
inline

Definition at line 452 of file G4SteppingManager.hh.

453 {
454 return fNavigator;
455 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfParticleChange()

G4VParticleChange * G4SteppingManager::GetfParticleChange ( )
inline

Definition at line 332 of file G4SteppingManager.hh.

333 {
334 return fParticleChange;
335 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepDoItProcTriggered()

size_t G4SteppingManager::GetfPostStepDoItProcTriggered ( )
inline

Definition at line 432 of file G4SteppingManager.hh.

433 {
434 return fPostStepDoItProcTriggered;
435 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfPostStepDoItVector ( )
inline

Definition at line 387 of file G4SteppingManager.hh.

388 {
389 return fPostStepDoItVector;
390 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfPostStepGetPhysIntVector ( )
inline

Definition at line 402 of file G4SteppingManager.hh.

403 {
404 return fPostStepGetPhysIntVector;
405 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPostStepPoint()

G4StepPoint * G4SteppingManager::GetfPostStepPoint ( )
inline

Definition at line 357 of file G4SteppingManager.hh.

358 {
359 return fPostStepPoint;
360 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPreStepPoint()

G4StepPoint * G4SteppingManager::GetfPreStepPoint ( )
inline

Definition at line 352 of file G4SteppingManager.hh.

353 {
354 return fPreStepPoint;
355 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfPreviousStepSize()

G4double G4SteppingManager::GetfPreviousStepSize ( )
inline

Definition at line 480 of file G4SteppingManager.hh.

481 {
482 return fPreviousStepSize;
483 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSecondary()

G4TrackVector * G4SteppingManager::GetfSecondary ( )
inline

Definition at line 342 of file G4SteppingManager.hh.

343 {
344 return fStep->GetfSecondary();
345 }
G4TrackVector * GetfSecondary()

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

◆ GetfSelectedAlongStepDoItVector()

G4SelectedAlongStepDoItVector * G4SteppingManager::GetfSelectedAlongStepDoItVector ( )
inline

Definition at line 469 of file G4SteppingManager.hh.

470 {
471 return fSelectedAlongStepDoItVector;
472 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSelectedAtRestDoItVector()

G4SelectedAtRestDoItVector * G4SteppingManager::GetfSelectedAtRestDoItVector ( )
inline

Definition at line 463 of file G4SteppingManager.hh.

464 {
465 return fSelectedAtRestDoItVector;
466 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSelectedPostStepDoItVector()

G4SelectedPostStepDoItVector * G4SteppingManager::GetfSelectedPostStepDoItVector ( )
inline

Definition at line 475 of file G4SteppingManager.hh.

476 {
477 return fSelectedPostStepDoItVector;
478 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfSensitive()

G4VSensitiveDetector * G4SteppingManager::GetfSensitive ( )
inline

Definition at line 367 of file G4SteppingManager.hh.

368 {
369 return fSensitive;
370 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfStep()

G4Step * G4SteppingManager::GetfStep ( )
inline

Definition at line 347 of file G4SteppingManager.hh.

348 {
349 return fStep;
350 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfStepStatus()

G4StepStatus G4SteppingManager::GetfStepStatus ( )
inline

Definition at line 307 of file G4SteppingManager.hh.

308 {
309 return fStepStatus;
310 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetfTrack()

G4Track * G4SteppingManager::GetfTrack ( )
inline

Definition at line 337 of file G4SteppingManager.hh.

338 {
339 return fTrack;
340 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetGeometricalStep()

G4double G4SteppingManager::GetGeometricalStep ( )
inline

Definition at line 287 of file G4SteppingManager.hh.

288 {
289 return GeometricalStep;
290 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMass()

G4double G4SteppingManager::GetMass ( )
inline

Definition at line 322 of file G4SteppingManager.hh.

323 {
324 return Mass;
325 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMAXofAlongStepLoops()

size_t G4SteppingManager::GetMAXofAlongStepLoops ( )
inline

Definition at line 412 of file G4SteppingManager.hh.

413 {
414 return MAXofAlongStepLoops;
415 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMAXofAtRestLoops()

size_t G4SteppingManager::GetMAXofAtRestLoops ( )
inline

Definition at line 407 of file G4SteppingManager.hh.

408 {
409 return MAXofAtRestLoops;
410 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetMAXofPostStepLoops()

size_t G4SteppingManager::GetMAXofPostStepLoops ( )
inline

Definition at line 417 of file G4SteppingManager.hh.

418 {
419 return MAXofPostStepLoops;
420 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetnumberOfInteractionLengthLeft()

G4double G4SteppingManager::GetnumberOfInteractionLengthLeft ( )

◆ GetPhysicalStep()

G4double G4SteppingManager::GetPhysicalStep ( )
inline

Definition at line 282 of file G4SteppingManager.hh.

283 {
284 return PhysicalStep;
285 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetphysIntLength()

G4double G4SteppingManager::GetphysIntLength ( )
inline

Definition at line 495 of file G4SteppingManager.hh.

496 {
497 return physIntLength;
498 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetPreStepPointIsGeom()

G4bool G4SteppingManager::GetPreStepPointIsGeom ( )
inline

Definition at line 297 of file G4SteppingManager.hh.

298 {
299 return PreStepPointIsGeom;
300 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetProcessNumber()

void G4SteppingManager::GetProcessNumber ( )

Definition at line 409 of file G4SteppingManager.cc.

411{
412 #ifdef debug
413 G4cout << "G4SteppingManager::GetProcessNumber: is called track="
414 << fTrack << G4endl;
415 #endif
416
418 if(pm == nullptr)
419 {
420 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
421 << " ProcessManager is NULL for particle = "
422 << fTrack->GetDefinition()->GetParticleName() << ", PDG_code = "
423 << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
424 G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011",
425 FatalException, "Process Manager is not found.");
426 return;
427 }
428
429 // AtRestDoits
430 //
431 MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
432 fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
433 fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
434
435 #ifdef debug
436 G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest="
437 << MAXofAtRestLoops << G4endl;
438 #endif
439
440 // AlongStepDoits
441 //
442 MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
443 fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
444 fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
445
446 #ifdef debug
447 G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp="
448 << MAXofAlongStepLoops << G4endl;
449 #endif
450
451 // PostStepDoits
452 //
453 MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
454 fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
455 fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
456
457 #ifdef debug
458 G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep="
459 << MAXofPostStepLoops << G4endl;
460 #endif
461
462 if (SizeOfSelectedDoItVector<MAXofAtRestLoops ||
463 SizeOfSelectedDoItVector<MAXofAlongStepLoops ||
464 SizeOfSelectedDoItVector<MAXofPostStepLoops )
465 {
466 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
467 << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
468 << " ; is smaller then one of MAXofAtRestLoops= "
469 << MAXofAtRestLoops << G4endl
470 << " or MAXofAlongStepLoops= " << MAXofAlongStepLoops
471 << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
472 G4Exception("G4SteppingManager::GetProcessNumber()",
473 "Tracking0012", FatalException,
474 "The array size is smaller than the actual No of processes.");
475 }
476}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
@ 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 510 of file G4SteppingManager.hh.

511 {
512 return fStep->GetSecondary();
513 }
const G4TrackVector * GetSecondary() const

◆ GetStep()

G4Step * G4SteppingManager::GetStep ( ) const
inline

Definition at line 545 of file G4SteppingManager.hh.

546 {
547 return fStep;
548 }

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

◆ GetStepControlFlag()

G4SteppingControl G4SteppingManager::GetStepControlFlag ( )
inline

Definition at line 490 of file G4SteppingManager.hh.

491 {
492 return StepControlFlag;
493 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetsumEnergyChange()

G4double G4SteppingManager::GetsumEnergyChange ( )
inline

Definition at line 327 of file G4SteppingManager.hh.

328 {
329 return sumEnergyChange;
330 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTempInitVelocity()

G4double G4SteppingManager::GetTempInitVelocity ( )
inline

Definition at line 312 of file G4SteppingManager.hh.

313 {
314 return TempInitVelocity;
315 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTempVelocity()

G4double G4SteppingManager::GetTempVelocity ( )
inline

Definition at line 317 of file G4SteppingManager.hh.

318 {
319 return TempVelocity;
320 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTouchableHandle()

const G4TouchableHandle & G4SteppingManager::GetTouchableHandle ( )
inline

Definition at line 485 of file G4SteppingManager.hh.

486 {
487 return fTouchableHandle;
488 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetTrack()

G4Track * G4SteppingManager::GetTrack ( ) const
inline

Definition at line 530 of file G4SteppingManager.hh.

531 {
532 return fTrack;
533 }

Referenced by G4TrackingMessenger::SetNewValue().

◆ GetUserAction()

G4UserSteppingAction * G4SteppingManager::GetUserAction ( )
inline

Definition at line 525 of file G4SteppingManager.hh.

526 {
527 return fUserSteppingAction;
528 }

Referenced by G4VSteppingVerbose::CopyState().

◆ GetverboseLevel()

G4int G4SteppingManager::GetverboseLevel ( )
inline

Definition at line 457 of file G4SteppingManager.hh.

458 {
459 return verboseLevel;
460 }

Referenced by G4VSteppingVerbose::CopyState().

◆ SetInitialStep()

void G4SteppingManager::SetInitialStep ( G4Track valueTrack)

Definition at line 288 of file G4SteppingManager.cc.

290{
291 // Set up several local variables
292 //
293 PreStepPointIsGeom = false;
294 FirstStep = true;
295 fParticleChange = nullptr;
296 fPreviousStepSize = 0.;
297 fStepStatus = fUndefined;
298
299 fTrack = valueTrack;
300 Mass = fTrack->GetDynamicParticle()->GetMass();
301
302 PhysicalStep = 0.;
303 GeometricalStep = 0.;
304 CorrectedStep = 0.;
305 PreStepPointIsGeom = false;
306 FirstStep = false;
307
308 TempInitVelocity = 0.;
309 TempVelocity = 0.;
310 sumEnergyChange = 0.;
311
312 // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
313 // set the track state to 'Alive'
314 //
315 if ( (fTrack->GetTrackStatus()==fSuspend)
316 || (fTrack->GetTrackStatus()==fPostponeToNextEvent) )
317 {
318 fTrack->SetTrackStatus(fAlive);
319 }
320
321 // If the primary track has 'zero' kinetic energy, set the track
322 // state to 'StopButAlive'
323 //
324 if(fTrack->GetKineticEnergy() <= 0.0)
325 {
326 fTrack->SetTrackStatus( fStopButAlive );
327 }
328
329 // Set Touchable to track and a private attribute of G4SteppingManager
330
331 if ( ! fTrack->GetTouchableHandle() )
332 {
333 G4ThreeVector direction= fTrack->GetMomentumDirection();
334 fNavigator->LocateGlobalPointAndSetup( fTrack->GetPosition(),
335 &direction, false, false );
336 fTouchableHandle = fNavigator->CreateTouchableHistory();
337 fTrack->SetTouchableHandle( fTouchableHandle );
338 fTrack->SetNextTouchableHandle( fTouchableHandle );
339 }
340 else
341 {
342 fTrack->SetNextTouchableHandle( fTouchableHandle = fTrack->GetTouchableHandle() );
343 G4VPhysicalVolume* oldTopVolume = fTrack->GetTouchableHandle()->GetVolume();
344 G4VPhysicalVolume* newTopVolume =
345 fNavigator->ResetHierarchyAndLocate( fTrack->GetPosition(),
346 fTrack->GetMomentumDirection(),
347 *((G4TouchableHistory*)fTrack->GetTouchableHandle()()) );
348 if ( newTopVolume != oldTopVolume
349 || oldTopVolume->GetRegularStructureId() == 1 )
350 {
351 fTouchableHandle = fNavigator->CreateTouchableHistory();
352 fTrack->SetTouchableHandle( fTouchableHandle );
353 fTrack->SetNextTouchableHandle( fTouchableHandle );
354 }
355 }
356
357 // Set OriginTouchableHandle for primary track
358 //
359 if(fTrack->GetParentID()==0)
360 {
362 }
363
364 // Set vertex information of G4Track at here
365 //
366 if ( fTrack->GetCurrentStepNumber() == 0 )
367 {
368 fTrack->SetVertexPosition( fTrack->GetPosition() );
370 fTrack->SetVertexKineticEnergy( fTrack->GetKineticEnergy() );
372 }
373
374 // Initial set up for attributes of 'G4SteppingManager'
375 fCurrentVolume = fTouchableHandle->GetVolume();
376
377 // If track is already outside the world boundary, kill it
378 //
379 if( fCurrentVolume==nullptr )
380 {
381 // If the track is a primary, stop processing
382 if(fTrack->GetParentID()==0)
383 {
384 G4cerr << "ERROR - G4SteppingManager::SetInitialStep()" << G4endl
385 << " Primary particle starting at - "
386 << fTrack->GetPosition()
387 << " - is outside of the world volume." << G4endl;
388 G4Exception("G4SteppingManager::SetInitialStep()", "Tracking0010",
389 FatalException, "Primary vertex outside of the world!");
390 }
391
392 fTrack->SetTrackStatus( fStopAndKill );
393 G4cout << "WARNING - G4SteppingManager::SetInitialStep()" << G4endl
394 << " Initial track position is outside world! - "
395 << fTrack->GetPosition() << G4endl;
396 }
397 else
398 {
399 // Initial set up for attributes of 'Step'
400 fStep->InitializeStep( fTrack );
401 }
402
403 #ifdef G4VERBOSE
404 if(verboseLevel>0) fVerbose->TrackingStarted();
405 #endif
406}
@ 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:132
virtual G4VPhysicalVolume * ResetHierarchyAndLocate(const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
Definition: G4Navigator.cc:102
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:34

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

◆ SetNavigator()

void G4SteppingManager::SetNavigator ( G4Navigator value)
inline

Definition at line 515 of file G4SteppingManager.hh.

516 {
517 fNavigator = value;
518 }

Referenced by G4SteppingManager().

◆ SetUserAction()

void G4SteppingManager::SetUserAction ( G4UserSteppingAction apAction)
inline

Definition at line 520 of file G4SteppingManager.hh.

521 {
522 fUserSteppingAction = apAction;
523 }

Referenced by G4TrackingManager::SetUserAction().

◆ SetVerbose()

void G4SteppingManager::SetVerbose ( G4VSteppingVerbose yourVerbose)
inline

Definition at line 540 of file G4SteppingManager.hh.

541 {
542 fVerbose = yourVerbose;
543 }

◆ SetVerboseLevel()

void G4SteppingManager::SetVerboseLevel ( G4int  vLevel)
inline

Definition at line 535 of file G4SteppingManager.hh.

536 {
537 verboseLevel = vLevel;
538 }

Referenced by G4ErrorPropagatorManager::SetSteppingManagerVerboseLevel().

◆ Stepping()

G4StepStatus G4SteppingManager::Stepping ( )

Definition at line 125 of file G4SteppingManager.cc.

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