Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/tracking/src/G4SteppingManager.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 // G4SteppingManager class implementation
 27 //
 28 // Contact:
 29 //   Questions and comments to this code should be sent to
 30 //     Katsuya Amako  (e-mail: Katsuya.Amako@kek.jp)
 31 //     Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
 32 // --------------------------------------------------------------------
 33 
 34 #include "G4SteppingManager.hh"
 35 
 36 #include "G4ForceCondition.hh"
 37 #include "G4GPILSelection.hh"
 38 #include "G4GeometryTolerance.hh"
 39 #include "G4ParticleTable.hh"
 40 #include "G4SteppingControl.hh"
 41 #include "G4SteppingVerbose.hh"
 42 #include "G4SteppingVerboseWithUnits.hh"
 43 #include "G4TransportationManager.hh"
 44 #include "G4UImanager.hh"
 45 #include "G4UserLimits.hh"
 46 #include "G4VSensitiveDetector.hh"  // Include from 'hits/digi'
 47 
 48 // #define debug
 49 
 50 //////////////////////////////////////
 51 G4SteppingManager::G4SteppingManager()
 52 //////////////////////////////////////
 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
 62   fVerbose = G4VSteppingVerbose::GetInstance();
 63   if (fVerbose == nullptr) {
 64     if (G4VSteppingVerbose::GetMasterInstance() == nullptr) {
 65       G4int prec = G4SteppingVerbose::BestUnitPrecision();
 66       if (prec > 0) {
 67         fVerbose = new G4SteppingVerboseWithUnits(prec);
 68       }
 69       else {
 70         fVerbose = new G4SteppingVerbose();
 71       }
 72     }
 73     else {
 74       fVerbose = G4VSteppingVerbose::GetMasterInstance()->Clone();
 75     }
 76     KillVerbose = true;
 77   }
 78   else {
 79     KillVerbose = false;
 80   }
 81   fVerbose->SetManager(this);
 82 #endif
 83 
 84   SetNavigator(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking());
 85 
 86   fSelectedAtRestDoItVector = new G4SelectedAtRestDoItVector(SizeOfSelectedDoItVector, 0);
 87   fSelectedAlongStepDoItVector = new G4SelectedAlongStepDoItVector(SizeOfSelectedDoItVector, 0);
 88   fSelectedPostStepDoItVector = new G4SelectedPostStepDoItVector(SizeOfSelectedDoItVector, 0);
 89 
 90   SetNavigator(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking());
 91 
 92   physIntLength = DBL_MAX;
 93   kCarTolerance = 0.5 * G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
 94 
 95   fNoProcess = new G4NoProcess;
 96 }
 97 
 98 ///////////////////////////////////////
 99 G4SteppingManager::~G4SteppingManager()
