Geant4  10.02.p03
G4SteppingManager Class Reference

#include <G4SteppingManager.hh>

Collaboration diagram for G4SteppingManager:

Public Member Functions

 G4SteppingManager ()
 
 ~G4SteppingManager ()
 
const G4TrackVector * GetSecondary () const
 
void SetUserAction (G4UserSteppingAction *apAction)
 
G4Track * GetTrack () const
 
void SetVerboseLevel (G4int vLevel)
 
void SetVerbose (G4VSteppingVerbose *)
 
G4Step * GetStep () 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 ()
 
G4VParticleChange * GetfParticleChange ()
 
G4Track * GetfTrack ()
 
G4TrackVector * GetfSecondary ()
 
G4Step * GetfStep ()
 
G4StepPoint * GetfPreStepPoint ()
 
G4StepPoint * GetfPostStepPoint ()
 
G4VPhysicalVolumeGetfCurrentVolume ()
 
G4VSensitiveDetectorGetfSensitive ()
 
G4VProcessGetfCurrentProcess ()
 
G4ProcessVectorGetfAtRestDoItVector ()
 
G4ProcessVectorGetfAlongStepDoItVector ()
 
G4ProcessVectorGetfPostStepDoItVector ()
 
G4ProcessVectorGetfAlongStepGetPhysIntVector ()
 
G4ProcessVectorGetfPostStepGetPhysIntVector ()
 
G4ProcessVectorGetfAtRestGetPhysIntVector ()
 
G4double GetcurrentMinimumStep ()
 
G4double GetnumberOfInteractionLengthLeft ()
 
size_t GetfAtRestDoItProcTriggered ()
 
size_t GetfAlongStepDoItProcTriggered ()
 
size_t GetfPostStepDoItProcTriggered ()
 
G4int GetfN2ndariesAtRestDoIt ()
 
G4int GetfN2ndariesAlongStepDoIt ()
 
G4int GetfN2ndariesPostStepDoIt ()
 
G4NavigatorGetfNavigator ()
 
G4int GetverboseLevel ()
 
size_t GetMAXofAtRestLoops ()
 
size_t GetMAXofAlongStepLoops ()
 
size_t GetMAXofPostStepLoops ()
 
G4SelectedAtRestDoItVectorGetfSelectedAtRestDoItVector ()
 
G4SelectedAlongStepDoItVectorGetfSelectedAlongStepDoItVector ()
 
G4SelectedPostStepDoItVectorGetfSelectedPostStepDoItVector ()
 
G4double GetfPreviousStepSize ()
 
const G4TouchableHandleGetTouchableHandle ()
 
G4SteppingControl GetStepControlFlag ()
 
G4UserSteppingActionGetUserAction ()
 
G4double GetphysIntLength ()
 
G4ForceCondition GetfCondition ()
 
G4GPILSelection GetfGPILSelection ()
 

Public Attributes

G4bool KillVerbose
 

Private Member Functions

void DefinePhysicalStepLength ()
 
void InvokeAtRestDoItProcs ()
 
void InvokeAlongStepDoItProcs ()
 
void InvokePostStepDoItProcs ()
 
void InvokePSDIP (size_t)
 
G4double CalculateSafety ()
 
void ApplyProductionCut (G4Track *)
 

Private Attributes

G4UserSteppingActionfUserSteppingAction
 
G4VSteppingVerbosefVerbose
 
G4double PhysicalStep
 
G4double GeometricalStep
 
G4double CorrectedStep
 
G4bool PreStepPointIsGeom
 
G4bool FirstStep
 
G4StepStatus fStepStatus
 
G4double TempInitVelocity
 
G4double TempVelocity
 
G4double Mass
 
G4double sumEnergyChange
 
G4VParticleChange * fParticleChange
 
G4Track * fTrack
 
G4TrackVector * fSecondary
 
G4Step * fStep
 
G4StepPoint * fPreStepPoint
 
G4StepPoint * fPostStepPoint
 
G4VPhysicalVolumefCurrentVolume
 
G4VSensitiveDetectorfSensitive
 
G4VProcessfCurrentProcess
 
G4ProcessVectorfAtRestDoItVector
 
G4ProcessVectorfAlongStepDoItVector
 
G4ProcessVectorfPostStepDoItVector
 
G4ProcessVectorfAtRestGetPhysIntVector
 
G4ProcessVectorfAlongStepGetPhysIntVector
 
G4ProcessVectorfPostStepGetPhysIntVector
 
size_t MAXofAtRestLoops
 
size_t MAXofAlongStepLoops
 
size_t MAXofPostStepLoops
 
size_t fAtRestDoItProcTriggered
 
size_t fAlongStepDoItProcTriggered
 
size_t fPostStepDoItProcTriggered
 
G4int fN2ndariesAtRestDoIt
 
G4int fN2ndariesAlongStepDoIt
 
G4int fN2ndariesPostStepDoIt
 
G4NavigatorfNavigator
 
G4int verboseLevel
 
G4SelectedAtRestDoItVectorfSelectedAtRestDoItVector
 
G4SelectedAlongStepDoItVectorfSelectedAlongStepDoItVector
 
G4SelectedPostStepDoItVectorfSelectedPostStepDoItVector
 
G4double fPreviousStepSize
 
G4TouchableHandle fTouchableHandle
 
G4SteppingControl StepControlFlag
 
G4double kCarTolerance
 
G4double proposedSafety
 
G4ThreeVector endpointSafOrigin
 
G4double endpointSafety
 
G4double physIntLength
 
G4ForceCondition fCondition
 
G4GPILSelection fGPILSelection
 

Detailed Description

Definition at line 92 of file G4SteppingManager.hh.

Constructor & Destructor Documentation

◆ G4SteppingManager()

G4SteppingManager::G4SteppingManager ( )

Definition at line 56 of file G4SteppingManager.cc.

59 {
60 
61 // Construct simple 'has-a' related objects
62  fStep = new G4Step();
63  fSecondary = fStep->NewSecondaryVector();
64  fPreStepPoint = fStep->GetPreStepPoint();
65  fPostStepPoint = fStep->GetPostStepPoint();
66 #ifdef G4VERBOSE
70  fVerbose -> SetManager(this);
71  KillVerbose = true;
72  }
73  else {
75  fVerbose -> SetManager(this);
76  KillVerbose = false;
77  }
78 #endif
80  ->GetNavigatorForTracking());
81 
88 
90  ->GetNavigatorForTracking());
91 
94 }
static G4VSteppingVerbose * GetInstance()
void SetNavigator(G4Navigator *value)
class std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
class std::vector< int, std::allocator< int > > G4SelectedAlongStepDoItVector
G4VSteppingVerbose * fVerbose
G4UserSteppingAction * fUserSteppingAction
G4double GetSurfaceTolerance() const
G4SelectedPostStepDoItVector * fSelectedPostStepDoItVector
G4SelectedAlongStepDoItVector * fSelectedAlongStepDoItVector
class std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
static G4TransportationManager * GetTransportationManager()
static void SetInstance(G4VSteppingVerbose *Instance)
G4StepPoint * fPostStepPoint
G4TrackVector * fSecondary
#define DBL_MAX
Definition: templates.hh:83
G4SelectedAtRestDoItVector * fSelectedAtRestDoItVector
static G4GeometryTolerance * GetInstance()
static const size_t SizeOfSelectedDoItVector
G4StepPoint * fPreStepPoint
Here is the call graph for this function:

