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

#include <G4SteppingManager.hh>

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 73 of file G4SteppingManager.hh.

Constructor & Destructor Documentation

◆ G4SteppingManager()

G4SteppingManager::G4SteppingManager ( )

Definition at line 51 of file G4SteppingManager.cc.

53{
54 // Construct simple 'has-a' related objects
55
56 fStep = new G4Step();
57 fSecondary = fStep->NewSecondaryVector();
58 fPreStepPoint = fStep->GetPreStepPoint();
59 fPostStepPoint = fStep->GetPostStepPoint();
60
61#ifdef G4VERBOSE
63 if (fVerbose == nullptr) {
66 if (prec > 0) {
67 fVerbose = new G4SteppingVerboseWithUnits(prec);
68 }
69 else {
70 fVerbose = new G4SteppingVerbose();
71 }
72 }
73 else {
75 }
76 KillVerbose = true;
77 }
78 else {
79 KillVerbose = false;
80 }
81 fVerbose->SetManager(this);
82#endif
83
85
86 fSelectedAtRestDoItVector = new G4SelectedAtRestDoItVector(SizeOfSelectedDoItVector, 0);
87 fSelectedAlongStepDoItVector = new G4SelectedAlongStepDoItVector(SizeOfSelectedDoItVector, 0);
88 fSelectedPostStepDoItVector = new G4SelectedPostStepDoItVector(SizeOfSelectedDoItVector, 0);
89
91
92 physIntLength = DBL_MAX;
94
95 fNoProcess = new G4NoProcess;
96}
std::vector< G4int > G4SelectedPostStepDoItVector
std::vector< G4int > G4SelectedAtRestDoItVector
std::vector< G4int > G4SelectedAlongStepDoItVector
int G4int
Definition G4Types.hh:85
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
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 99 of file G4SteppingManager.cc.

101{
102 fTouchableHandle = nullptr;
103
104 // Destruct simple 'has-a' objects
105 //
106 fStep->DeleteSecondaryVector();
107
108 // delete fSecondary;
109 delete fStep;
110 delete fSelectedAtRestDoItVector;
111 delete fSelectedAlongStepDoItVector;
112 delete fSelectedPostStepDoItVector;
113 delete fUserSteppingAction;
114#ifdef G4VERBOSE
115 if (KillVerbose) delete fVerbose;
116#endif
117}

Member Function Documentation

◆ GetCorrectedStep()

G4double G4SteppingManager::GetCorrectedStep ( )
inline

Definition at line 273 of file G4SteppingManager.hh.

273{ return CorrectedStep; }

◆ GetcurrentMinimumStep()

G4double G4SteppingManager::GetcurrentMinimumStep ( )

◆ GetfAlongStepDoItProcTriggered()

size_t G4SteppingManager::GetfAlongStepDoItProcTriggered ( )
inline

Definition at line 339 of file G4SteppingManager.hh.

340{
341 return fAtRestDoItProcTriggered;
342}

◆ GetfAlongStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepDoItVector ( )
inline

Definition at line 309 of file G4SteppingManager.hh.

310{
311 return fAlongStepDoItVector;
312}

◆ GetfAlongStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepGetPhysIntVector ( )
inline

Definition at line 321 of file G4SteppingManager.hh.

322{
323 return fAlongStepGetPhysIntVector;
324}

◆ GetfAtRestDoItProcTriggered()

size_t G4SteppingManager::GetfAtRestDoItProcTriggered ( )
inline

Definition at line 337 of file G4SteppingManager.hh.

337{ return fAtRestDoItProcTriggered; }

◆ GetfAtRestDoItVector()

G4ProcessVector * G4SteppingManager::GetfAtRestDoItVector ( )
inline

Definition at line 307 of file G4SteppingManager.hh.

307{ return fAtRestDoItVector; }

◆ GetfAtRestGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAtRestGetPhysIntVector ( )
inline

Definition at line 316 of file G4SteppingManager.hh.

317{
318 return fAtRestGetPhysIntVector;
319}

◆ GetfCondition()

G4ForceCondition G4SteppingManager::GetfCondition ( )
inline

