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 ]

  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 "G4Profiler.hh"
 41 #include "G4SteppingControl.hh"
 42 #include "G4SteppingVerbose.hh"
 43 #include "G4SteppingVerboseWithUnits.hh"
 44 #include "G4TiMemory.hh"
 45 #include "G4TransportationManager.hh"
 46 #include "G4UImanager.hh"
 47 #include "G4UserLimits.hh"
 48 #include "G4VSensitiveDetector.hh"  // Include from 'hits/digi'
 49 
 50 // #define debug
 51 
 52 //////////////////////////////////////
 53 G4SteppingManager::G4SteppingManager()
 54 //////////////////////////////////////
 55 {
 56   // Construct simple 'has-a' related objects
 57 
 58   fStep = new G4Step();
 59   fSecondary = fStep->NewSecondaryVector();
 60   fPreStepPoint = fStep->GetPreStepPoint();
 61   fPostStepPoint = fStep->GetPostStepPoint();
 62 
 63 #ifdef G4VERBOSE
 64   fVerbose = G4VSteppingVerbose::GetInstance();
 65   if (fVerbose == nullptr) {
 66     if (G4VSteppingVerbose::GetMasterInstance() == nullptr) {
 67       G4int prec = G4SteppingVerbose::BestUnitPrecision();
 68       if (prec > 0) {
 69         fVerbose = new G4SteppingVerboseWithUnits(prec);
 70       }
 71       else {
 72         fVerbose = new G4SteppingVerbose();
 73       }
 74     }
 75     else {
 76       fVerbose = G4VSteppingVerbose::GetMasterInstance()->Clone();
 77     }
 78     KillVerbose = true;
 79   }
 80   else {
 81     KillVerbose = false;
 82   }
 83   fVerbose->SetManager(this);
 84 #endif
 85 
 86   SetNavigator(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking());
 87 
 88   fSelectedAtRestDoItVector = new G4SelectedAtRestDoItVector(SizeOfSelectedDoItVector, 0);
 89   fSelectedAlongStepDoItVector = new G4SelectedAlongStepDoItVector(SizeOfSelectedDoItVector, 0);
 90   fSelectedPostStepDoItVector = new G4SelectedPostStepDoItVector(SizeOfSelectedDoItVector, 0);
 91 
 92   SetNavigator(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking());
 93 
 94   physIntLength = DBL_MAX;
 95   kCarTolerance = 0.5 * G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
 96 
 97   fNoProcess = new G4NoProcess;
 98 }
 99 