◆ ~G4SteppingManager()

G4SteppingManager::~G4SteppingManager ( )

Definition at line 97 of file G4SteppingManager.cc.

99 {
100  fTouchableHandle = 0;
101 // Destruct simple 'has-a' objects
102  fStep->DeleteSecondaryVector();
104  delete fStep;
109 #ifdef G4VERBOSE
110  if(KillVerbose) delete fVerbose;
111 #endif
112 }
G4VSteppingVerbose * fVerbose
G4UserSteppingAction * fUserSteppingAction
G4SelectedPostStepDoItVector * fSelectedPostStepDoItVector
G4SelectedAlongStepDoItVector * fSelectedAlongStepDoItVector
G4SelectedAtRestDoItVector * fSelectedAtRestDoItVector
G4TouchableHandle fTouchableHandle

Member Function Documentation

◆ ApplyProductionCut()

void G4SteppingManager::ApplyProductionCut ( G4Track *  aSecondary)
private

Definition at line 592 of file G4SteppingManager2.cc.

593 {
594  G4bool tBelowCutEnergyAndSafety = false;
595  G4int tPtclIdx
596  = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
597  if (tPtclIdx<0) { return; }
598  G4ProductionCutsTable* tCutsTbl
600  G4int tCoupleIdx
601  = tCutsTbl->GetCoupleIndex(fPreStepPoint->GetMaterialCutsCouple());
602  G4double tProdThreshold
603  = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
604  if( aSecondary->GetKineticEnergy()<tProdThreshold )
605  {
606  tBelowCutEnergyAndSafety = true;
607  if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
608  {
609  G4double currentRange
610  = G4LossTableManager::Instance()->GetRange(aSecondary->GetDefinition(),
611  aSecondary->GetKineticEnergy(),
612  fPreStepPoint->GetMaterialCutsCouple());
613  tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
614  }
615  }
616 
617  if( tBelowCutEnergyAndSafety )
618  {
619  if( !(aSecondary->IsGoodForTracking()) )
620  {
621 //G4cout << "!! Warning - G4SteppingManager:" << G4endl
622 //<< " This physics process generated a secondary"
623 //<< " of which energy is below cut but"
624 //<< " GoodForTracking is off !!!!!" << G4endl;
625  // Add kinetic energy to the total energy deposit
626  fStep->AddTotalEnergyDeposit(
627  aSecondary->GetKineticEnergy() );
628  aSecondary->SetKineticEnergy(0.0);
629  }
630  }
631 }
static G4LossTableManager * Instance()
static G4int GetIndex(const G4String &name)
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
static G4ProductionCutsTable * GetProductionCutsTable()
#define DBL_MIN
Definition: templates.hh:75
G4int GetCoupleIndex(const G4MaterialCutsCouple *aCouple) const
double G4double
Definition: G4Types.hh:76
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4StepPoint * fPreStepPoint
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalculateSafety()

G4double G4SteppingManager::CalculateSafety ( )
inlineprivate

Definition at line 490 of file G4SteppingManager.hh.

490  {
491  return std::max( endpointSafety -
492  (endpointSafOrigin - fPostStepPoint->GetPosition()).mag(),
493  kCarTolerance );
494  }
G4ThreeVector endpointSafOrigin
G4StepPoint * fPostStepPoint
Here is the caller graph for this function:

◆ DefinePhysicalStepLength()

void G4SteppingManager::DefinePhysicalStepLength ( )
private

Definition at line 128 of file G4SteppingManager2.cc.

130 {
131 
132 // ReSet the counter etc.
133  PhysicalStep = DBL_MAX; // Initialize by a huge number
134  physIntLength = DBL_MAX; // Initialize by a huge number
135 #ifdef G4VERBOSE
136  // !!!!! Verbose
138 #endif
139 
140 // Obtain the user defined maximum allowed Step in the volume
141 // 1997.12.13 adds argument for GetMaxAllowedStep by K.Kurashige
142 // 2004.01.20 This block will be removed by Geant4 7.0
143 // G4UserLimits* ul= fCurrentVolume->GetLogicalVolume()->GetUserLimits();
144 // if (ul) {
145 // physIntLength = ul->GetMaxAllowedStep(*fTrack);
146 //#ifdef G4VERBOSE
147 // // !!!!! Verbose
148 // if(verboseLevel>0) fVerbose->DPSLUserLimit();
149 //#endif
150 // }
151 //
152 // if(physIntLength < PhysicalStep ){
153 // PhysicalStep = physIntLength;
154 // fStepStatus = fUserDefinedLimit;
155 // fStep->GetPostStepPoint()
156 // ->SetProcessDefinedStep(0);
157 // // Take note that the process pointer is 'NULL' if the Step
158 // // is defined by the user defined limit.
159 // }
160 // 2004.01.20 This block will be removed by Geant4 7.0
161 
162 // GPIL for PostStep
164 
165  for(size_t np=0; np < MAXofPostStepLoops; np++){
166  fCurrentProcess = (*fPostStepGetPhysIntVector)(np);
167  if (fCurrentProcess== 0) {
168  (*fSelectedPostStepDoItVector)[np] = InActivated;
169  continue;
170  } // NULL means the process is inactivated by a user on fly.
171 
173  PostStepGPIL( *fTrack,
175  &fCondition );
176 #ifdef G4VERBOSE
177  // !!!!! Verbose
179 #endif
180 
181  switch (fCondition) {
182  case ExclusivelyForced:
183  (*fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
184  fStepStatus = fExclusivelyForcedProc;
185  fStep->GetPostStepPoint()
186  ->SetProcessDefinedStep(fCurrentProcess);
187  break;
188  case Conditionally:
189  // (*fSelectedPostStepDoItVector)[np] = Conditionally;
190  G4Exception("G4SteppingManager::DefinePhysicalStepLength()", "Tracking1001", FatalException, "This feature no more supported");
191 
192  break;
193  case Forced:
194  (*fSelectedPostStepDoItVector)[np] = Forced;
195  break;
196  case StronglyForced:
197  (*fSelectedPostStepDoItVector)[np] = StronglyForced;
198  break;
199  default:
200  (*fSelectedPostStepDoItVector)[np] = InActivated;
201  break;
202  }
203 
204 
205 
206  if (fCondition==ExclusivelyForced) {
207  for(size_t nrest=np+1; nrest < MAXofPostStepLoops; nrest++){
208  (*fSelectedPostStepDoItVector)[nrest] = InActivated;
209  }
210  return; // Take note the 'return' at here !!!
211  }
212  else{
215  fStepStatus = fPostStepDoItProc;
217  fStep->GetPostStepPoint()
218  ->SetProcessDefinedStep(fCurrentProcess);
219  }
220  }
221 
222 
223  }
224 
225  if (fPostStepDoItProcTriggered<MAXofPostStepLoops) {
227  InActivated) {
228  (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] =
229  NotForced;
230  }
231  }
232 
233 // GPIL for AlongStep
235  G4double safetyProposedToAndByProcess = proposedSafety;
236 
237  for(size_t kp=0; kp < MAXofAlongStepLoops; kp++){
238  fCurrentProcess = (*fAlongStepGetPhysIntVector)[kp];
239  if (fCurrentProcess== 0) continue;
240  // NULL means the process is inactivated by a user on fly.
241 
243  AlongStepGPIL( *fTrack, fPreviousStepSize,
244  PhysicalStep,
245  safetyProposedToAndByProcess,
246  &fGPILSelection );
247 #ifdef G4VERBOSE
248  // !!!!! Verbose
250 #endif
253 
254  // Check if the process wants to be the GPIL winner. For example,
255  // multi-scattering proposes Step limit, but won't be the winner.
256  if(fGPILSelection==CandidateForSelection){
257  fStepStatus = fAlongStepDoItProc;
258  fStep->GetPostStepPoint()
259  ->SetProcessDefinedStep(fCurrentProcess);
260  }
261 
262  // Transportation is assumed to be the last process in the vector
263  if(kp == MAXofAlongStepLoops-1)
264  fStepStatus = fGeomBoundary;
265  }
266 
267  // Make sure to check the safety, even if Step is not limited
268  // by this process. J. Apostolakis, June 20, 1998
269  //
270  if (safetyProposedToAndByProcess < proposedSafety)
271  // proposedSafety keeps the smallest value:
272  proposedSafety = safetyProposedToAndByProcess;
273  else
274  // safetyProposedToAndByProcess always proposes a valid safety:
275  safetyProposedToAndByProcess = proposedSafety;
276 
277  }
278 } // void G4SteppingManager::DefinePhysicalStepLength() //
G4VSteppingVerbose * fVerbose
G4SelectedPostStepDoItVector * fSelectedPostStepDoItVector
int G4int
Definition: G4Types.hh:78
G4VProcess * fCurrentProcess
virtual void DPSLAlongStep()=0
virtual void DPSLPostStep()=0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ForceCondition fCondition
G4GPILSelection fGPILSelection
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
virtual void DPSLStarted()=0
G4StepStatus fStepStatus
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCorrectedStep()