100 ///////////////////////////////////////
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 }
118 
119 //////////////////////////////////////////
120 G4StepStatus G4SteppingManager::Stepping()
121 //////////////////////////////////////////
122 {
123 //--------
124 // Prelude
125 //--------
126 #ifdef G4VERBOSE
127   if (verboseLevel > 0) {
128     fVerbose->NewStep();
129   }
130   else if (verboseLevel == -1) {
131     G4VSteppingVerbose::SetSilent(1);
132   }
133   else {
134     G4VSteppingVerbose::SetSilent(0);
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 }
266 
267 ///////////////////////////////////////////////////////////
268 void G4SteppingManager::SetInitialStep(G4Track* valueTrack)
269 ///////////////////////////////////////////////////////////
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 }
372 
373 /////////////////////////////////////////////////
374 void G4SteppingManager::GetProcessNumber()
375 /////////////////////////////////////////////////
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 }
435 
436 // ************************************************************************
437 //
438 //  Private Member Functions
439 //
440 // ************************************************************************
441 
442 /////////////////////////////////////////////////////////
443 void G4SteppingManager::DefinePhysicalStepLength()
444 /////////////////////////////////////////////////////////
445 {
446   // ReSet the counter etc.
447   //
448   PhysicalStep = DBL_MAX;  // Initialize by a huge number
449   physIntLength = DBL_MAX;  // Initialize by a huge number
450 
451 #ifdef G4VERBOSE
452   if (verboseLevel > 0) fVerbose->DPSLStarted();
453 #endif
454 
455   // GPIL for PostStep
456   //
457   fPostStepDoItProcTriggered = MAXofPostStepLoops;
458 
459   for (std::size_t np = 0; np < MAXofPostStepLoops; ++np) {
460     fCurrentProcess = (*fPostStepGetPhysIntVector)((G4int)np);
461     if (fCurrentProcess == nullptr) {
462       (*fSelectedPostStepDoItVector)[np] = InActivated;
463       continue;
464     }  // NULL means the process is inactivated by a user on fly
465 
466     physIntLength = fCurrentProcess->PostStepGPIL(*fTrack, fPreviousStepSize, &fCondition);
467 #ifdef G4VERBOSE
468     if (verboseLevel > 0) fVerbose->DPSLPostStep();
469 #endif
470 
471     switch (fCondition) {
472       case ExclusivelyForced:
473         (*fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
474         fStepStatus = fExclusivelyForcedProc;
475         fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
476         break;
477       case Conditionally:
478         // (*fSelectedPostStepDoItVector)[np] = Conditionally;
479         G4Exception("G4SteppingManager::DefinePhysicalStepLength()", "Tracking1001", FatalException,
480           "This feature no more supported");
481         break;
482       case Forced:
483         (*fSelectedPostStepDoItVector)[np] = Forced;
484         break;
485       case StronglyForced:
486         (*fSelectedPostStepDoItVector)[np] = StronglyForced;
487         break;
488       default:
489         (*fSelectedPostStepDoItVector)[np] = InActivated;
490         break;
491     }
492 
493     if (fCondition == ExclusivelyForced) {
494       for (std::size_t nrest = np + 1; nrest < MAXofPostStepLoops; ++nrest) {
495         (*fSelectedPostStepDoItVector)[nrest] = InActivated;
496       }
497       return;  // Take note the 'return' at here !!!
498     }
499 
500     if (physIntLength < PhysicalStep) {
501       PhysicalStep = physIntLength;
502       fStepStatus = fPostStepDoItProc;
503       fPostStepDoItProcTriggered = G4int(np);
504       fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
505     }
506   }
507 
508   if (fPostStepDoItProcTriggered < MAXofPostStepLoops) {
509     if ((*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] == InActivated) {
510       (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = NotForced;
511     }
512   }
513 
514   // GPIL for AlongStep
515   //
516   proposedSafety = DBL_MAX;
517   G4double safetyProposedToAndByProcess = proposedSafety;
518   G4bool delegateToTransportation = false;
519 
520   for (std::size_t kp = 0; kp < MAXofAlongStepLoops; ++kp) {
521     fCurrentProcess = (*fAlongStepGetPhysIntVector)[(G4int)kp];
522     if (fCurrentProcess == nullptr) continue;
523     // NULL means the process is inactivated by a user on fly
524 
525     physIntLength = fCurrentProcess->AlongStepGPIL(
526       *fTrack, fPreviousStepSize, PhysicalStep, safetyProposedToAndByProcess, &fGPILSelection);
527 #ifdef G4VERBOSE
528     if (verboseLevel > 0) fVerbose->DPSLAlongStep();
529 #endif
530 
531     if (physIntLength < PhysicalStep) {
532       PhysicalStep = physIntLength;
533 
534       // Check if the process wants to be the GPIL winner. For example,
535       // multi-scattering proposes Step limit, but won't be the winner
536       //
537       if (fGPILSelection == CandidateForSelection) {
538         fStepStatus = fAlongStepDoItProc;
539         fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
540       }
541       else if (fCurrentProcess->GetProcessType() == fParallel)
542       {  // a parallel world is proposing the shortest but expecting Transportation
543         // to win.
544         delegateToTransportation = true;
545       }
546 
547       // Transportation is assumed to be the last process in the vector
548       // Transportation is winning
549       if (kp == MAXofAlongStepLoops - 1) {
550         // This used to set fStepStatus = fGeomBoundary, but it was moved to
551         // G4Transportation::AlongStepDoIt where the process can actually
552         // decide if there is a volume boundary.
553         delegateToTransportation = false;
554       }
555     }
556 
557     // Make sure to check the safety, even if Step is not limited
558     // by this process
559     //
560     if (safetyProposedToAndByProcess < proposedSafety) {
561       // proposedSafety keeps the smallest value
562       //
563       proposedSafety = safetyProposedToAndByProcess;
564     }
565     else {
566       // safetyProposedToAndByProcess always proposes a valid safety
567       //
568       safetyProposedToAndByProcess = proposedSafety;
569     }
570   }
571   if (delegateToTransportation) {
572     fStepStatus = fGeomBoundary;
573     fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
574   }
575 }
576 
577 //////////////////////////////////////////////////////
578 G4int G4SteppingManager::ProcessSecondariesFromParticleChange()
579 //////////////////////////////////////////////////////
580 {
581   G4Track* tempSecondaryTrack;
582   G4int num2ndaries;
583   G4int pushedSecondaries = 0;
584 
585   num2ndaries = fParticleChange->GetNumberOfSecondaries();
586   if (num2ndaries == 0) {
587     return 0;
588   }
589 
590   // Get the creator process. This may be different from fCurrentProcess for a
591   // "combined" process such as G4GammaGeneralProcess.
592   const G4VProcess* creatorProcess = fCurrentProcess->GetCreatorProcess();
593 
594   for (G4int DSecLoop = 0; DSecLoop < num2ndaries; ++DSecLoop) {
595     tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
596 
597     // Set parentID
598     tempSecondaryTrack->SetParentID(fTrack->GetTrackID());
599 
600     // Set the process pointer which created this track
601     tempSecondaryTrack->SetCreatorProcess(creatorProcess);
602 
603     // If this 2ndry particle has 'zero' kinetic energy, make sure
604     // it invokes a rest process at the beginning of the tracking
605     //
606     if (tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN) {
607       G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
608       if (pm == nullptr) {
609         G4ExceptionDescription ED;
610         ED << "A track without proper process manager is pushed\n"
611            << "into the track stack.\n"
612            << " Particle name : " << tempSecondaryTrack->GetDefinition()->GetParticleName()
613            << " -- created by " << creatorProcess->GetProcessName() << ".";
614         G4Exception("G4SteppingManager::ProcessSecondariesFromParticleChange()", "Tracking10051",
615           FatalException, ED);
616       }
617       if (pm->GetAtRestProcessVector()->entries() > 0) {
618         tempSecondaryTrack->SetTrackStatus(fStopButAlive);
619         fSecondary->push_back(tempSecondaryTrack);
620         ++pushedSecondaries;
621       }
622       else {
623         delete tempSecondaryTrack;
624       }
625     }
626     else {
627       fSecondary->push_back(tempSecondaryTrack);
628       ++pushedSecondaries;
629     }
630   }  // end of loop on secondary
631 
632   return pushedSecondaries;
633 }
634 
635 //////////////////////////////////////////////////////
636 void G4SteppingManager::InvokeAtRestDoItProcs()
637 //////////////////////////////////////////////////////
638 {
639   // Select the rest process which has the shortest time before
640   // it is invoked. In rest processes, GPIL()
641   // returns the time before a process occurs
642 
643   G4double lifeTime, shortestLifeTime;
644 
645   fAtRestDoItProcTriggered = 0;
646   shortestLifeTime = DBL_MAX;
647 
648   for (std::size_t ri = 0; ri < MAXofAtRestLoops; ++ri) {
649     fCurrentProcess = (*fAtRestGetPhysIntVector)[(G4int)ri];
650     if (fCurrentProcess == nullptr) {
651       (*fSelectedAtRestDoItVector)[ri] = InActivated;
652       continue;
653     }  // nullptr means the process is inactivated by a user on fly
654 
655     lifeTime = fCurrentProcess->AtRestGPIL(*fTrack, &fCondition);
656 
657     if (fCondition == Forced) {
658       (*fSelectedAtRestDoItVector)[ri] = Forced;
659     }
660     else {
661       (*fSelectedAtRestDoItVector)[ri] = InActivated;
662       if (lifeTime < shortestLifeTime) {
663         shortestLifeTime = lifeTime;
664         fAtRestDoItProcTriggered = G4int(ri);
665         fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
666       }
667     }
668   }
669 
670   (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
671 
672   fStep->SetStepLength(0.);  // the particle has stopped
673   fTrack->SetStepLength(0.);
674 
675   // Condition to avoid that stable ions are handled by Radioactive Decay.
676   // We use a very large time threshold (many orders of magnitude bigger than
677   // the universe's age) but not DBL_MAX because shortestLifeTime can be
678   // sometimes slightly smaller for stable ions.
679   if (shortestLifeTime < 1.0e+100)  // Unstable ion at rest: Radioactive Decay will decay it
680   {
681     // invoke selected process
682     //
683     for (std::size_t np = 0; np < MAXofAtRestLoops; ++np) {
684       //
685       // Note: DoItVector has inverse order against GetPhysIntVector
686       //       and SelectedAtRestDoItVector.
687       //
688       if ((*fSelectedAtRestDoItVector)[MAXofAtRestLoops - np - 1] != InActivated) {
689         fCurrentProcess = (*fAtRestDoItVector)[(G4int)np];
690         fParticleChange = fCurrentProcess->AtRestDoIt(*fTrack, *fStep);
691 
692         // Update Step
693         //
694         fParticleChange->UpdateStepForAtRest(fStep);
695 
696         // Now Store the secondaries from ParticleChange to SecondaryList
697         fN2ndariesAtRestDoIt += ProcessSecondariesFromParticleChange();
698 
699         // clear ParticleChange
700         fParticleChange->Clear();
701 
702       }  // if(fSelectedAtRestDoItVector[np] != InActivated){
703     }  // for(std::size_t np=0; np<MAXofAtRestLoops; ++np){
704   }
705   else  // Stable ion at rest
706   {
707     fStep->GetPostStepPoint()->SetProcessDefinedStep(fNoProcess);
708   }  // if(shortestLifeTime < 1.0e+100)
709 
710   fStep->UpdateTrack();
711 
712   fTrack->SetTrackStatus(fStopAndKill);
713 }
714 
715 /////////////////////////////////////////////////////////
716 void G4SteppingManager::InvokeAlongStepDoItProcs()
717 /////////////////////////////////////////////////////////
718 {
719   // If the current Step is defined by a 'ExclusivelyForced'
720   // PostStepDoIt, then don't invoke any AlongStepDoIt
721   //
722   if (fStepStatus == fExclusivelyForcedProc) {
723     return;  // Take note 'return' is here !!!
724   }
725 
726   // Invoke all active continuous processes
727   //
728   for (std::size_t ci = 0; ci < MAXofAlongStepLoops; ++ci) {
729     fCurrentProcess = (*fAlongStepDoItVector)[(G4int)ci];
730     if (fCurrentProcess == nullptr) continue;
731     // NULL means the process is inactivated by a user on fly.
732 
733     fParticleChange = fCurrentProcess->AlongStepDoIt(*fTrack, *fStep);
734 
735     // Update the PostStepPoint of Step according to ParticleChange
736     fParticleChange->UpdateStepForAlongStep(fStep);
737 
738 #ifdef G4VERBOSE
739     if (verboseLevel > 0) fVerbose->AlongStepDoItOneByOne();
740 #endif
741 
742     // Now Store the secondaries from ParticleChange to SecondaryList
743     fN2ndariesAlongStepDoIt += ProcessSecondariesFromParticleChange();
744 
745     // Set the track status according to what the process defined
746     // if kinetic energy >0, otherwise set  fStopButAlive
747     //
748     fTrack->SetTrackStatus(fParticleChange->GetTrackStatus());
749 
750     // clear ParticleChange
751     fParticleChange->Clear();
752   }
753 
754   fStep->UpdateTrack();
755   G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
756 
757   if (fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN) {
758     if (MAXofAtRestLoops > 0)
759       fNewStatus = fStopButAlive;
760     else
761       fNewStatus = fStopAndKill;
762     fTrack->SetTrackStatus(fNewStatus);
763   }
764 }
765 
766 ////////////////////////////////////////////////////////
767 void G4SteppingManager::InvokePostStepDoItProcs()
768 ////////////////////////////////////////////////////////
769 {
770   // Invoke the specified discrete processes
771   //
772   for (std::size_t np = 0; np < MAXofPostStepLoops; ++np) {
773     //
774     // Note: DoItVector has inverse order against GetPhysIntVector
775     //       and SelectedPostStepDoItVector.
776     //
777     G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops - np - 1];
778     if (Cond != InActivated) {
779       if (((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
780           ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
781           ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) ||
782           ((Cond == StronglyForced)))
783       {
784         InvokePSDIP(np);
785         if ((np == 0) && (fTrack->GetNextVolume() == nullptr)) {
786           fStepStatus = fWorldBoundary;
787           fStep->GetPostStepPoint()->SetStepStatus(fStepStatus);
788         }
789       }
790     }
791 
792     // Exit from PostStepLoop if the track has been killed,
793     // but extra treatment for processes with Strongly Forced flag
794     //
795     if (fTrack->GetTrackStatus() == fStopAndKill) {
796       for (std::size_t np1 = np + 1; np1 < MAXofPostStepLoops; ++np1) {
797         G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops - np1 - 1];
798         if (Cond2 == StronglyForced) {
799           InvokePSDIP(np1);
800         }
801       }
802       break;
803     }
804   }
805 }
806 
807 ////////////////////////////////////////////////////////
808 void G4SteppingManager::InvokePSDIP(size_t np)
809 ////////////////////////////////////////////////////////
810 {
811   fCurrentProcess = (*fPostStepDoItVector)[(G4int)np];
812   fParticleChange = fCurrentProcess->PostStepDoIt(*fTrack, *fStep);
813 
814   // Update PostStepPoint of Step according to ParticleChange
815   fParticleChange->UpdateStepForPostStep(fStep);
816 
817 #ifdef G4VERBOSE
818   if (verboseLevel > 0) fVerbose->PostStepDoItOneByOne();
819 #endif
820 
821   // Update G4Track according to ParticleChange after each PostStepDoIt
822   fStep->UpdateTrack();
823 
824   // Update safety after each invocation of PostStepDoIts
825   fStep->GetPostStepPoint()->SetSafety(CalculateSafety());
826 
827   // Now Store the secondaries from ParticleChange to SecondaryList
828   fN2ndariesPostStepDoIt += ProcessSecondariesFromParticleChange();
829 
830   // Set the track status according to what the process defined
831   fTrack->SetTrackStatus(fParticleChange->GetTrackStatus());
832 
833   // clear ParticleChange
834   fParticleChange->Clear();
835 }
836