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