G4double G4SteppingManager::GetCorrectedStep ( )
inline

Definition at line 307 of file G4SteppingManager.hh.

307  {
308  return CorrectedStep;
309  }
Here is the caller graph for this function:

◆ GetcurrentMinimumStep()

G4double G4SteppingManager::GetcurrentMinimumStep ( )

◆ GetfAlongStepDoItProcTriggered()

size_t G4SteppingManager::GetfAlongStepDoItProcTriggered ( )
inline

Definition at line 403 of file G4SteppingManager.hh.

403  {
405  }
Here is the caller graph for this function:

◆ GetfAlongStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepDoItVector ( )
inline

Definition at line 371 of file G4SteppingManager.hh.

371  {
372  return fAlongStepDoItVector;
373  }
G4ProcessVector * fAlongStepDoItVector
Here is the caller graph for this function:

◆ GetfAlongStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAlongStepGetPhysIntVector ( )
inline

Definition at line 382 of file G4SteppingManager.hh.

382  {
384  }
G4ProcessVector * fAlongStepGetPhysIntVector
Here is the caller graph for this function:

◆ GetfAtRestDoItProcTriggered()

size_t G4SteppingManager::GetfAtRestDoItProcTriggered ( )
inline

Definition at line 400 of file G4SteppingManager.hh.

400  {
402  }
Here is the caller graph for this function:

◆ GetfAtRestDoItVector()

G4ProcessVector * G4SteppingManager::GetfAtRestDoItVector ( )
inline

Definition at line 368 of file G4SteppingManager.hh.

368  {
369  return fAtRestDoItVector;
370  }
G4ProcessVector * fAtRestDoItVector
Here is the caller graph for this function:

◆ GetfAtRestGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfAtRestGetPhysIntVector ( )
inline

Definition at line 378 of file G4SteppingManager.hh.

378  {
380  }
G4ProcessVector * fAtRestGetPhysIntVector
Here is the caller graph for this function:

◆ GetfCondition()

G4ForceCondition G4SteppingManager::GetfCondition ( )
inline

Definition at line 452 of file G4SteppingManager.hh.

452  {
453  return fCondition;
454  }
G4ForceCondition fCondition
Here is the caller graph for this function:

◆ GetfCurrentProcess()

G4VProcess * G4SteppingManager::GetfCurrentProcess ( )
inline

Definition at line 364 of file G4SteppingManager.hh.

364  {
365  return fCurrentProcess;
366  }
G4VProcess * fCurrentProcess
Here is the caller graph for this function:

◆ GetfCurrentVolume()

G4VPhysicalVolume * G4SteppingManager::GetfCurrentVolume ( )
inline

Definition at line 358 of file G4SteppingManager.hh.

358  {
359  return fCurrentVolume;
360  }
G4VPhysicalVolume * fCurrentVolume
Here is the caller graph for this function:

◆ GetfGPILSelection()

G4GPILSelection G4SteppingManager::GetfGPILSelection ( )
inline

Definition at line 455 of file G4SteppingManager.hh.

455  {
456  return fGPILSelection;
457  }
G4GPILSelection fGPILSelection
Here is the caller graph for this function:

◆ GetFirstStep()

G4bool G4SteppingManager::GetFirstStep ( )
inline

Definition at line 315 of file G4SteppingManager.hh.

315  {
316  return FirstStep;
317  }
Here is the caller graph for this function:

◆ GetfN2ndariesAlongStepDoIt()

G4int G4SteppingManager::GetfN2ndariesAlongStepDoIt ( )
inline

Definition at line 412 of file G4SteppingManager.hh.

412  {
414  }
Here is the caller graph for this function:

◆ GetfN2ndariesAtRestDoIt()

G4int G4SteppingManager::GetfN2ndariesAtRestDoIt ( )
inline

Definition at line 409 of file G4SteppingManager.hh.

409  {
410  return fN2ndariesAtRestDoIt;
411  }
Here is the caller graph for this function:

◆ GetfN2ndariesPostStepDoIt()

G4int G4SteppingManager::GetfN2ndariesPostStepDoIt ( )
inline

Definition at line 415 of file G4SteppingManager.hh.

415  {
416  return fN2ndariesPostStepDoIt;
417  }
Here is the caller graph for this function:

◆ GetfNavigator()

G4Navigator * G4SteppingManager::GetfNavigator ( )
inline

Definition at line 419 of file G4SteppingManager.hh.

419  {
420  return fNavigator;
421  }
G4Navigator * fNavigator
Here is the caller graph for this function:

◆ GetfParticleChange()

G4VParticleChange * G4SteppingManager::GetfParticleChange ( )
inline

Definition at line 337 of file G4SteppingManager.hh.

337  {
338  return fParticleChange;
339  }
G4VParticleChange * fParticleChange
Here is the caller graph for this function:

◆ GetfPostStepDoItProcTriggered()

size_t G4SteppingManager::GetfPostStepDoItProcTriggered ( )
inline

Definition at line 406 of file G4SteppingManager.hh.

406  {
408  }
Here is the caller graph for this function:

◆ GetfPostStepDoItVector()

G4ProcessVector * G4SteppingManager::GetfPostStepDoItVector ( )
inline

Definition at line 374 of file G4SteppingManager.hh.

374  {
375  return fPostStepDoItVector;
376  }
G4ProcessVector * fPostStepDoItVector
Here is the caller graph for this function:

◆ GetfPostStepGetPhysIntVector()

G4ProcessVector * G4SteppingManager::GetfPostStepGetPhysIntVector ( )
inline

Definition at line 386 of file G4SteppingManager.hh.

386  {
388  }
G4ProcessVector * fPostStepGetPhysIntVector
Here is the caller graph for this function:

◆ GetfPostStepPoint()

G4StepPoint * G4SteppingManager::GetfPostStepPoint ( )
inline

Definition at line 354 of file G4SteppingManager.hh.

354  {
355  return fPostStepPoint;
356  }
G4StepPoint * fPostStepPoint
Here is the caller graph for this function:

◆ GetfPreStepPoint()

G4StepPoint * G4SteppingManager::GetfPreStepPoint ( )
inline

Definition at line 351 of file G4SteppingManager.hh.

351  {
352  return fPreStepPoint;
353  }
G4StepPoint * fPreStepPoint
Here is the caller graph for this function:

◆ GetfPreviousStepSize()

G4double G4SteppingManager::GetfPreviousStepSize ( )
inline

Definition at line 438 of file G4SteppingManager.hh.

438  {
439  return fPreviousStepSize;
440  }
Here is the caller graph for this function:

◆ GetfSecondary()

G4TrackVector * G4SteppingManager::GetfSecondary ( )
inline

Definition at line 345 of file G4SteppingManager.hh.

345  {
346  return fStep->GetfSecondary();
347  }
Here is the caller graph for this function:

◆ GetfSelectedAlongStepDoItVector()

G4SelectedAlongStepDoItVector * G4SteppingManager::GetfSelectedAlongStepDoItVector ( )
inline

Definition at line 430 of file G4SteppingManager.hh.

430  {
432  }
G4SelectedAlongStepDoItVector * fSelectedAlongStepDoItVector
Here is the caller graph for this function:

◆ GetfSelectedAtRestDoItVector()

G4SelectedAtRestDoItVector * G4SteppingManager::GetfSelectedAtRestDoItVector ( )
inline

Definition at line 426 of file G4SteppingManager.hh.

426  {
428  }
G4SelectedAtRestDoItVector * fSelectedAtRestDoItVector
Here is the caller graph for this function:

◆ GetfSelectedPostStepDoItVector()

G4SelectedPostStepDoItVector * G4SteppingManager::GetfSelectedPostStepDoItVector ( )
inline

Definition at line 434 of file G4SteppingManager.hh.

434  {
436  }
G4SelectedPostStepDoItVector * fSelectedPostStepDoItVector
Here is the caller graph for this function:

◆ GetfSensitive()

G4VSensitiveDetector * G4SteppingManager::GetfSensitive ( )
inline

Definition at line 361 of file G4SteppingManager.hh.

361  {
362  return fSensitive;
363  }
G4VSensitiveDetector * fSensitive
Here is the caller graph for this function:

◆ GetfStep()

G4Step * G4SteppingManager::GetfStep ( )
inline

Definition at line 348 of file G4SteppingManager.hh.

348  {
349  return fStep;
350  }
Here is the caller graph for this function:

◆ GetfStepStatus()

G4StepStatus G4SteppingManager::GetfStepStatus ( )
inline

Definition at line 319 of file G4SteppingManager.hh.

319  {
320  return fStepStatus;
321  }
G4StepStatus fStepStatus
Here is the caller graph for this function:

◆ GetfTrack()

G4Track * G4SteppingManager::GetfTrack ( )
inline

Definition at line 341 of file G4SteppingManager.hh.

341  {
342  return fTrack;
343  }
Here is the caller graph for this function:

◆ GetGeometricalStep()

G4double G4SteppingManager::GetGeometricalStep ( )
inline

Definition at line 303 of file G4SteppingManager.hh.

303  {
304  return GeometricalStep;
305  }
Here is the caller graph for this function:

◆ GetMass()

G4double G4SteppingManager::GetMass ( )
inline

Definition at line 329 of file G4SteppingManager.hh.

329  {
330  return Mass;
331  }
Here is the caller graph for this function:

◆ GetMAXofAlongStepLoops()

size_t G4SteppingManager::GetMAXofAlongStepLoops ( )
inline

Definition at line 393 of file G4SteppingManager.hh.

393  {
394  return MAXofAlongStepLoops;
395  }
Here is the caller graph for this function:

◆ GetMAXofAtRestLoops()

size_t G4SteppingManager::GetMAXofAtRestLoops ( )
inline

Definition at line 390 of file G4SteppingManager.hh.

390  {
391  return MAXofAtRestLoops;
392  }
Here is the caller graph for this function:

◆ GetMAXofPostStepLoops()

size_t G4SteppingManager::GetMAXofPostStepLoops ( )
inline

Definition at line 396 of file G4SteppingManager.hh.

396  {
397  return MAXofPostStepLoops;
398  }
Here is the caller graph for this function:

◆ GetnumberOfInteractionLengthLeft()

G4double G4SteppingManager::GetnumberOfInteractionLengthLeft ( )

◆ GetPhysicalStep()

G4double G4SteppingManager::GetPhysicalStep ( )
inline

Definition at line 299 of file G4SteppingManager.hh.

299  {
300  return PhysicalStep;
301  }
Here is the caller graph for this function:

◆ GetphysIntLength()

G4double G4SteppingManager::GetphysIntLength ( )
inline

Definition at line 449 of file G4SteppingManager.hh.

449  {
450  return physIntLength;
451  }
Here is the caller graph for this function:

◆ GetPreStepPointIsGeom()

G4bool G4SteppingManager::GetPreStepPointIsGeom ( )
inline

Definition at line 311 of file G4SteppingManager.hh.

311  {
312  return PreStepPointIsGeom;
313  }
Here is the caller graph for this function:

◆ GetProcessNumber()

void G4SteppingManager::GetProcessNumber ( )

Definition at line 56 of file G4SteppingManager2.cc.

58 {
59 #ifdef debug
60  G4cout<<"G4SteppingManager::GetProcessNumber: is called track="
61  <<fTrack<<G4endl;
62 #endif
63 
64  G4ProcessManager* pm= fTrack->GetDefinition()->GetProcessManager();
65  if(!pm)
66  {
67  G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
68  << " ProcessManager is NULL for particle = "
69  << fTrack->GetDefinition()->GetParticleName() << ", PDG_code = "
70  << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
71  G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011",
72  FatalException, "Process Manager is not found.");
73  return;
74  }
75 
76 // AtRestDoits
80 #ifdef debug
81  G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest="
83 #endif
84 
85 // AlongStepDoits
89 #ifdef debug
90  G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp="
92 #endif
93 
94 // PostStepDoits
98 #ifdef debug
99  G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep="
101 #endif
102 
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 }
G4ProcessVector * fPostStepGetPhysIntVector
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * fPostStepDoItVector
G4ProcessVector * fAlongStepGetPhysIntVector
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4int entries() const
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * fAlongStepDoItVector
#define G4endl
Definition: G4ios.hh:61
G4ProcessVector * fAtRestGetPhysIntVector
G4ProcessVector * fAtRestDoItVector
static const size_t SizeOfSelectedDoItVector
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetSecondary()

const G4TrackVector * G4SteppingManager::GetSecondary ( ) const
inline