Definition at line 382 of file G4SteppingManager.hh.

382{ return fCondition; }

◆ GetfCurrentProcess()

G4VProcess * G4SteppingManager::GetfCurrentProcess ( )
inline

Definition at line 305 of file G4SteppingManager.hh.

305{ return fCurrentProcess; }

◆ GetfCurrentVolume()

G4VPhysicalVolume * G4SteppingManager::GetfCurrentVolume ( )
inline

Definition at line 301 of file G4SteppingManager.hh.

301{ return fCurrentVolume; }

◆ GetfGPILSelection()

G4GPILSelection G4SteppingManager::GetfGPILSelection ( )
inline

Definition at line 384 of file G4SteppingManager.hh.

384{ return fGPILSelection; }

◆ GetFirstStep()

G4bool G4SteppingManager::GetFirstStep ( )
inline

Definition at line 277 of file G4SteppingManager.hh.

277{ return FirstStep; }

◆ GetfN2ndariesAlongStepDoIt()

G4int G4SteppingManager::GetfN2ndariesAlongStepDoIt ( )
inline

Definition at line 351 of file G4SteppingManager.hh.

351{ return fN2ndariesAlongStepDoIt; }

◆ GetfN2ndariesAtRestDoIt()

G4int G4SteppingManager::GetfN2ndariesAtRestDoIt ( )
inline

Definition at line 349 of file G4SteppingManager.hh.

349{ return fN2ndariesAtRestDoIt; }

◆ GetfN2ndariesPostStepDoIt()

G4int G4SteppingManager::GetfN2ndariesPostStepDoIt ( )
inline

Definition at line 353 of file G4SteppingManager.hh.

353{ return fN2ndariesPostStepDoIt; }

◆ GetfNavigator()

G4Navigator * G4SteppingManager::GetfNavigator ( )
inline

Definition at line 355 of file G4SteppingManager.hh.

355{ return fNavigator; }

◆ GetfParticleChange()

G4VParticleChange * G4SteppingManager::GetfParticleChange ( )
inline

Definition at line 289 of file G4SteppingManager.hh.

289{ return fParticleChange; }

◆ GetfPostStepDoItProcTriggered()

size_t G4SteppingManager::GetfPostStepDoItProcTriggered ( )
inline

Definition at line 344 of file G4SteppingManager.hh.

345{
346 return fPostStepDoItProcTriggered;
347}

◆ GetfPostStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfPostStepDoItVector ( )
inline

Definition at line 314 of file G4SteppingManager.hh.

314{ return fPostStepDoItVector; }

◆ GetfPostStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfPostStepGetPhysIntVector ( )
inline

Definition at line 326 of file G4SteppingManager.hh.

327{
328 return fPostStepGetPhysIntVector;
329}

◆ GetfPostStepPoint()

G4StepPoint * G4SteppingManager::GetfPostStepPoint ( )
inline

Definition at line 299 of file G4SteppingManager.hh.

299{ return fPostStepPoint; }

◆ GetfPreStepPoint()

G4StepPoint * G4SteppingManager::GetfPreStepPoint ( )
inline

Definition at line 297 of file G4SteppingManager.hh.

297{ return fPreStepPoint; }

◆ GetfPreviousStepSize()

G4double G4SteppingManager::GetfPreviousStepSize ( )
inline

Definition at line 374 of file G4SteppingManager.hh.

374{ return fPreviousStepSize; }

◆ GetfSecondary()

G4TrackVector * G4SteppingManager::GetfSecondary ( )
inline

Definition at line 293 of file G4SteppingManager.hh.

293{ return fStep->GetfSecondary(); }

◆ GetfSelectedAlongStepDoItVector()

G4SelectedAlongStepDoItVector * G4SteppingManager::GetfSelectedAlongStepDoItVector ( )
inline

Definition at line 364 of file G4SteppingManager.hh.

365{
366 return fSelectedAlongStepDoItVector;
367}

◆ GetfSelectedAtRestDoItVector()

G4SelectedAtRestDoItVector * G4SteppingManager::GetfSelectedAtRestDoItVector ( )
inline