100 ///////////////////////////////////////
101 G4SteppingManager::~G4SteppingManager()
102 ///////////////////////////////////////
103 {
104   fTouchableHandle = nullptr;
105 
106   // Destruct simple 'has-a' objects
107   //
108   fStep->DeleteSecondaryVector();
109 
110   // delete fSecondary;
111   delete fStep;
112   delete fSelectedAtRestDoItVector;
113   delete fSelectedAlongStepDoItVector;
114   delete fSelectedPostStepDoItVector;
115   delete fUserSteppingAction;
116 #ifdef G4VERBOSE
117   if (KillVerbose) delete fVerbose;
118 #endif
119 }
120 
121 //////////////////////////////////////////
122 G4StepStatus G4SteppingManager::Stepping()
123 //////////////////////////////////////////
124 {
125 #ifdef GEANT4_USE_TIMEMORY
126   ProfilerConfig profiler{fStep};
127 #endif
128 
129 //--------
130 // Prelude
131 //--------
132 #ifdef G4VERBOSE
133   if (verboseLevel > 0) {
134     fVerbose->NewStep();
135   }
136   else if (verboseLevel == -1) {
137     G4VSteppingVerbose::SetSilent(1);
138   }
139   else {
140     G4VSteppingVerbose::SetSilent(0);
141   }
142 #endif
143 
144   // Store last PostStepPoint to PreStepPoint, and swap current and nex
145   // volume information of G4Track. Reset total energy deposit in one Step.
146   //
147   fStep->CopyPostToPreStepPoint();
148   fStep->ResetTotalEnergyDeposit();
149 
150   // Switch next touchable in track to current one
151   //
152   fTrack->SetTouchableHandle(fTrack->GetNextTouchableHandle());
153 
154   // Reset the secondary particles
155   //
156   fN2ndariesAtRestDoIt = 0;
157   fN2ndariesAlongStepDoIt = 0;
158   fN2ndariesPostStepDoIt = 0;
159 
160   // Set the volume before it is used (in DefineStepLength() for User Limit)
161   //
162   fCurrentVolume = fStep->GetPreStepPoint()->GetPhysicalVolume();
163 
164   // Reset the step's auxiliary points vector pointer
165   //
166   fStep->SetPointerToVectorOfAuxiliaryPoints(nullptr);
167 
168   //-----------------
169   // AtRest Processes
170   //-----------------
171 
172   if (fTrack->GetTrackStatus() == fStopButAlive) {
173     if (MAXofAtRestLoops > 0) {
174       InvokeAtRestDoItProcs();
175       fStepStatus = fAtRestDoItProc;
176       fStep->GetPostStepPoint()->SetStepStatus(fStepStatus);
177 
178 #ifdef G4VERBOSE
179       if (verboseLevel > 0) fVerbose->AtRestDoItInvoked();
180 #endif
181     }
182     // Make sure the track is killed
183     //
184     fTrack->SetTrackStatus(fStopAndKill);
185   }
186 
187   //---------------------------------
188   // AlongStep and PostStep Processes
189   //---------------------------------
190 
191   else {
192     // Find minimum Step length demanded by active disc./cont. processes
193     DefinePhysicalStepLength();
194 
195     // Store the Step length (geometrical length) to G4Step and G4Track
196     fStep->SetStepLength(PhysicalStep);
197     fTrack->SetStepLength(PhysicalStep);
198     G4double GeomStepLength = PhysicalStep;
199 
200     // Store StepStatus to PostStepPoint
201     fStep->GetPostStepPoint()->SetStepStatus(fStepStatus);
202 
203     // Invoke AlongStepDoIt
204     InvokeAlongStepDoItProcs();
205 
206     // Get StepStatus from PostStepPoint - a process such as transportation
207     // might have changed it.
208     fStepStatus = fStep->GetPostStepPoint()->GetStepStatus();
209 
210     // Update track by taking into account all changes by AlongStepDoIt
211     fStep->UpdateTrack();
212 
213     // Update safety after invocation of all AlongStepDoIts
214     endpointSafOrigin = fPostStepPoint->GetPosition();
215     // endpointSafety=  std::max( proposedSafety - GeomStepLength, 0.);
216     endpointSafety = std::max(proposedSafety - GeomStepLength, kCarTolerance);
217 
218     fStep->GetPostStepPoint()->SetSafety(endpointSafety);
219 
220 #ifdef G4VERBOSE
221     if (verboseLevel > 0) fVerbose->AlongStepDoItAllDone();
222 #endif
223 
224     // Invoke PostStepDoIt
225     InvokePostStepDoItProcs();
226 
227 #ifdef G4VERBOSE
228     if (verboseLevel > 0) fVerbose->PostStepDoItAllDone();
229 #endif
230   }
231 
232   //-------
233   // Finale
234   //-------
235 
236   // Update 'TrackLength' and remeber the Step length of the current Step
237   //
238   fTrack->AddTrackLength(fStep->GetStepLength());
239   fPreviousStepSize = fStep->GetStepLength();
240   fStep->SetTrack(fTrack);
241 
242 #ifdef G4VERBOSE
243   if (verboseLevel > 0) fVerbose->StepInfo();
244 #endif
245 
246   // Send G4Step information to Hit/Dig if the volume is sensitive
247   //
248   fCurrentVolume = fStep->GetPreStepPoint()->GetPhysicalVolume();
249   StepControlFlag = fStep->GetControlFlag();
250   if (fCurrentVolume != nullptr && StepControlFlag != AvoidHitInvocation) {
251     fSensitive = fStep->GetPreStepPoint()->GetSensitiveDetector();
252     if (fSensitive != nullptr) {
253       fSensitive->Hit(fStep);
254     }
255   }
256 
257   // User intervention process
258   //
259   if (fUserSteppingAction != nullptr) {
260     fUserSteppingAction->UserSteppingAction(fStep);
261   }
262 
263   G4UserSteppingAction* regionalAction =
264     fCurrentVolume->GetLogicalVolume()->GetRegion()->GetRegionalSteppingAction();
265 
266   if (regionalAction != nullptr) regionalAction->UserSteppingAction(fStep);
267 
268   // Stepping process finish. Return the value of the StepStatus
269   //
270   return fStepStatus;
271 }
272 
273 ///////////////////////////////////////////////////////////
274 void G4SteppingManager::SetInitialStep(G4Track* valueTrack)
275 ///////////////////////////////////////////////////////////
276 {
277   // Set up several local variables
278   //
279   PreStepPointIsGeom = false;
280   FirstStep = true;
281   fParticleChange = nullptr;
282   fPreviousStepSize = 0.;
283   fStepStatus = fUndefined;
284 
285   fTrack = valueTrack;
286   Mass = fTrack->GetDynamicParticle()->GetMass();
287 
288   PhysicalStep = 0.;
289   GeometricalStep = 0.;
290   CorrectedStep = 0.;
291   PreStepPointIsGeom = false;
292   FirstStep = false;
293 
294   TempInitVelocity = 0.;
295   TempVelocity = 0.;
296   sumEnergyChange = 0.;
297 
298   // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
299   // set the track state to 'Alive'
300   //
301   if ((fTrack->GetTrackStatus() == fSuspend) || (fTrack->GetTrackStatus() == fPostponeToNextEvent))
302   {
303     fTrack->SetTrackStatus(fAlive);
304   }
305 
306   // If the primary track has 'zero' kinetic energy, set the track
307   // state to 'StopButAlive'
308   //
309   if (fTrack->GetKineticEnergy() <= 0.0) {
310     fTrack->SetTrackStatus(fStopButAlive);
311   }
312 
313   // Set Touchable to track and a private attribute of G4SteppingManager
314 
315   if (! fTrack->GetTouchableHandle()) {
316     G4ThreeVector direction = fTrack->GetMomentumDirection();
317     fNavigator->LocateGlobalPointAndSetup(fTrack->GetPosition(), &direction, false, false);
318     fTouchableHandle = fNavigator->CreateTouchableHistory();
319     fTrack->SetTouchableHandle(fTouchableHandle);
320     fTrack->SetNextTouchableHandle(fTouchableHandle);
321   }
322   else {
323     fTrack->SetNextTouchableHandle(fTouchableHandle = fTrack->GetTouchableHandle());
324     G4VPhysicalVolume* oldTopVolume = fTrack->GetTouchableHandle()->GetVolume();
325     G4VPhysicalVolume* newTopVolume = fNavigator->ResetHierarchyAndLocate(fTrack->GetPosition(),
326       fTrack->GetMomentumDirection(), *((G4TouchableHistory*)fTrack->GetTouchableHandle()()));
327     if (newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId() == 1) {
328       fTouchableHandle = fNavigator->CreateTouchableHistory();
329       fTrack->SetTouchableHandle(fTouchableHandle);
330       fTrack->SetNextTouchableHandle(fTouchableHandle);
331     }
332   }
333 
334   // Set OriginTouchableHandle for primary track
335   //
336   if (fTrack->GetParentID() == 0) {
337     fTrack->SetOriginTouchableHandle(fTrack->GetTouchableHandle());
338   }
339 
340   // Set vertex information of G4Track at here
341   //
342   if (fTrack->GetCurrentStepNumber() == 0) {
343     fTrack->SetVertexPosition(fTrack->GetPosition());
344     fTrack->SetVertexMomentumDirection(fTrack->GetMomentumDirection());
345     fTrack->SetVertexKineticEnergy(fTrack->GetKineticEnergy());
346     fTrack->SetLogicalVolumeAtVertex(fTrack->GetVolume()->GetLogicalVolume());
347   }
348 
349   // Initial set up for attributes of 'G4SteppingManager'
350   fCurrentVolume = fTouchableHandle->GetVolume();
351 
352   // If track is already outside the world boundary, kill it
353   //
354   if (fCurrentVolume == nullptr) {
355     // If the track is a primary, stop processing
356     if (fTrack->GetParentID() == 0) {
357       G4cerr << "ERROR - G4SteppingManager::SetInitialStep()" << G4endl
358              << "        Primary particle starting at - " << fTrack->GetPosition()
359              << " - is outside of the world volume." << G4endl;
360       G4Exception("G4SteppingManager::SetInitialStep()", "Tracking0010", FatalException,
361         "Primary vertex outside of the world!");
362     }
363 
364     fTrack->SetTrackStatus(fStopAndKill);
365     G4cout << "WARNING - G4SteppingManager::SetInitialStep()" << G4endl
366            << "          Initial track position is outside world! - " << fTrack->GetPosition()
367            << G4endl;
368   }
369   else {
370     // Initial set up for attributes of 'Step'
371     fStep->InitializeStep(fTrack);
372   }
373 
374 #ifdef G4VERBOSE
375   if (verboseLevel > 0) fVerbose->TrackingStarted();
376 #endif
377 }
378 
379 /////////////////////////////////////////////////
380 void G4SteppingManager::GetProcessNumber()
381 /////////////////////////////////////////////////
382 {
383 #ifdef debug
384   G4cout << "G4SteppingManager::GetProcessNumber: is called track=" << fTrack << G4endl;
385 #endif
386 
387   G4ProcessManager* pm = fTrack->GetDefinition()->GetProcessManager();
388   if (pm == nullptr) {
389     G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
390            << "        ProcessManager is NULL for particle = "
391            << fTrack->GetDefinition()->GetParticleName()
392            << ", PDG_code = " << fTrack->GetDefinition()->GetPDGEncoding() << G4endl;
393     G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0011", FatalException,
394       "Process Manager is not found.");
395     return;
396   }
397 
398   // AtRestDoits
399   //
400   MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
401   fAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
402   fAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
403 
404 #ifdef debug
405   G4cout << "G4SteppingManager::GetProcessNumber: #ofAtRest=" << MAXofAtRestLoops << G4endl;
406 #endif
407 
408   // AlongStepDoits
409   //
410   MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
411   fAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
412   fAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
413 
414 #ifdef debug
415   G4cout << "G4SteppingManager::GetProcessNumber:#ofAlongStp=" << MAXofAlongStepLoops << G4endl;
416 #endif
417 
418   // PostStepDoits
419   //
420   MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
421   fPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
422   fPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
423 
424 #ifdef debug
425   G4cout << "G4SteppingManager::GetProcessNumber: #ofPostStep=" << MAXofPostStepLoops << G4endl;
426 #endif
427 
428   if (SizeOfSelectedDoItVector < MAXofAtRestLoops ||
429       SizeOfSelectedDoItVector < MAXofAlongStepLoops ||
430       SizeOfSelectedDoItVector < MAXofPostStepLoops)
431   {
432     G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
433            << "        SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
434            << " ; is smaller then one of MAXofAtRestLoops= " << MAXofAtRestLoops << G4endl
435            << "        or MAXofAlongStepLoops= " << MAXofAlongStepLoops
436            << " or MAXofPostStepLoops= " << MAXofPostStepLoops << G4endl;
437     G4Exception("G4SteppingManager::GetProcessNumber()", "Tracking0012", FatalException,
438       "The array size is smaller than the actual No of processes.");
439   }
440 }
441 
442 // ************************************************************************
443 //
444 //  Private Member Functions
445 //
446 // ************************************************************************
447 
448 /////////////////////////////////////////////////////////
449 void G4SteppingManager::DefinePhysicalStepLength()
450 /////////////////////////////////////////////////////////
451 {
452   // ReSet the counter etc.
453   //
454   PhysicalStep = DBL_MAX;  // Initialize by a huge number
455   physIntLength = DBL_MAX;  // Initialize by a huge number
456 
457 #ifdef G4VERBOSE
458   if (verboseLevel > 0) fVerbose->DPSLStarted();
459 #endif
460 
461   // GPIL for PostStep
462   //
463   fPostStepDoItProcTriggered = MAXofPostStepLoops;
464 
465   for (std::size_t np = 0; np < MAXofPostStepLoops; ++np) {
466     fCurrentProcess = (*fPostStepGetPhysIntVector)((G4int)np);
467     if (fCurrentProcess == nullptr) {
468       (*fSelectedPostStepDoItVector)[np] = InActivated;
469       continue;
470     }  // NULL means the process is inactivated by a user on fly
471 
472     physIntLength = fCurrentProcess->PostStepGPIL(*fTrack, fPreviousStepSize, &fCondition);
473 #ifdef G4VERBOSE
474     if (verboseLevel > 0) fVerbose->DPSLPostStep();
475 #endif
476 
477     switch (fCondition) {
478       case ExclusivelyForced:
479         (*fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
480         fStepStatus = fExclusivelyForcedProc;
481         fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
482         break;
483       case Conditionally:
484         // (*fSelectedPostStepDoItVector)[np] = Conditionally;
485         G4Exception("G4SteppingManager::DefinePhysicalStepLength()", "Tracking1001", FatalException,
486           "This feature no more supported");
487         break;
488       case Forced:
489         (*fSelectedPostStepDoItVector)[np] = Forced;
490         break;
491       case StronglyForced:
492         (*fSelectedPostStepDoItVector)[np] = StronglyForced;
493         break;
494       default:
495         (*fSelectedPostStepDoItVector)[np] = InActivated;
496         break;
497     }
498 
499     if (fCondition == ExclusivelyForced) {
500       for (std::size_t nrest = np + 1; nrest < MAXofPostStepLoops; ++nrest) {
501         (*fSelectedPostStepDoItVector)[nrest] = InActivated;
502       }
503       return;  // Take note the 'return' at here !!!
504     }
505 
506     if (physIntLength < PhysicalStep) {
507       PhysicalStep = physIntLength;
508       fStepStatus = fPostStepDoItProc;
509       fPostStepDoItProcTriggered = G4int(np);
510       fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
511     }
512   }
513 
514   if (fPostStepDoItProcTriggered < MAXofPostStepLoops) {
515     if ((*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] == InActivated) {
516       (*fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = NotForced;
517     }
518   }
519 
520   // GPIL for AlongStep
521   //
522   proposedSafety = DBL_MAX;
523   G4double safetyProposedToAndByProcess = proposedSafety;
524   G4bool delegateToTransportation = false;
525 
526   for (std::size_t kp = 0; kp < MAXofAlongStepLoops; ++kp) {
527     fCurrentProcess = (*fAlongStepGetPhysIntVector)[(G4int)kp];
528     if (fCurrentProcess == nullptr) continue;
529     // NULL means the process is inactivated by a user on fly
530 
531     physIntLength = fCurrentProcess->AlongStepGPIL(
532       *fTrack, fPreviousStepSize, PhysicalStep, safetyProposedToAndByProcess, &fGPILSelection);
533 #ifdef G4VERBOSE
534     if (verboseLevel > 0) fVerbose->DPSLAlongStep();
535 #endif
536 
537     if (physIntLength < PhysicalStep) {
538       PhysicalStep = physIntLength;
539 
540       // Check if the process wants to be the GPIL winner. For example,
541       // multi-scattering proposes Step limit, but won't be the winner
542       //
543       if (fGPILSelection == CandidateForSelection) {
544         fStepStatus = fAlongStepDoItProc;
545         fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
546       }
547       else if (fCurrentProcess->GetProcessType() == fParallel)
548       {  // a parallel world is proposing the shortest but expecting Transportation
549         // to win.
550         delegateToTransportation = true;
551       }
552 
553       // Transportation is assumed to be the last process in the vector
554       // Transportation is winning
555       if (kp == MAXofAlongStepLoops - 1) {
556         // This used to set fStepStatus = fGeomBoundary, but it was moved to
557         // G4Transportation::AlongStepDoIt where the process can actually
558         // decide if there is a volume boundary.
559         delegateToTransportation = false;
560       }
561     }
562 
563     // Make sure to check the safety, even if Step is not limited
564     // by this process
565     //
566     if (safetyProposedToAndByProcess < proposedSafety) {
567       // proposedSafety keeps the smallest value
568       //
569       proposedSafety = safetyProposedToAndByProcess;
570     }
571     else {
572       // safetyProposedToAndByProcess always proposes a valid safety
573       //
574       safetyProposedToAndByProcess = proposedSafety;
575     }
576   }
577   if (delegateToTransportation) {
578     fStepStatus = fGeomBoundary;
579     fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
580   }
581 }
582 
583 //////////////////////////////////////////////////////
584 G4int G4SteppingManager::ProcessSecondariesFromParticleChange()
585 //////////////////////////////////////////////////////
586 {
587   G4Track* tempSecondaryTrack;
588   G4int num2ndaries;
589   G4int pushedSecondaries = 0;
590 
591   num2ndaries = fParticleChange->GetNumberOfSecondaries();
592   if (num2ndaries == 0) {
593     return 0;
594   }
595 
596   // Get the creator process. This may be different from fCurrentProcess for a
597   // "combined" process such as G4GammaGeneralProcess.
598   const G4VProcess* creatorProcess = fCurrentProcess->GetCreatorProcess();
599 
600   for (G4int DSecLoop = 0; DSecLoop < num2ndaries; ++DSecLoop) {
601     tempSecondaryTrack = fParticleChange->GetSecondary(DSecLoop);
602 
603     // Set parentID
604     tempSecondaryTrack->SetParentID(fTrack->GetTrackID());
605 
606     // Set the process pointer which created this track
607     tempSecondaryTrack->SetCreatorProcess(creatorProcess);
608 
609     // If this 2ndry particle has 'zero' kinetic energy, make sure
610     // it invokes a rest process at the beginning of the tracking
611     //
612     if (tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN) {
613       G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
614       if (pm == nullptr) {
615         G4ExceptionDescription ED;
616         ED << "A track without proper process manager is pushed\n"
617            << "into the track stack.\n"
618            << " Particle name : " << tempSecondaryTrack->GetDefinition()->GetParticleName()
619            << " -- created by " << creatorProcess->GetProcessName() << ".";
620         G4Exception("G4SteppingManager::ProcessSecondariesFromParticleChange()", "Tracking10051",
621           FatalException, ED);
622       }
623       if (pm->GetAtRestProcessVector()->entries() > 0) {
624         tempSecondaryTrack->SetTrackStatus(fStopButAlive);
625         fSecondary->push_back(tempSecondaryTrack);
626         ++pushedSecondaries;
627       }
628       else {
629         delete tempSecondaryTrack;
630       }
631     }
632     else {
633       fSecondary->push_back(tempSecondaryTrack);
634       ++pushedSecondaries;
635     }
636   }  // end of loop on secondary
637 
638   return pushedSecondaries;
639 }
640 
641 //////////////////////////////////////////////////////
642 void G4SteppingManager::InvokeAtRestDoItProcs()
643 //////////////////////////////////////////////////////
644 {
645   // Select the rest process which has the shortest time before
646   // it is invoked. In rest processes, GPIL()
647   // returns the time before a process occurs
648 
649   G4double lifeTime, shortestLifeTime;
650 
651   fAtRestDoItProcTriggered = 0;
652   shortestLifeTime = DBL_MAX;
653 
654   for (std::size_t ri = 0; ri < MAXofAtRestLoops; ++ri) {
655     fCurrentProcess = (*fAtRestGetPhysIntVector)[(G4int)ri];
656     if (fCurrentProcess == nullptr) {
657       (*fSelectedAtRestDoItVector)[ri] = InActivated;
658       continue;
659     }  // nullptr means the process is inactivated by a user on fly
660 
661     lifeTime = fCurrentProcess->AtRestGPIL(*fTrack, &fCondition);
662 
663     if (fCondition == Forced) {
664       (*fSelectedAtRestDoItVector)[ri] = Forced;
665     }
666     else {
667       (*fSelectedAtRestDoItVector)[ri] = InActivated;
668       if (lifeTime < shortestLifeTime) {
669         shortestLifeTime = lifeTime;
670         fAtRestDoItProcTriggered = G4int(ri);
671         fStep->GetPostStepPoint()->SetProcessDefinedStep(fCurrentProcess);
672       }
673     }
674   }
675 
676   (*fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
677 
678   fStep->SetStepLength(0.);  // the particle has stopped
679   fTrack->SetStepLength(0.);
680 
681   // Condition to avoid that stable ions are handled by Radioactive Decay.
682   // We use a very large time threshold (many orders of magnitude bigger than
683   // the universe's age) but not DBL_MAX because shortestLifeTime can be
684   // sometimes slightly smaller for stable ions.
685   if (shortestLifeTime < 1.0e+100)  // Unstable ion at rest: Radioactive Decay will decay it
686   {
687     // invoke selected process
688     //
689     for (std::size_t np = 0; np < MAXofAtRestLoops; ++np) {
690       //
691       // Note: DoItVector has inverse order against GetPhysIntVector
692       //       and SelectedAtRestDoItVector.
693       //
694       if ((*fSelectedAtRestDoItVector)[MAXofAtRestLoops - np - 1] != InActivated) {
695         fCurrentProcess = (*fAtRestDoItVector)[(G4int)np];
696         fParticleChange = fCurrentProcess->AtRestDoIt(*fTrack, *fStep);
697 
698         // Update Step
699         //
700         fParticleChange->UpdateStepForAtRest(fStep);
701 
702         // Now Store the secondaries from ParticleChange to SecondaryList
703         fN2ndariesAtRestDoIt += ProcessSecondariesFromParticleChange();
704 
705         // clear ParticleChange
706         fParticleChange->Clear();
707 
708       }  // if(fSelectedAtRestDoItVector[np] != InActivated){
709     }  // for(std::size_t np=0; np<MAXofAtRestLoops; ++np){
710   }
711   else  // Stable ion at rest
712   {
713     fStep->GetPostStepPoint()->SetProcessDefinedStep(fNoProcess);
714   }  // if(shortestLifeTime < 1.0e+100)
715 
716   fStep->UpdateTrack();
717 
718   fTrack->SetTrackStatus(fStopAndKill);
719 }
720 
721 /////////////////////////////////////////////////////////
722 void G4SteppingManager::InvokeAlongStepDoItProcs()
723 /////////////////////////////////////////////////////////
724 {
725   // If the current Step is defined by a 'ExclusivelyForced'
726   // PostStepDoIt, then don't invoke any AlongStepDoIt
727   //
728   if (fStepStatus == fExclusivelyForcedProc) {
729     return;  // Take note 'return' is here !!!
730   }
731 
732   // Invoke all active continuous processes
733   //
734   for (std::size_t ci = 0; ci < MAXofAlongStepLoops; ++ci) {
735     fCurrentProcess = (*fAlongStepDoItVector)[(G4int)ci];
736     if (fCurrentProcess == nullptr) continue;
737     // NULL means the process is inactivated by a user on fly.
738 
739     fParticleChange = fCurrentProcess->AlongStepDoIt(*fTrack, *fStep);
740 
741     // Update the PostStepPoint of Step according to ParticleChange
742     fParticleChange->UpdateStepForAlongStep(fStep);
743 
744 #ifdef G4VERBOSE
745     if (verboseLevel > 0) fVerbose->AlongStepDoItOneByOne();
746 #endif
747 
748     // Now Store the secondaries from ParticleChange to SecondaryList
749     fN2ndariesAlongStepDoIt += ProcessSecondariesFromParticleChange();
750 
751     // Set the track status according to what the process defined
752     // if kinetic energy >0, otherwise set  fStopButAlive
753     //
754     fTrack->SetTrackStatus(fParticleChange->GetTrackStatus());
755 
756     // clear ParticleChange
757     fParticleChange->Clear();
758   }
759 
760   fStep->UpdateTrack();
761   G4TrackStatus fNewStatus = fTrack->GetTrackStatus();
762 
763   if (fNewStatus == fAlive && fTrack->GetKineticEnergy() <= DBL_MIN) {
764     if (MAXofAtRestLoops > 0)
765       fNewStatus = fStopButAlive;
766     else
767       fNewStatus = fStopAndKill;
768     fTrack->SetTrackStatus(fNewStatus);
769   }
770 }
771 
772 ////////////////////////////////////////////////////////
773 void G4SteppingManager::InvokePostStepDoItProcs()
774 ////////////////////////////////////////////////////////
775 {
776   // Invoke the specified discrete processes
777   //
778   for (std::size_t np = 0; np < MAXofPostStepLoops; ++np) {
779     //
780     // Note: DoItVector has inverse order against GetPhysIntVector
781     //       and SelectedPostStepDoItVector.
782     //
783     G4int Cond = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops - np - 1];
784     if (Cond != InActivated) {
785       if (((Cond == NotForced) && (fStepStatus == fPostStepDoItProc)) ||
786           ((Cond == Forced) && (fStepStatus != fExclusivelyForcedProc)) ||
787           ((Cond == ExclusivelyForced) && (fStepStatus == fExclusivelyForcedProc)) ||
788           ((Cond == StronglyForced)))
789       {
790         InvokePSDIP(np);
791         if ((np == 0) && (fTrack->GetNextVolume() == nullptr)) {
792           fStepStatus = fWorldBoundary;
793           fStep->GetPostStepPoint()->SetStepStatus(fStepStatus);
794         }
795       }
796     }
797 
798     // Exit from PostStepLoop if the track has been killed,
799     // but extra treatment for processes with Strongly Forced flag
800     //
801     if (fTrack->GetTrackStatus() == fStopAndKill) {
802       for (std::size_t np1 = np + 1; np1 < MAXofPostStepLoops; ++np1) {
803         G4int Cond2 = (*fSelectedPostStepDoItVector)[MAXofPostStepLoops - np1 - 1];
804         if (Cond2 == StronglyForced) {
805           InvokePSDIP(np1);
806         }
807       }
808       break;
809     }
810   }
811 }
812 
813 ////////////////////////////////////////////////////////
814 void G4SteppingManager::InvokePSDIP(size_t np)
815 ////////////////////////////////////////////////////////
816 {
817   fCurrentProcess = (*fPostStepDoItVector)[(G4int)np];
818   fParticleChange = fCurrentProcess->PostStepDoIt(*fTrack, *fStep);
819 
820   // Update PostStepPoint of Step according to ParticleChange
821   fParticleChange->UpdateStepForPostStep(fStep);
822 
823 #ifdef G4VERBOSE
824   if (verboseLevel > 0) fVerbose->PostStepDoItOneByOne();
825 #endif
826 
827   // Update G4Track according to ParticleChange after each PostStepDoIt
828   fStep->UpdateTrack();
829 
830   // Update safety after each invocation of PostStepDoIts
831   fStep->GetPostStepPoint()->SetSafety(CalculateSafety());
832 
833   // Now Store the secondaries from ParticleChange to SecondaryList
834   fN2ndariesPostStepDoIt += ProcessSecondariesFromParticleChange();
835 
836   // Set the track status according to what the process defined
837   fTrack->SetTrackStatus(fParticleChange->GetTrackStatus());
838 
839   // clear ParticleChange
840   fParticleChange->Clear();
841 }
842