Definition at line 459 of file G4SteppingManager.hh.

459  {
460  return fStep->GetSecondary();
461  }
Here is the caller graph for this function:

◆ GetStep()

G4Step * G4SteppingManager::GetStep ( ) const
inline

Definition at line 486 of file G4SteppingManager.hh.

486  {
487  return fStep;
488  }
Here is the caller graph for this function:

◆ GetStepControlFlag()

G4SteppingControl G4SteppingManager::GetStepControlFlag ( )
inline

Definition at line 446 of file G4SteppingManager.hh.

446  {
447  return StepControlFlag;
448  }
G4SteppingControl StepControlFlag
Here is the caller graph for this function:

◆ GetsumEnergyChange()

G4double G4SteppingManager::GetsumEnergyChange ( )
inline

Definition at line 333 of file G4SteppingManager.hh.

333  {
334  return sumEnergyChange;
335  }
Here is the caller graph for this function:

◆ GetTempInitVelocity()

G4double G4SteppingManager::GetTempInitVelocity ( )
inline

Definition at line 323 of file G4SteppingManager.hh.

323  {
324  return TempInitVelocity;
325  }
Here is the caller graph for this function:

◆ GetTempVelocity()

G4double G4SteppingManager::GetTempVelocity ( )
inline

Definition at line 326 of file G4SteppingManager.hh.

326  {
327  return TempVelocity;
328  }
Here is the caller graph for this function:

◆ GetTouchableHandle()

const G4TouchableHandle & G4SteppingManager::GetTouchableHandle ( )
inline

Definition at line 442 of file G4SteppingManager.hh.

442  {
443  return fTouchableHandle;
444  }
G4TouchableHandle fTouchableHandle
Here is the caller graph for this function:

◆ GetTrack()

G4Track * G4SteppingManager::GetTrack ( ) const
inline

Definition at line 474 of file G4SteppingManager.hh.

474  {
475  return fTrack;
476  }
Here is the caller graph for this function:

◆ GetUserAction()

G4UserSteppingAction * G4SteppingManager::GetUserAction ( )
inline

Definition at line 470 of file G4SteppingManager.hh.

470  {
471  return fUserSteppingAction;
472  }
G4UserSteppingAction * fUserSteppingAction
Here is the caller graph for this function:

◆ GetverboseLevel()

G4int G4SteppingManager::GetverboseLevel ( )
inline

Definition at line 422 of file G4SteppingManager.hh.

422  {
423  return verboseLevel;
424  }
Here is the caller graph for this function:

◆ InvokeAlongStepDoItProcs()

void G4SteppingManager::InvokeAlongStepDoItProcs ( )
private

Definition at line 400 of file G4SteppingManager2.cc.

402 {
403 
404 // If the current Step is defined by a 'ExclusivelyForced'
405 // PostStepDoIt, then don't invoke any AlongStepDoIt
406  if(fStepStatus == fExclusivelyForcedProc){
407  return; // Take note 'return' at here !!!
408  }
409 
410 // Invoke the all active continuous processes
411  for( size_t ci=0 ; ci<MAXofAlongStepLoops ; ci++ ){
412  fCurrentProcess = (*fAlongStepDoItVector)[ci];
413  if (fCurrentProcess== 0) continue;
414  // NULL means the process is inactivated by a user on fly.
415 
418 
419  // Update the PostStepPoint of Step according to ParticleChange
420  fParticleChange->UpdateStepForAlongStep(fStep);
421 #ifdef G4VERBOSE
422  // !!!!! Verbose
424 #endif
425 
426  // Now Store the secondaries from ParticleChange to SecondaryList
427  G4Track* tempSecondaryTrack;
428  G4int num2ndaries;
429 
430  num2ndaries = fParticleChange->GetNumberOfSecondaries();
431 
432  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
433  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
434 
435  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
436  { ApplyProductionCut(tempSecondaryTrack); }
437 
438  // Set parentID
439  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
440 
441  // Set the process pointer which created this track
442  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
443 
444  // If this 2ndry particle has 'zero' kinetic energy, make sure
445  // it invokes a rest process at the beginning of the tracking
446  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
447  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
448  if(!pm && tempSecondaryTrack->GetDefinition()->IsGeneralIon())
450  if (pm->GetAtRestProcessVector()->entries()>0){
451  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
452  fSecondary->push_back( tempSecondaryTrack );
454  } else {
455  delete tempSecondaryTrack;
456  }
457  } else {
458  fSecondary->push_back( tempSecondaryTrack );
460  }
461  } //end of loop on secondary
462 
463  // Set the track status according to what the process defined
464  // if kinetic energy >0, otherwise set fStopButAlive
465  fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
466 
467  // clear ParticleChange
468  fParticleChange->Clear();
469  }
470 
471  fStep->UpdateTrack();
472  G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
473 
474  if ( fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN ) {
475  if(MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
476  else fNewStatus = fStopAndKill;
477  fTrack->SetTrackStatus( fNewStatus );
478  }
479 
480 }
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4VSteppingVerbose * fVerbose
G4ProcessManager * GetProcessManager() const
G4ParticleDefinition * GetGenericIon() const
int G4int
Definition: G4Types.hh:78
G4int entries() const
G4VProcess * fCurrentProcess
void ApplyProductionCut(G4Track *)
static G4ParticleTable * GetParticleTable()
G4TrackVector * fSecondary
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
virtual void AlongStepDoItOneByOne()=0
#define DBL_MIN
Definition: templates.hh:75
G4VParticleChange * fParticleChange
G4StepStatus fStepStatus
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InvokeAtRestDoItProcs()

void G4SteppingManager::InvokeAtRestDoItProcs ( )
private

Definition at line 282 of file G4SteppingManager2.cc.