Definition at line 359 of file G4SteppingManager.hh.

360{
361 return fSelectedAtRestDoItVector;
362}

◆ GetfSelectedPostStepDoItVector()

G4SelectedPostStepDoItVector * G4SteppingManager::GetfSelectedPostStepDoItVector ( )
inline

Definition at line 369 of file G4SteppingManager.hh.

370{
371 return fSelectedPostStepDoItVector;
372}

◆ GetfSensitive()

G4VSensitiveDetector * G4SteppingManager::GetfSensitive ( )
inline

Definition at line 303 of file G4SteppingManager.hh.

303{ return fSensitive; }

◆ GetfStep()

G4Step * G4SteppingManager::GetfStep ( )
inline

Definition at line 295 of file G4SteppingManager.hh.

295{ return fStep; }

◆ GetfStepStatus()

G4StepStatus G4SteppingManager::GetfStepStatus ( )
inline

Definition at line 279 of file G4SteppingManager.hh.

279{ return fStepStatus; }

◆ GetfTrack()

G4Track * G4SteppingManager::GetfTrack ( )
inline

Definition at line 291 of file G4SteppingManager.hh.

291{ return fTrack; }

◆ GetGeometricalStep()

G4double G4SteppingManager::GetGeometricalStep ( )
inline

Definition at line 271 of file G4SteppingManager.hh.

271{ return GeometricalStep; }

◆ GetMass()

G4double G4SteppingManager::GetMass ( )
inline

Definition at line 285 of file G4SteppingManager.hh.

285{ return Mass; }

◆ GetMAXofAlongStepLoops()

size_t G4SteppingManager::GetMAXofAlongStepLoops ( )
inline

Definition at line 333 of file G4SteppingManager.hh.

333{ return MAXofAlongStepLoops; }

◆ GetMAXofAtRestLoops()

size_t G4SteppingManager::GetMAXofAtRestLoops ( )
inline

Definition at line 331 of file G4SteppingManager.hh.

331{ return MAXofAtRestLoops; }

◆ GetMAXofPostStepLoops()

size_t G4SteppingManager::GetMAXofPostStepLoops ( )
inline

Definition at line 335 of file G4SteppingManager.hh.

335{ return MAXofPostStepLoops; }

◆ GetnumberOfInteractionLengthLeft()

G4double G4SteppingManager::GetnumberOfInteractionLengthLeft ( )

◆ GetPhysicalStep()

G4double G4SteppingManager::GetPhysicalStep ( )
inline

Definition at line 269 of file G4SteppingManager.hh.

269{ return PhysicalStep; }

◆ GetphysIntLength()

G4double G4SteppingManager::GetphysIntLength ( )
inline

Definition at line 380 of file G4SteppingManager.hh.

380{ return physIntLength; }

◆ GetPreStepPointIsGeom()

G4bool G4SteppingManager::GetPreStepPointIsGeom ( )
inline

Definition at line 275 of file G4SteppingManager.hh.

275{ return PreStepPointIsGeom; }

◆ GetProcessNumber()

void G4SteppingManager::GetProcessNumber ( )

Definition at line 374 of file G4SteppingManager.cc.

376{
377#ifdef debug
378 G4cout << "G4SteppingManager::GetProcessNumber: is called track=" << fTrack << G4endl;
379#endif
380
381 G4ProcessManager* pm = fTrack->GetDefinition()->GetProcessManager();
382 if (pm == nullptr) {
383 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
384 << " ProcessManager is NULL for particle = "
385 << fTrack->GetDefinition()->GetParticleName()
386 << ", PDG_code = " << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
387 G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011", FatalException,
388 "Process Manager is not found.");
389 return;
390 }
391
392 // AtRestDoits
393 //
394 MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
395 fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
396 fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
397
398#ifdef debug
399 G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest=" << MAXofAtRestLoops << G4endl;
400#endif
401
402 // AlongStepDoits
403 //
404 MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
405 fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
406 fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
407
408#ifdef debug
409 G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp=" << MAXofAlongStepLoops << G4endl;
410#endif
411
412 // PostStepDoits
413 //
414 MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
415 fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
416 fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
417
418#ifdef debug
419 G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep=" << MAXofPostStepLoops << G4endl;
420#endif
421
422 if (SizeOfSelectedDoItVector < MAXofAtRestLoops ||
423 SizeOfSelectedDoItVector < MAXofAlongStepLoops ||
424 SizeOfSelectedDoItVector < MAXofPostStepLoops)
425 {
426 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
427 << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
428 << " ; is smaller then one of MAXofAtRestLoops= " << MAXofAtRestLoops << G4endl
429 << " or MAXofAlongStepLoops= " << MAXofAlongStepLoops
430 << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
431 G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0012", FatalException,
432 "The array size is smaller than the actual No of processes.");
433 }
434}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ typeGPIL
@ typeDoIt
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t entries() const