284 {
285 // Select the rest process which has the shortest time before
286 // it is invoked. In rest processes, GPIL()
287 // returns the time before a process occurs.
288  G4double lifeTime, shortestLifeTime;
289 
291  shortestLifeTime = DBL_MAX;
292 
293  unsigned int NofInactiveProc=0;
294  for( size_t ri=0 ; ri < MAXofAtRestLoops ; ri++ ){
295  fCurrentProcess = (*fAtRestGetPhysIntVector)[ri];
296  if (fCurrentProcess== 0) {
297  (*fSelectedAtRestDoItVector)[ri] = InActivated;
298  NofInactiveProc++;
299  continue;
300  } // NULL means the process is inactivated by a user on fly.
301 
302  lifeTime =
304 
305 
306 
307  if(fCondition==Forced && fCurrentProcess){
308  (*fSelectedAtRestDoItVector)[ri] = Forced;
309  }
310  else{
311  (*fSelectedAtRestDoItVector)[ri] = InActivated;
312  if(lifeTime < shortestLifeTime ){
313  shortestLifeTime = lifeTime;
315  }
316  }
317  }
318 
319  (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
320 
321 // at least one process is necessary to destroy the particle
322 // exit with warning
323  if(NofInactiveProc==MAXofAtRestLoops){
324  G4Exception("G4SteppingManager::InvokeAtRestDoItProcs()", "Tracking0013",
325  JustWarning, "No AtRestDoIt process is active!" );
326  }
327 
328  fStep->SetStepLength( 0. ); //the particle has stopped
329  fTrack->SetStepLength( 0. );
330 
331 // invoke selected process
332  for(size_t np=0; np < MAXofAtRestLoops; np++){
333  //
334  // Note: DoItVector has inverse order against GetPhysIntVector
335  // and SelectedAtRestDoItVector.
336  //
337  if( (*fSelectedAtRestDoItVector)[MAXofAtRestLoops-np-1] != InActivated){
338 
339  fCurrentProcess = (*fAtRestDoItVector)[np];
342 
343  // Set the current process as a process which defined this Step length
344  fStep->GetPostStepPoint()
345  ->SetProcessDefinedStep(fCurrentProcess);
346 
347  // Update Step
348  fParticleChange->UpdateStepForAtRest(fStep);
349 
350  // Now Store the secondaries from ParticleChange to SecondaryList
351  G4Track* tempSecondaryTrack;
352  G4int num2ndaries;
353 
354  num2ndaries = fParticleChange->GetNumberOfSecondaries();
355 
356  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
357  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
358 
359  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
360  { ApplyProductionCut(tempSecondaryTrack); }
361 
362  // Set parentID
363  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
364 
365  // Set the process pointer which created this track
366  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
367 
368  // If this 2ndry particle has 'zero' kinetic energy, make sure
369  // it invokes a rest process at the beginning of the tracking
370  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
371  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
372  if(!pm && tempSecondaryTrack->GetDefinition()->IsGeneralIon())
374  if (pm->GetAtRestProcessVector()->entries()>0){
375  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
376  fSecondary->push_back( tempSecondaryTrack );
378  } else {
379  delete tempSecondaryTrack;
380  }
381  } else {
382  fSecondary->push_back( tempSecondaryTrack );
384  }
385  } //end of loop on secondary
386 
387 
388  // clear ParticleChange
389  fParticleChange->Clear();
390 
391  } //if(fSelectedAtRestDoItVector[np] != InActivated){
392  } //for(size_t np=0; np < MAXofAtRestLoops; np++){
393  fStep->UpdateTrack();
394 
395  fTrack->SetTrackStatus( fStopAndKill );
396 }
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessManager * GetProcessManager() const
G4ParticleDefinition * GetGenericIon() const
int G4int
Definition: G4Types.hh:78
G4int entries() const
G4VProcess * fCurrentProcess
G4double AtRestGPIL(const G4Track &track, G4ForceCondition *condition)
Definition: G4VProcess.hh:490
void ApplyProductionCut(G4Track *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ParticleTable * GetParticleTable()
G4TrackVector * fSecondary
#define DBL_MIN
Definition: templates.hh:75
G4ForceCondition fCondition
double G4double
Definition: G4Types.hh:76
G4VParticleChange * fParticleChange
#define DBL_MAX
Definition: templates.hh:83
G4SelectedAtRestDoItVector * fSelectedAtRestDoItVector
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InvokePostStepDoItProcs()

void G4SteppingManager::InvokePostStepDoItProcs ( )
private

Definition at line 483 of file G4SteppingManager2.cc.

485 {
486 
487 // Invoke the specified discrete processes
488  for(size_t np=0; np < MAXofPostStepLoops; np++){
489  //
490  // Note: DoItVector has inverse order against GetPhysIntVector
491  // and SelectedPostStepDoItVector.
492  //
493  G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np-1];
494  if(Cond != InActivated){
495  if( ((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
496  ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
497  // ((Cond == Conditionally) && (fStepStatus == fAlongStepDoItProc)) ||
498  ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) ||
499  ((Cond == StronglyForced) )
500  ) {
501 
502  InvokePSDIP(np);
503  if ((np==0) && (fTrack->GetNextVolume() == 0)){
504  fStepStatus = fWorldBoundary;
505  fStep->GetPostStepPoint()->SetStepStatus( fStepStatus );
506  }
507  }
508  } //if(*fSelectedPostStepDoItVector(np)........
509 
510  // Exit from PostStepLoop if the track has been killed,
511  // but extra treatment for processes with Strongly Forced flag
512  if(fTrack->GetTrackStatus() == fStopAndKill) {
513  for(size_t np1=np+1; np1 < MAXofPostStepLoops; np1++){
514  G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops-np1-1];
515  if (Cond2 == StronglyForced) {
516  InvokePSDIP(np1);
517  }
518  }
519  break;
520  }
521  } //for(size_t np=0; np < MAXofPostStepLoops; np++){
522 }
int G4int
Definition: G4Types.hh:78
G4StepStatus fStepStatus
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InvokePSDIP()

void G4SteppingManager::InvokePSDIP ( size_t  np)
private

Definition at line 526 of file G4SteppingManager2.cc.

527 {
528  fCurrentProcess = (*fPostStepDoItVector)[np];
531 
532  // Update PostStepPoint of Step according to ParticleChange
533  fParticleChange->UpdateStepForPostStep(fStep);
534 #ifdef G4VERBOSE
535  // !!!!! Verbose
537 #endif
538  // Update G4Track according to ParticleChange after each PostStepDoIt
539  fStep->UpdateTrack();
540 
541  // Update safety after each invocation of PostStepDoIts
542  fStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
543 
544  // Now Store the secondaries from ParticleChange to SecondaryList
545  G4Track* tempSecondaryTrack;
546  G4int num2ndaries;
547 
548  num2ndaries = fParticleChange->GetNumberOfSecondaries();
549 
550  for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
551  tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
552 
553  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
554  { ApplyProductionCut(tempSecondaryTrack); }
555 
556  // Set parentID
557  tempSecondaryTrack->SetParentID( fTrack->GetTrackID() );
558 
559  // Set the process pointer which created this track
560  tempSecondaryTrack->SetCreatorProcess( fCurrentProcess );
561 
562  // If this 2ndry particle has 'zero' kinetic energy, make sure
563  // it invokes a rest process at the beginning of the tracking
564  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN){
565  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
566  if(!pm && tempSecondaryTrack->GetDefinition()->IsGeneralIon())
568  if (pm->GetAtRestProcessVector()->entries()>0){
569  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
570  fSecondary->push_back( tempSecondaryTrack );
572  } else {
573  delete tempSecondaryTrack;
574  }
575  } else {
576  fSecondary->push_back( tempSecondaryTrack );
578  }
579  } //end of loop on secondary
580 
581  // Set the track status according to what the process defined
582  fTrack->SetTrackStatus( fParticleChange->GetTrackStatus() );
583 
584  // clear ParticleChange
585  fParticleChange->Clear();
586 }
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
virtual void PostStepDoItOneByOne()=0
G4VSteppingVerbose * fVerbose
G4ProcessManager * GetProcessManager() const
G4ParticleDefinition * GetGenericIon() const
int G4int
Definition: G4Types.hh:78
G4int entries() const
G4VProcess * fCurrentProcess
void ApplyProductionCut(G4Track *)
static G4ParticleTable * GetParticleTable()
G4TrackVector * fSecondary
#define DBL_MIN
Definition: templates.hh:75
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4VParticleChange * fParticleChange
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetInitialStep()

void G4SteppingManager::SetInitialStep ( G4Track *  valueTrack)

Definition at line 256 of file G4SteppingManager.cc.

258 {
259 
260 // Set up several local variables.
261  PreStepPointIsGeom = false;
262  FirstStep = true;
263  fParticleChange = 0;
264  fPreviousStepSize = 0.;
265  fStepStatus = fUndefined;
266 
267  fTrack = valueTrack;
268  Mass = fTrack->GetDynamicParticle()->GetMass();
269 
270  PhysicalStep = 0.;
271  GeometricalStep = 0.;
272  CorrectedStep = 0.;
273  PreStepPointIsGeom = false;
274  FirstStep = false;
275  fStepStatus = fUndefined;
276 
277  TempInitVelocity = 0.;
278  TempVelocity = 0.;
279  sumEnergyChange = 0.;
280 
281 
282 // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
283 // set the track state to 'Alive'.
284  if( (fTrack->GetTrackStatus()==fSuspend) ||
285  (fTrack->GetTrackStatus()==fPostponeToNextEvent) ){
286  fTrack->SetTrackStatus(fAlive);
287  }
288 
289 // If the primary track has 'zero' kinetic energy, set the track
290 // state to 'StopButAlive'.
291  if(fTrack->GetKineticEnergy() <= 0.0){
292  fTrack->SetTrackStatus( fStopButAlive );
293  }
294 
295 
296 // Set Touchable to track and a private attribute of G4SteppingManager
297 
298 
299  if ( ! fTrack->GetTouchableHandle() ) {
300  G4ThreeVector direction= fTrack->GetMomentumDirection();
302  &direction, false, false );
304 
305  fTrack->SetTouchableHandle( fTouchableHandle );
306  fTrack->SetNextTouchableHandle( fTouchableHandle );
307  }else{
308  fTrack->SetNextTouchableHandle( fTouchableHandle = fTrack->GetTouchableHandle() );
309  G4VPhysicalVolume* oldTopVolume= fTrack->GetTouchableHandle()->GetVolume();
310  G4VPhysicalVolume* newTopVolume=
311  fNavigator->ResetHierarchyAndLocate( fTrack->GetPosition(),
312  fTrack->GetMomentumDirection(),
313  *((G4TouchableHistory*)fTrack->GetTouchableHandle()()) );
314 // if(newTopVolume != oldTopVolume ){
315  if(newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId() == 1 ) {
317  fTrack->SetTouchableHandle( fTouchableHandle );
318  fTrack->SetNextTouchableHandle( fTouchableHandle );
319  }
320  }
321 // Set OriginTouchableHandle for primary track
322  if(fTrack->GetParentID()==0){
323  fTrack->SetOriginTouchableHandle(fTrack->GetTouchableHandle());
324  }
325 
326 // Set vertex information of G4Track at here
327  if ( fTrack->GetCurrentStepNumber() == 0 ) {
328  fTrack->SetVertexPosition( fTrack->GetPosition() );
329  fTrack->SetVertexMomentumDirection( fTrack->GetMomentumDirection() );
330  fTrack->SetVertexKineticEnergy( fTrack->GetKineticEnergy() );
331  fTrack->SetLogicalVolumeAtVertex( fTrack->GetVolume()->GetLogicalVolume() );
332  }
333 // Initial set up for attributes of 'G4SteppingManager'
335 
336 // If track is already outside the world boundary, kill it
337  if( fCurrentVolume==0 ){
338  // If the track is a primary, stop processing
339  if(fTrack->GetParentID()==0)
340  {
341  G4cerr << "ERROR - G4SteppingManager::SetInitialStep()" << G4endl
342  << " Primary particle starting at - "
343  << fTrack->GetPosition()
344  << " - is outside of the world volume." << G4endl;
345  G4Exception("G4SteppingManager::SetInitialStep()", "Tracking0010",
346  FatalException, "Primary vertex outside of the world!");
347  }
348 
349  fTrack->SetTrackStatus( fStopAndKill );
350  G4cout << "WARNING - G4SteppingManager::SetInitialStep()" << G4endl
351  << " Initial track position is outside world! - "
352  << fTrack->GetPosition() << G4endl;
353  }
354  else {
355 // Initial set up for attribues of 'Step'
356  fStep->InitializeStep( fTrack );
357  }
358 #ifdef G4VERBOSE
359  // !!!!! Verbose
361 #endif
362 }
G4Navigator * fNavigator
G4VSteppingVerbose * fVerbose
virtual G4VPhysicalVolume * ResetHierarchyAndLocate(const G4ThreeVector &point, const G4ThreeVector &direction, const G4TouchableHistory &h)
Definition: G4Navigator.cc:95
G4GLOB_DLL std::ostream G4cout
virtual G4int GetRegularStructureId() const =0
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4VPhysicalVolume * LocateGlobalPointAndSetup(const G4ThreeVector &point, const G4ThreeVector *direction=0, const G4bool pRelativeSearch=true, const G4bool ignoreDirection=true)
Definition: G4Navigator.cc:125
G4VPhysicalVolume * fCurrentVolume
#define G4endl
Definition: G4ios.hh:61
G4TouchableHistory * CreateTouchableHistory() const
G4VParticleChange * fParticleChange
G4GLOB_DLL std::ostream G4cerr
G4StepStatus fStepStatus
G4TouchableHandle fTouchableHandle
virtual void TrackingStarted()=0
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetNavigator()

void G4SteppingManager::SetNavigator ( G4Navigator value)
inline

Definition at line 463 of file G4SteppingManager.hh.

463  {
464  fNavigator = value;
465  }
G4Navigator * fNavigator
Here is the caller graph for this function:

◆ SetUserAction()

void G4SteppingManager::SetUserAction ( G4UserSteppingAction apAction)
inline

Definition at line 467 of file G4SteppingManager.hh.

467  {
468  fUserSteppingAction = apAction;
469  }
G4UserSteppingAction * fUserSteppingAction
Here is the caller graph for this function:

◆ SetVerbose()

void G4SteppingManager::SetVerbose ( G4VSteppingVerbose yourVerbose)
inline

Definition at line 482 of file G4SteppingManager.hh.

482  {
483  fVerbose = yourVerbose;
484  }
G4VSteppingVerbose * fVerbose

◆ SetVerboseLevel()

void G4SteppingManager::SetVerboseLevel ( G4int  vLevel)
inline

Definition at line 478 of file G4SteppingManager.hh.

478  {
479  verboseLevel = vLevel;
480  }
Here is the caller graph for this function:

◆ Stepping()

G4StepStatus G4SteppingManager::Stepping ( )

Definition at line 116 of file G4SteppingManager.cc.

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

Member Data Documentation

◆ CorrectedStep

G4double G4SteppingManager::CorrectedStep
private

Definition at line 213 of file G4SteppingManager.hh.

◆ endpointSafety

G4double G4SteppingManager::endpointSafety
private

Definition at line 279 of file G4SteppingManager.hh.

◆ endpointSafOrigin

G4ThreeVector G4SteppingManager::endpointSafOrigin
private

Definition at line 278 of file G4SteppingManager.hh.

◆ fAlongStepDoItProcTriggered

size_t G4SteppingManager::fAlongStepDoItProcTriggered
private

Definition at line 251 of file G4SteppingManager.hh.

◆ fAlongStepDoItVector

G4ProcessVector* G4SteppingManager::fAlongStepDoItVector
private

Definition at line 239 of file G4SteppingManager.hh.

◆ fAlongStepGetPhysIntVector

G4ProcessVector* G4SteppingManager::fAlongStepGetPhysIntVector
private