◆ GetSecondary()

const G4TrackVector * G4SteppingManager::GetSecondary ( ) const
inline

Definition at line 386 of file G4SteppingManager.hh.

387{
388 return fStep->GetSecondary();
389}

◆ GetStep()

G4Step * G4SteppingManager::GetStep ( ) const
inline

Definition at line 409 of file G4SteppingManager.hh.

409{ return fStep; }

◆ GetStepControlFlag()

G4SteppingControl G4SteppingManager::GetStepControlFlag ( )
inline

Definition at line 378 of file G4SteppingManager.hh.

378{ return StepControlFlag; }

◆ GetsumEnergyChange()

G4double G4SteppingManager::GetsumEnergyChange ( )
inline

Definition at line 287 of file G4SteppingManager.hh.

287{ return sumEnergyChange; }

◆ GetTempInitVelocity()

G4double G4SteppingManager::GetTempInitVelocity ( )
inline

Definition at line 281 of file G4SteppingManager.hh.

281{ return TempInitVelocity; }

◆ GetTempVelocity()

G4double G4SteppingManager::GetTempVelocity ( )
inline

Definition at line 283 of file G4SteppingManager.hh.

283{ return TempVelocity; }

◆ GetTouchableHandle()

const G4TouchableHandle & G4SteppingManager::GetTouchableHandle ( )
inline

Definition at line 376 of file G4SteppingManager.hh.

376{ return fTouchableHandle; }

◆ GetTrack()

G4Track * G4SteppingManager::GetTrack ( ) const
inline

Definition at line 400 of file G4SteppingManager.hh.

400{ return fTrack; }

◆ GetUserAction()

G4UserSteppingAction * G4SteppingManager::GetUserAction ( )
inline

Definition at line 398 of file G4SteppingManager.hh.

398{ return fUserSteppingAction; }

◆ GetverboseLevel()

G4int G4SteppingManager::GetverboseLevel ( )
inline

Definition at line 357 of file G4SteppingManager.hh.

357{ return verboseLevel; }

◆ SetInitialStep()

void G4SteppingManager::SetInitialStep ( G4Track * valueTrack)

Definition at line 268 of file G4SteppingManager.cc.

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

◆ SetNavigator()

void G4SteppingManager::SetNavigator ( G4Navigator * value)
inline

Definition at line 391 of file G4SteppingManager.hh.

391{ fNavigator = value; }

Referenced by G4SteppingManager().

◆ SetUserAction()

void G4SteppingManager::SetUserAction ( G4UserSteppingAction * apAction)
inline

Definition at line 393 of file G4SteppingManager.hh.

394{
395 fUserSteppingAction = apAction;
396}

◆ SetVerbose()

void G4SteppingManager::SetVerbose ( G4VSteppingVerbose * yourVerbose)
inline

Definition at line 404 of file G4SteppingManager.hh.

405{
406 fVerbose = yourVerbose;
407}

◆ SetVerboseLevel()

void G4SteppingManager::SetVerboseLevel ( G4int vLevel)
inline

Definition at line 402 of file G4SteppingManager.hh.

402{ verboseLevel = vLevel; }

Referenced by G4ErrorPropagatorManager::SetSteppingManagerVerboseLevel().

◆ Stepping()

G4StepStatus G4SteppingManager::Stepping ( )

Definition at line 120 of file G4SteppingManager.cc.

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

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