Definition at line 243 of file G4SteppingManager.hh.

◆ fAtRestDoItProcTriggered

size_t G4SteppingManager::fAtRestDoItProcTriggered
private

Definition at line 250 of file G4SteppingManager.hh.

◆ fAtRestDoItVector

G4ProcessVector* G4SteppingManager::fAtRestDoItVector
private

Definition at line 238 of file G4SteppingManager.hh.

◆ fAtRestGetPhysIntVector

G4ProcessVector* G4SteppingManager::fAtRestGetPhysIntVector
private

Definition at line 242 of file G4SteppingManager.hh.

◆ fCondition

G4ForceCondition G4SteppingManager::fCondition
private

Definition at line 283 of file G4SteppingManager.hh.

◆ fCurrentProcess

G4VProcess* G4SteppingManager::fCurrentProcess
private

Definition at line 233 of file G4SteppingManager.hh.

◆ fCurrentVolume

G4VPhysicalVolume* G4SteppingManager::fCurrentVolume
private

Definition at line 231 of file G4SteppingManager.hh.

◆ fGPILSelection

G4GPILSelection G4SteppingManager::fGPILSelection
private

Definition at line 284 of file G4SteppingManager.hh.

◆ FirstStep

G4bool G4SteppingManager::FirstStep
private

Definition at line 215 of file G4SteppingManager.hh.

◆ fN2ndariesAlongStepDoIt

G4int G4SteppingManager::fN2ndariesAlongStepDoIt
private

Definition at line 255 of file G4SteppingManager.hh.

◆ fN2ndariesAtRestDoIt

G4int G4SteppingManager::fN2ndariesAtRestDoIt
private

Definition at line 254 of file G4SteppingManager.hh.

◆ fN2ndariesPostStepDoIt

G4int G4SteppingManager::fN2ndariesPostStepDoIt
private

Definition at line 256 of file G4SteppingManager.hh.

◆ fNavigator

G4Navigator* G4SteppingManager::fNavigator
private

Definition at line 260 of file G4SteppingManager.hh.

◆ fParticleChange

G4VParticleChange* G4SteppingManager::fParticleChange
private

Definition at line 224 of file G4SteppingManager.hh.

◆ fPostStepDoItProcTriggered

size_t G4SteppingManager::fPostStepDoItProcTriggered
private

Definition at line 252 of file G4SteppingManager.hh.

◆ fPostStepDoItVector

G4ProcessVector* G4SteppingManager::fPostStepDoItVector
private

Definition at line 240 of file G4SteppingManager.hh.

◆ fPostStepGetPhysIntVector

G4ProcessVector* G4SteppingManager::fPostStepGetPhysIntVector
private

Definition at line 244 of file G4SteppingManager.hh.

◆ fPostStepPoint

G4StepPoint* G4SteppingManager::fPostStepPoint
private

Definition at line 229 of file G4SteppingManager.hh.

◆ fPreStepPoint

G4StepPoint* G4SteppingManager::fPreStepPoint
private

Definition at line 228 of file G4SteppingManager.hh.

◆ fPreviousStepSize

G4double G4SteppingManager::fPreviousStepSize
private

Definition at line 268 of file G4SteppingManager.hh.

◆ fSecondary

G4TrackVector* G4SteppingManager::fSecondary
private

Definition at line 226 of file G4SteppingManager.hh.

◆ fSelectedAlongStepDoItVector

G4SelectedAlongStepDoItVector* G4SteppingManager::fSelectedAlongStepDoItVector
private

Definition at line 265 of file G4SteppingManager.hh.

◆ fSelectedAtRestDoItVector

G4SelectedAtRestDoItVector* G4SteppingManager::fSelectedAtRestDoItVector
private

Definition at line 264 of file G4SteppingManager.hh.

◆ fSelectedPostStepDoItVector

G4SelectedPostStepDoItVector* G4SteppingManager::fSelectedPostStepDoItVector
private

Definition at line 266 of file G4SteppingManager.hh.

◆ fSensitive

G4VSensitiveDetector* G4SteppingManager::fSensitive
private

Definition at line 232 of file G4SteppingManager.hh.

◆ fStep

G4Step* G4SteppingManager::fStep
private

Definition at line 227 of file G4SteppingManager.hh.

◆ fStepStatus

G4StepStatus G4SteppingManager::fStepStatus
private

Definition at line 216 of file G4SteppingManager.hh.

◆ fTouchableHandle

G4TouchableHandle G4SteppingManager::fTouchableHandle
private

Definition at line 270 of file G4SteppingManager.hh.

◆ fTrack

G4Track* G4SteppingManager::fTrack
private

Definition at line 225 of file G4SteppingManager.hh.

◆ fUserSteppingAction

G4UserSteppingAction* G4SteppingManager::fUserSteppingAction
private

Definition at line 207 of file G4SteppingManager.hh.

◆ fVerbose

G4VSteppingVerbose* G4SteppingManager::fVerbose
private

Definition at line 209 of file G4SteppingManager.hh.

◆ GeometricalStep

G4double G4SteppingManager::GeometricalStep
private

Definition at line 212 of file G4SteppingManager.hh.

◆ kCarTolerance

G4double G4SteppingManager::kCarTolerance
private

Definition at line 274 of file G4SteppingManager.hh.

◆ KillVerbose

G4bool G4SteppingManager::KillVerbose

Definition at line 186 of file G4SteppingManager.hh.

◆ Mass

G4double G4SteppingManager::Mass
private

Definition at line 220 of file G4SteppingManager.hh.

◆ MAXofAlongStepLoops

size_t G4SteppingManager::MAXofAlongStepLoops
private

Definition at line 247 of file G4SteppingManager.hh.

◆ MAXofAtRestLoops

size_t G4SteppingManager::MAXofAtRestLoops
private

Definition at line 246 of file G4SteppingManager.hh.

◆ MAXofPostStepLoops

size_t G4SteppingManager::MAXofPostStepLoops
private

Definition at line 248 of file G4SteppingManager.hh.

◆ PhysicalStep

G4double G4SteppingManager::PhysicalStep
private

Definition at line 211 of file G4SteppingManager.hh.

◆ physIntLength

G4double G4SteppingManager::physIntLength
private

Definition at line 282 of file G4SteppingManager.hh.

◆ PreStepPointIsGeom

G4bool G4SteppingManager::PreStepPointIsGeom
private

Definition at line 214 of file G4SteppingManager.hh.

◆ proposedSafety

G4double G4SteppingManager::proposedSafety
private

Definition at line 276 of file G4SteppingManager.hh.

◆ StepControlFlag

G4SteppingControl G4SteppingManager::StepControlFlag
private

Definition at line 272 of file G4SteppingManager.hh.

◆ sumEnergyChange

G4double G4SteppingManager::sumEnergyChange
private

Definition at line 222 of file G4SteppingManager.hh.

◆ TempInitVelocity

G4double G4SteppingManager::TempInitVelocity
private

Definition at line 218 of file G4SteppingManager.hh.

◆ TempVelocity

G4double G4SteppingManager::TempVelocity
private

Definition at line 219 of file G4SteppingManager.hh.

◆ verboseLevel

G4int G4SteppingManager::verboseLevel
private

Definition at line 262 of file G4SteppingManager.hh.


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