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

#include <G4ITStepProcessor.hh>

Public Member Functions

 G4ITStepProcessor ()
 
virtual ~G4ITStepProcessor ()
 
void SetPreviousStepTime (G4double)
 
G4TrackGetTrack ()
 
G4StepGetStep ()
 
const G4StepGetStep () const
 
void SetStep (G4Step *val)
 
G4TrackVectorGetSecondaries () const
 
void SetTrackingManager (G4ITTrackingManager *trackMan)
 
G4ITTrackingManagerGetTrackingManager ()
 
virtual void Initialize ()
 
void ForceReInitialization ()
 
void ResetLeadingTracks ()
 
void PrepareLeadingTracks ()
 
G4double ComputeInteractionLength (double previousTimeStep)
 
void DefinePhysicalStepLength (G4Track *)
 
G4double GetILTimeStep ()
 
void DoIt (double timeStep)
 
void ExtractDoItData ()
 
void Stepping (G4Track *, const double &)
 
void FindTransportationStep ()
 
double GetInteractionTime ()
 
const G4TrackGetTrack () const
 
void CleanProcessor ()
 
std::size_t GetAtRestDoItProcTriggered () const
 
G4GPILSelection GetGPILSelection () const
 
G4int GetN2ndariesAlongStepDoIt () const
 
G4int GetN2ndariesAtRestDoIt () const
 
G4int GetN2ndariesPostStepDoIt () const
 
const G4VITProcessGetCurrentProcess () const
 
G4double GetPhysIntLength () const
 
std::size_t GetPostStepAtTimeDoItProcTriggered () const
 
std::size_t GetPostStepDoItProcTriggered () const
 
const ProcessGeneralInfoGetCurrentProcessInfo () const
 
const G4ITStepProcessorStateGetProcessorState () const
 
const G4VParticleChangeGetParticleChange () const
 
const G4VPhysicalVolumeGetCurrentVolume () const
 
G4ForceCondition GetCondition () const
 

Protected Member Functions

void ExtractILData ()
 
void SetupGeneralProcessInfo (G4ParticleDefinition *, G4ProcessManager *)
 
void ClearProcessInfo ()
 
void SetTrack (G4Track *)
 
void GetProcessInfo ()
 
void SetupMembers ()
 
void ResetSecondaries ()
 
void InitDefineStep ()
 
void SetInitialStep ()
 
void GetAtRestIL ()
 
void DoDefinePhysicalStepLength ()
 
void DoStepping ()
 
void PushSecondaries ()
 
void ActiveOnlyITProcess ()
 
void ActiveOnlyITProcess (G4ProcessManager *)
 
void DealWithSecondaries (G4int &)
 
void InvokeAtRestDoItProcs ()
 
void InvokeAlongStepDoItProcs ()
 
void InvokePostStepDoItProcs ()
 
void InvokePSDIP (std::size_t)
 
void InvokeTransportationProc ()
 
void SetNavigator (G4ITNavigator *value)
 
G4double CalculateSafety ()
 
void ApplyProductionCut (G4Track *)
 
 G4ITStepProcessor (const G4ITStepProcessor &other)
 
G4ITStepProcessoroperator= (const G4ITStepProcessor &other)
 

Friends

class G4Scheduler
 

Detailed Description

Its role is the same as G4StepManager :

  • Find the minimum physical length and corresponding time step
  • Step one track BUT on a given time step.

Definition at line 152 of file G4ITStepProcessor.hh.

Constructor & Destructor Documentation

◆ G4ITStepProcessor() [1/2]

G4ITStepProcessor::G4ITStepProcessor ( )

Definition at line 71 of file G4ITStepProcessor.cc.

72{
73 fpVerbose = nullptr;
74 // fpUserSteppingAction = 0 ;
75 fStoreTrajectory = 0;
76 fpTrackingManager = nullptr;
77 fpNavigator = nullptr;
78 kCarTolerance = -1.;
79 fInitialized = false;
80 fPreviousTimeStep = DBL_MAX;
81 fILTimeStep = DBL_MAX;
82 fpTrackContainer = nullptr;
83
86}
#define DBL_MAX
Definition templates.hh:62

◆ ~G4ITStepProcessor()

G4ITStepProcessor::~G4ITStepProcessor ( )
virtual

Definition at line 223 of file G4ITStepProcessor.cc.

224{
225 if(fpStep != nullptr)
226 {
227 fpStep->DeleteSecondaryVector();
228 delete fpStep;
229 }
230
231 delete fpSecondary;
233 //G4ITTransportationManager::DeleteInstance();
234
235 // if(fpUserSteppingAction) delete fpUserSteppingAction;
236}
void DeleteSecondaryVector()

◆ G4ITStepProcessor() [2/2]

G4ITStepProcessor::G4ITStepProcessor ( const G4ITStepProcessor & other)
protected

Definition at line 239 of file G4ITStepProcessor.cc.

240{
241 fpVerbose = rhs.fpVerbose;
242 fStoreTrajectory = rhs.fStoreTrajectory;
243
244 // fpUserSteppingAction = 0 ;
245 fpTrackingManager = nullptr;
246 fpNavigator = nullptr;
247 fInitialized = false;
248
249 kCarTolerance = rhs.kCarTolerance;
250 fInitialized = false;
251 fPreviousTimeStep = DBL_MAX;
252
255 fpTrackContainer = nullptr;
256 fILTimeStep = DBL_MAX;
257}

Member Function Documentation

◆ ActiveOnlyITProcess() [1/2]

void G4ITStepProcessor::ActiveOnlyITProcess ( )
protected

Definition at line 284 of file G4ITStepProcessor.cc.

285{
286 // Method not used for the time being
287#ifdef debug
288 G4cout<<"G4ITStepProcessor::CloneProcesses: is called"<<G4endl;
289#endif
290
293 ->GetIterator();
294
296 // TODO : Ne faire la boucle que sur les IT **** !!!
297 while((*theParticleIterator)())
298 {
299 G4ParticleDefinition* particle = theParticleIterator->value();
300 G4ProcessManager* pm = particle->GetProcessManager();
301
302 if(pm == nullptr)
303 {
304 G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl<< " ProcessManager is NULL for particle = "
305 << particle->GetParticleName() << ", PDG_code = "
306 << particle->GetPDGEncoding() << G4endl;
307 G4Exception("G4ITStepProcessor::GetProcessNumber()", "ITStepProcessor0001",
308 FatalException, "Process Manager is not found.");
309 return;
310 }
311
313 }
314}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#define theParticleIterator
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
void reset(G4bool ifSkipIon=true)
G4PTblDicIterator * GetIterator() const
static G4ParticleTable * GetParticleTable()

Referenced by ActiveOnlyITProcess().

◆ ActiveOnlyITProcess() [2/2]

void G4ITStepProcessor::ActiveOnlyITProcess ( G4ProcessManager * processManager)
protected

Definition at line 318 of file G4ITStepProcessor.cc.

319{
320 // Method not used for the time being
321 G4ProcessVector* processVector = processManager->GetProcessList();
322
323 G4VITProcess* itProcess = nullptr;
324 for(G4int i = 0; i < (G4int)processVector->size(); ++i)
325 {
326 G4VProcess* base_process = (*processVector)[i];
327 itProcess = dynamic_cast<G4VITProcess*>(base_process);
328
329 if(itProcess == nullptr)
330 {
331 processManager->SetProcessActivation(base_process, false);
332 }
333 }
334}
int G4int
Definition G4Types.hh:85
G4VProcess * SetProcessActivation(G4VProcess *aProcess, G4bool fActive)
G4ProcessVector * GetProcessList() const
std::size_t size() const

◆ ApplyProductionCut()

void G4ITStepProcessor::ApplyProductionCut ( G4Track * aSecondary)
protected

Definition at line 912 of file G4ITStepProcessor2.cc.

913{
914 G4bool tBelowCutEnergyAndSafety = false;
915 G4int tPtclIdx = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
916 if(tPtclIdx < 0)
917 {
918 return;
919 }
920 G4ProductionCutsTable* tCutsTbl =
922 G4int tCoupleIdx = tCutsTbl->GetCoupleIndex(fpPreStepPoint
923 ->GetMaterialCutsCouple());
924 G4double tProdThreshold =
925 (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
926 if(aSecondary->GetKineticEnergy() < tProdThreshold)
927 {
928 tBelowCutEnergyAndSafety = true;
929 if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
930 {
931 G4double currentRange
933 aSecondary->GetKineticEnergy(),
934 fpPreStepPoint->GetMaterialCutsCouple());
935 tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
936 }
937 }
938
939 if(tBelowCutEnergyAndSafety)
940 {
941 if(!(aSecondary->IsGoodForTracking()))
942 {
943 // Add kinetic energy to the total energy deposit
944 fpStep->AddTotalEnergyDeposit(aSecondary->GetKineticEnergy());
945 aSecondary->SetKineticEnergy(0.0);
946 }
947 }
948}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4double GetCharge() const
static G4LossTableManager * Instance()
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4int GetCoupleIndex(const G4MaterialCutsCouple *aCouple) const
static G4int GetIndex(const G4String &name)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void AddTotalEnergyDeposit(G4double value)
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
void SetKineticEnergy(const G4double aValue)
G4bool IsGoodForTracking() const
#define DBL_MIN
Definition templates.hh:54

Referenced by DealWithSecondaries().

◆ CalculateSafety()

G4double G4ITStepProcessor::CalculateSafety ( )
inlineprotected

Definition at line 440 of file G4ITStepProcessor.hh.

441{
442 return std::max(fpState->fEndpointSafety - (fpState->fEndpointSafOrigin
443 - fpPostStepPoint->GetPosition()).mag(),
444 kCarTolerance);
445}
const G4ThreeVector & GetPosition() const

Referenced by ApplyProductionCut(), and InvokePSDIP().

◆ CleanProcessor()

void G4ITStepProcessor::CleanProcessor ( )
inline

Definition at line 456 of file G4ITStepProcessor.hh.

457{
458 fTimeStep = DBL_MAX;
459 fPhysIntLength = DBL_MAX;
460
461 fpState = nullptr;
462 fpTrack = nullptr;
463 fpTrackingInfo = nullptr;
464 fpITrack = nullptr;
465 fpStep = nullptr;
466 fpPreStepPoint = nullptr;
467 fpPostStepPoint = nullptr;
468
469 fpParticleChange = nullptr;
470
471 fpCurrentVolume = nullptr;
472 // fpSensitive = 0;
473
474 fpSecondary = nullptr;
475
476 fpTransportation = nullptr;
477
478 fpCurrentProcess= nullptr;
479 fpProcessInfo = nullptr;
480
481 fAtRestDoItProcTriggered = INT_MAX;
482 fPostStepDoItProcTriggered = INT_MAX;
483 fPostStepAtTimeDoItProcTriggered = INT_MAX;
484 fGPILSelection = NotCandidateForSelection;
485 fCondition = NotForced;
486}
@ NotForced
@ NotCandidateForSelection
#define INT_MAX
Definition templates.hh:90

Referenced by ExtractDoItData(), ExtractILData(), G4ITStepProcessor(), G4ITStepProcessor(), Initialize(), and Stepping().

◆ ClearProcessInfo()

void G4ITStepProcessor::ClearProcessInfo ( )
protected

Definition at line 165 of file G4ITStepProcessor.cc.

166{
167 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it;
168
169 for(it = fProcessGeneralInfoMap.begin(); it != fProcessGeneralInfoMap.end();
170 it++)
171 {
172 if(it->second != nullptr)
173 {
174 delete it->second;
175 it->second = 0;
176 }
177 }
178
179 fProcessGeneralInfoMap.clear();
180}

Referenced by ForceReInitialization(), and ~G4ITStepProcessor().

◆ ComputeInteractionLength()

G4double G4ITStepProcessor::ComputeInteractionLength ( double previousTimeStep)

Definition at line 598 of file G4ITStepProcessor.cc.

599{
600 G4TrackManyList* mainList = fpTrackContainer->GetMainList();
601 G4TrackManyList::iterator it = mainList ->begin();
602 G4TrackManyList::iterator end = mainList ->end();
603
604 SetPreviousStepTime(previousTimeStep);
605
606 fILTimeStep = DBL_MAX;
607
608 for (; it != end; )
609 {
610 G4Track * track = *it;
611
612#ifdef DEBUG
613 G4cout << "*CIL* " << GetIT(track)->GetName()
614 << " ID: " << track->GetTrackID()
615 << " at time : " << track->GetGlobalTime()
616 << G4endl;
617#endif
618
619 ++it;
621
623 }
624
625 return fILTimeStep;
626}
G4IT * GetIT(const G4Track *track)
Definition G4IT.cc:48
void SetPreviousStepTime(G4double)
void DefinePhysicalStepLength(G4Track *)
G4TrackList * GetMainList(Key)
virtual const G4String & GetName() const =0
G4int GetTrackID() const
G4double GetGlobalTime() const

Referenced by G4Scheduler::Stepping().

◆ DealWithSecondaries()

void G4ITStepProcessor::DealWithSecondaries ( G4int & counter)
protected

Definition at line 63 of file G4ITStepProcessor2.cc.

64{
65 // Now Store the secondaries from ParticleChange to SecondaryList
66 G4Track* tempSecondaryTrack;
67
68 for(G4int DSecLoop = 0; DSecLoop < fpParticleChange->GetNumberOfSecondaries();
69 DSecLoop++)
70 {
71 tempSecondaryTrack = fpParticleChange->GetSecondary(DSecLoop);
72
73 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
74 {
75 ApplyProductionCut(tempSecondaryTrack);
76 }
77
78 // Set parentID
79 tempSecondaryTrack->SetParentID(fpTrack->GetTrackID());
80
81 // Set the process pointer which created this track
82 tempSecondaryTrack->SetCreatorProcess(fpCurrentProcess);
83
84 // If this 2ndry particle has 'zero' kinetic energy, make sure
85 // it invokes a rest process at the beginning of the tracking
86 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN)
87 {
88 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
89 if (pm->GetAtRestProcessVector()->entries()>0)
90 {
91 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
92 fpSecondary->push_back( tempSecondaryTrack );
93 fN2ndariesAtRestDoIt++;
94 }
95 else
96 {
97 delete tempSecondaryTrack;
98 }
99 }
100 else
101 {
102 fpSecondary->push_back( tempSecondaryTrack );
103 counter++;
104 }
105 } //end of loop on secondary
106}
@ fStopButAlive
void ApplyProductionCut(G4Track *)
G4bool GetApplyCutsFlag() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t entries() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
void SetParentID(const G4int aValue)
void SetCreatorProcess(const G4VProcess *aValue)
G4int GetNumberOfSecondaries() const
G4Track * GetSecondary(G4int anIndex) const

Referenced by InvokeAlongStepDoItProcs(), InvokeAtRestDoItProcs(), and InvokePSDIP().

◆ DefinePhysicalStepLength()

void G4ITStepProcessor::DefinePhysicalStepLength ( G4Track * track)

Definition at line 693 of file G4ITStepProcessor.cc.

694{
695 SetTrack(track);
697}
void SetTrack(G4Track *)

Referenced by ComputeInteractionLength().

◆ DoDefinePhysicalStepLength()

void G4ITStepProcessor::DoDefinePhysicalStepLength ( )
protected

Definition at line 951 of file G4ITStepProcessor.cc.

952{
953
955
956#ifdef G4VERBOSE
957 // !!!!! Verbose
958 if(fpVerbose != nullptr) fpVerbose->DPSLStarted();
959#endif
960
961 G4TrackStatus trackStatus = fpTrack->GetTrackStatus();
962
963 if(trackStatus == fStopAndKill)
964 {
965 return;
966 }
967
968 if(trackStatus == fStopButAlive)
969 {
970 fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator
971 ->GetNavigatorState());
972 fpNavigator->ResetNavigatorState();
973 return GetAtRestIL();
974 }
975
976 // Find minimum Step length and corresponding time
977 // demanded by active disc./cont. processes
978
979 // ReSet the counter etc.
980 fpState->fPhysicalStep = DBL_MAX; // Initialize by a huge number
981 fPhysIntLength = DBL_MAX; // Initialize by a huge number
982
983 G4double proposedTimeStep = DBL_MAX;
984 G4VProcess* processWithPostStepGivenByTimeStep(nullptr);
985
986 // GPIL for PostStep
987 fPostStepDoItProcTriggered = fpProcessInfo->MAXofPostStepLoops;
988 fPostStepAtTimeDoItProcTriggered = fpProcessInfo->MAXofPostStepLoops;
989
990 // G4cout << "fpProcessInfo->MAXofPostStepLoops : "
991 // << fpProcessInfo->MAXofPostStepLoops
992 // << " mol : " << fpITrack -> GetName()
993 // << " id : " << fpTrack->GetTrackID()
994 // << G4endl;
995
996 for(std::size_t np = 0; np < fpProcessInfo->MAXofPostStepLoops; ++np)
997 {
998 fpCurrentProcess = dynamic_cast<G4VITProcess*>((*fpProcessInfo
1000 if(fpCurrentProcess == nullptr)
1001 {
1003 continue;
1004 } // NULL means the process is inactivated by a user on fly.
1005
1006 fCondition = NotForced;
1007 fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
1008 ->GetProcessID()));
1009
1010 // G4cout << "Is going to call : "
1011 // << fpCurrentProcess -> GetProcessName()
1012 // << G4endl;
1013 fPhysIntLength = fpCurrentProcess->PostStepGPIL(*fpTrack,
1014 fpState->fPreviousStepSize,
1015 &fCondition);
1016
1017#ifdef G4VERBOSE
1018 // !!!!! Verbose
1019 if(fpVerbose != nullptr) fpVerbose->DPSLPostStep();
1020#endif
1021
1022 fpCurrentProcess->ResetProcessState();
1023 //fpCurrentProcess->SetProcessState(0);
1024
1025 switch(fCondition)
1026 {
1027 case ExclusivelyForced: // Will need special treatment
1030 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpCurrentProcess);
1031 break;
1032
1033 case Conditionally:
1034 // (fpState->fSelectedPostStepDoItVector)[np] = Conditionally;
1035 G4Exception("G4ITStepProcessor::DefinePhysicalStepLength()",
1036 "ITStepProcessor0008",
1038 "This feature is no more supported");
1039 break;
1040
1041 case Forced:
1042 (fpState->fSelectedPostStepDoItVector)[np] = Forced;
1043 break;
1044
1045 case StronglyForced:
1047 break;
1048
1049 default:
1051 break;
1052 }
1053
1054 if(fCondition == ExclusivelyForced)
1055 {
1056 for(std::size_t nrest = np + 1; nrest < fpProcessInfo->MAXofPostStepLoops;
1057 ++nrest)
1058 {
1059 (fpState->fSelectedPostStepDoItVector)[nrest] = InActivated;
1060 }
1061 return; // Please note the 'return' at here !!!
1062 }
1063
1064 if(fPhysIntLength < fpState->fPhysicalStep)
1065 {
1066 // To avoid checking whether the process is actually
1067 // proposing a time step, the returned time steps are
1068 // negative (just for tagging)
1069 if(fpCurrentProcess->ProposesTimeStep())
1070 {
1071 fPhysIntLength *= -1;
1072 if(fPhysIntLength < proposedTimeStep)
1073 {
1074 proposedTimeStep = fPhysIntLength;
1075 fPostStepAtTimeDoItProcTriggered = np;
1076 processWithPostStepGivenByTimeStep = fpCurrentProcess;
1077 }
1078 }
1079 else
1080 {
1081 fpState->fPhysicalStep = fPhysIntLength;
1082 fpState->fStepStatus = fPostStepDoItProc;
1083 fPostStepDoItProcTriggered = G4int(np);
1084 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpCurrentProcess);
1085 }
1086 }
1087 }
1088
1089 // GPIL for AlongStep
1090 fpState->fProposedSafety = DBL_MAX;
1091 G4double safetyProposedToAndByProcess = fpState->fProposedSafety;
1092
1093 for(std::size_t kp = 0; kp < fpProcessInfo->MAXofAlongStepLoops; ++kp)
1094 {
1095 fpCurrentProcess = dynamic_cast<G4VITProcess*>((*fpProcessInfo
1097 if(fpCurrentProcess == nullptr) continue;
1098 // NULL means the process is inactivated by a user on fly.
1099
1100 fpCurrentProcess->SetProcessState(
1101 fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
1102 fPhysIntLength =
1103 fpCurrentProcess->AlongStepGPIL(*fpTrack,
1104 fpState->fPreviousStepSize,
1105 fpState->fPhysicalStep,
1106 safetyProposedToAndByProcess,
1107 &fGPILSelection);
1108
1109#ifdef G4VERBOSE
1110 // !!!!! Verbose
1111 if(fpVerbose != nullptr) fpVerbose->DPSLAlongStep();
1112#endif
1113
1114 if(fPhysIntLength < fpState->fPhysicalStep)
1115 {
1116 fpState->fPhysicalStep = fPhysIntLength;
1117 // Should save PS and TS in IT
1118
1119 // Check if the process wants to be the GPIL winner. For example,
1120 // multi-scattering proposes Step limit, but won't be the winner.
1121 if(fGPILSelection == CandidateForSelection)
1122 {
1124 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpCurrentProcess);
1125 }
1126
1127 // Transportation is assumed to be the last process in the vector
1128 if(kp == fpProcessInfo->MAXofAlongStepLoops - 1)
1129 {
1130 fpTransportation = dynamic_cast<G4ITTransportation*>(fpCurrentProcess);
1131
1132 if(fpTransportation == nullptr)
1133 {
1134 G4ExceptionDescription exceptionDescription;
1135 exceptionDescription << "No transportation process found ";
1136 G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength",
1137 "ITStepProcessor0009",
1139 exceptionDescription);
1140 }
1141
1142 fTimeStep = fpTransportation->GetInteractionTimeLeft();
1143
1144 if(fpTrack->GetNextVolume() != nullptr) fpState->fStepStatus = fGeomBoundary;
1145 else fpState->fStepStatus = fWorldBoundary;
1146 }
1147 }
1148 else
1149 {
1150 if(kp == fpProcessInfo->MAXofAlongStepLoops - 1)
1151 {
1152 fpTransportation = dynamic_cast<G4ITTransportation*>(fpCurrentProcess);
1153
1154 if(fpTransportation == nullptr)
1155 {
1156 G4ExceptionDescription exceptionDescription;
1157 exceptionDescription << "No transportation process found ";
1158 G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength",
1159 "ITStepProcessor0010",
1161 exceptionDescription);
1162 }
1163
1164 fTimeStep = fpTransportation->GetInteractionTimeLeft();
1165 }
1166 }
1167
1168 // Handle PostStep processes sending back time steps rather than space length
1169 if(proposedTimeStep < fTimeStep)
1170 {
1171 if(fPostStepAtTimeDoItProcTriggered < fpProcessInfo->MAXofPostStepLoops)
1172 {
1173 if((fpState->fSelectedPostStepDoItVector)[fPostStepAtTimeDoItProcTriggered] == InActivated)
1174 {
1175 (fpState->fSelectedPostStepDoItVector)[fPostStepAtTimeDoItProcTriggered] =
1176 NotForced;
1177 // (fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = InActivated;
1178
1179 fpState->fStepStatus = fPostStepDoItProc;
1180 fpStep->GetPostStepPoint()->SetProcessDefinedStep(processWithPostStepGivenByTimeStep);
1181
1182 fTimeStep = proposedTimeStep;
1183
1184 fpTransportation->ComputeStep(*fpTrack,
1185 *fpStep,
1186 fTimeStep,
1187 fpState->fPhysicalStep);
1188 }
1189 }
1190 }
1191 else
1192 {
1193 if(fPostStepDoItProcTriggered < fpProcessInfo->MAXofPostStepLoops)
1194 {
1195 if((fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] == InActivated)
1196 {
1197 (fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] =
1198 NotForced;
1199 }
1200 }
1201 }
1202
1203// fpCurrentProcess->SetProcessState(0);
1204 fpCurrentProcess->ResetProcessState();
1205
1206 // Make sure to check the safety, even if Step is not limited
1207 // by this process. J. Apostolakis, June 20, 1998
1208 //
1209 if(safetyProposedToAndByProcess < fpState->fProposedSafety)
1210 // proposedSafety keeps the smallest value:
1211 fpState->fProposedSafety = safetyProposedToAndByProcess;
1212 else
1213 // safetyProposedToAndByProcess always proposes a valid safety:
1214 safetyProposedToAndByProcess = fpState->fProposedSafety;
1215
1216 }
1217
1218 fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator->GetNavigatorState());
1219 fpNavigator->ResetNavigatorState();
1220}
@ FatalErrorInArgument
std::ostringstream G4ExceptionDescription
@ InActivated
@ StronglyForced
@ Conditionally
@ ExclusivelyForced
@ Forced
@ CandidateForSelection
@ fGeomBoundary
@ fWorldBoundary
@ fPostStepDoItProc
@ fAlongStepDoItProc
@ fExclusivelyForcedProc
G4TrackStatus
@ fStopAndKill
G4SelectedPostStepDoItVector fSelectedPostStepDoItVector
virtual void ComputeStep(const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
G4TrackingInformation * GetTrackingInfo()
Definition G4IT.hh:143
void SetProcessDefinedStep(const G4VProcess *aValue)
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetNextVolume() const
G4shared_ptr< G4ProcessState_Lock > GetProcessState(size_t index)
void SetNavigatorState(G4ITNavigatorState_Lock *)
G4bool ProposesTimeStep() const
size_t GetProcessID() const
void ResetProcessState()
void SetProcessState(G4shared_ptr< G4ProcessState_Lock > aProcInfo)
G4double GetInteractionTimeLeft()
virtual void DPSLAlongStep()=0
virtual void DPSLStarted()=0
virtual void DPSLPostStep()=0
G4double PostStepGPIL(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double AlongStepGPIL(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
G4ProcessVector * fpPostStepGetPhysIntVector
std::size_t MAXofAlongStepLoops
G4ProcessVector * fpAlongStepGetPhysIntVector

Referenced by DefinePhysicalStepLength().

◆ DoIt()

void G4ITStepProcessor::DoIt ( double timeStep)

Definition at line 110 of file G4ITStepProcessor2.cc.

119{
120 if(fpVerbose != nullptr) fpVerbose->DoItStarted();
121
122 G4TrackManyList* mainList = fpTrackContainer->GetMainList();
123 G4TrackManyList::iterator it = mainList->end();
124 it--;
125 std::size_t initialSize = mainList->size();
126
127// G4cout << "initialSize = " << initialSize << G4endl;
128
129 for(std::size_t i = 0 ; i < initialSize ; ++i)
130 {
131
132// G4cout << "i = " << i << G4endl;
133
134 G4Track* track = *it;
135 if (track == nullptr)
136 {
137 G4ExceptionDescription exceptionDescription;
138 exceptionDescription << "No track was pop back the main track list.";
139 G4Exception("G4ITStepProcessor::DoIt", "NO_TRACK",
140 FatalException, exceptionDescription);
141 }
142 // G4TrackManyList::iterator next_it (it);
143 // next_it--;
144 // it = next_it;
145
146 it--;
147 // Must be called before EndTracking(track)
148 // Otherwise the iterator will point to the list of killed tracks
149
150 if(track->GetTrackStatus() == fStopAndKill)
151 {
152 fpTrackingManager->EndTracking(track);
153// G4cout << GetIT(track)->GetName() << G4endl;
154// G4cout << " ************************ CONTINUE ********************" << G4endl;
155 continue;
156 }
157
158#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
159 MemStat mem_first, mem_second, mem_diff;
160 mem_first = MemoryUsage();
161#endif
162
163 Stepping(track, timeStep);
164
165#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
166 MemStat mem_intermediaire = MemoryUsage();
167 mem_diff = mem_intermediaire-mem_first;
168 G4cout << "\t\t >> || MEM || In DoIT with track "
169 << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
170#endif
171
173
174#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
175 mem_second = MemoryUsage();
176 mem_diff = mem_second-mem_first;
177 G4cout << "\t >> || MEM || In DoIT with track "
178 << track->GetTrackID()
179 << ", diff is : " << mem_diff << G4endl;
180#endif
181 }
182
183
184 fpTrackContainer->MergeSecondariesWithMainList();
185 fpTrackContainer->KillTracks(); // (18-06-15 : MK) Remove it ?
186 fLeadingTracks.Reset();
187}
void Stepping(G4Track *, const double &)
void MergeSecondariesWithMainList()
void EndTracking(G4Track *)
size_t size() const
virtual void DoItStarted()=0
MemStat MemoryUsage()
Definition G4MemStat.cc:55

Referenced by G4Scheduler::Stepping().

◆ DoStepping()

void G4ITStepProcessor::DoStepping ( )
protected

Definition at line 292 of file G4ITStepProcessor2.cc.

293{
294 SetupMembers();
295
296#ifdef DEBUG_MEM
297 MemStat mem_first, mem_second, mem_diff;
298#endif
299
300#ifdef DEBUG_MEM
301 mem_first = MemoryUsage();
302#endif
303
304#ifdef G4VERBOSE
305 if(fpVerbose != nullptr) fpVerbose->PreStepVerbose(fpTrack);
306#endif
307
308 if(fpProcessInfo == nullptr)
309 {
310 G4ExceptionDescription exceptionDescription;
311 exceptionDescription << "No process info found for particle :"
312 << fpTrack->GetDefinition()->GetParticleName();
313 G4Exception("G4ITStepProcessor::DoStepping",
314 "ITStepProcessor0012",
316 exceptionDescription);
317 return;
318 }
319// else if(fpTrack->GetTrackStatus() == fStopAndKill)
320// {
321// fpState->fStepStatus = fUndefined;
322// return;
323// }
324
325 if(fpProcessInfo->MAXofPostStepLoops == 0 &&
326 fpProcessInfo->MAXofAlongStepLoops == 0
327 && fpProcessInfo->MAXofAtRestLoops == 0)
328 {/*
329 G4ExceptionDescription exceptionDescription;
330 exceptionDescription << "No process was found for particle :"
331 << fpTrack->GetDefinition()->GetParticleName();
332 G4Exception("G4ITStepProcessor::DoStepping",
333 "ITStepProcessorNoProcess",
334 JustWarning,
335 exceptionDescription);
336
337 fpTrack->SetTrackStatus(fStopAndKill);
338 fpState->fStepStatus = fUndefined;*/
339 return;
340 }
341
342 //--------
343 // Prelude
344 //--------
345#ifdef G4VERBOSE
346 // !!!!! Verbose
347 if(fpVerbose != nullptr) fpVerbose->NewStep();
348#endif
349
350 //---------------------------------
351 // AtRestStep, AlongStep and PostStep Processes
352 //---------------------------------
353
354 fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
355// fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
356// fpTrack->GetMomentumDirection(),
357// *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
358// fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
359 // We reset the navigator state before checking for AtRest
360 // in case a AtRest processe would use a navigator info
361
362#ifdef DEBUG_MEM
363 MemStat mem_intermediaire = MemoryUsage();
364 mem_diff = mem_intermediaire-mem_first;
365 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After dealing with navigator with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
366#endif
367
368 if(fpTrack->GetTrackStatus() == fStopButAlive)
369 {
370 if(fpProcessInfo->MAXofAtRestLoops > 0 && fpProcessInfo->fpAtRestDoItVector
371 != nullptr) // second condition to make coverity happy
372 {
373 //-----------------
374 // AtRestStepDoIt
375 //-----------------
377 fpState->fStepStatus = fAtRestDoItProc;
378 fpStep->GetPostStepPoint()->SetStepStatus(fpState->fStepStatus);
379
380#ifdef G4VERBOSE
381 // !!!!! Verbose
382 if(fpVerbose != nullptr) fpVerbose->AtRestDoItInvoked();
383#endif
384
385 }
386 // Make sure the track is killed
387 // fpTrack->SetTrackStatus(fStopAndKill);
388 }
389 else // if(fTimeStep > 0.) // Bye, because PostStepIL can return 0 => time =0
390 {
391 if(fpITrack == nullptr)
392 {
393 G4ExceptionDescription exceptionDescription;
394 exceptionDescription << " !!! TrackID : " << fpTrack->GetTrackID()
395 << G4endl<< " !!! Track status : "<< fpTrack->GetTrackStatus() << G4endl
396 << " !!! Particle Name : "<< fpTrack -> GetDefinition() -> GetParticleName() << G4endl
397 << "No G4ITStepProcessor::fpITrack found" << G4endl;
398
399 G4Exception("G4ITStepProcessor::DoStepping",
400 "ITStepProcessor0013",
402 exceptionDescription);
403 return; // to make coverity happy
404 }
405
406 if(!fpITrack->GetTrackingInfo()->IsLeadingStep())
407 {
408 // In case the track has NOT the minimum step length
409 // Given the final step time, the transportation
410 // will compute the final position of the particle
412 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpTransportation);
414 }
415
416#ifdef DEBUG_MEM
417 mem_intermediaire = MemoryUsage();
418 mem_diff = mem_intermediaire-mem_first;
419 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After FindTransportationStep() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
420#endif
421
422 // Store the Step length (geometrical length) to G4Step and G4Track
423 fpTrack->SetStepLength(fpState->fPhysicalStep);
424 fpStep->SetStepLength(fpState->fPhysicalStep);
425
426 G4double GeomStepLength = fpState->fPhysicalStep;
427
428 // Store StepStatus to PostStepPoint
429 fpStep->GetPostStepPoint()->SetStepStatus(fpState->fStepStatus);
430
431 // Invoke AlongStepDoIt
433
434#ifdef DEBUG_MEM
435 mem_intermediaire = MemoryUsage();
436 mem_diff = mem_intermediaire-mem_first;
437 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeAlongStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
438#endif
439
440#ifdef G4VERBOSE
441 // !!!!! Verbose
442 if(fpVerbose != nullptr) fpVerbose->AlongStepDoItAllDone();
443#endif
444
445 // Update track by taking into account all changes by AlongStepDoIt
446 // fpStep->UpdateTrack(); // done in InvokeAlongStepDoItProcs
447
448 // Update safety after invocation of all AlongStepDoIts
449 fpState->fEndpointSafOrigin = fpPostStepPoint->GetPosition();
450
451 fpState->fEndpointSafety =
452 std::max(fpState->fProposedSafety - GeomStepLength, kCarTolerance);
453
454 fpStep->GetPostStepPoint()->SetSafety(fpState->fEndpointSafety);
455
456 if(GetIT(fpTrack)->GetTrackingInfo()->IsLeadingStep())
457 {
458 // Invoke PostStepDoIt including G4ITTransportation::PSDI
460
461#ifdef DEBUG_MEM
462 mem_intermediaire = MemoryUsage();
463 mem_diff = mem_intermediaire-mem_first;
464 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokePostStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
465#endif
466#ifdef G4VERBOSE
467 // !!!!! Verbose
468 if(fpVerbose != nullptr) fpVerbose->StepInfoForLeadingTrack();
469#endif
470 }
471 else
472 {
473 // Only invoke transportation and all other forced processes
475 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpTransportation);
476
477#ifdef DEBUG_MEM
478 mem_intermediaire = MemoryUsage();
479 mem_diff = mem_intermediaire-mem_first;
480 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeTransportationProc() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
481#endif
482 }
483
484#ifdef G4VERBOSE
485 // !!!!! Verbose
486 if(fpVerbose != nullptr) fpVerbose->PostStepDoItAllDone();
487#endif
488 }
489
490 fpNavigator->ResetNavigatorState();
491
492#ifdef DEBUG_MEM
493 mem_intermediaire = MemoryUsage();
494 mem_diff = mem_intermediaire-mem_first;
495 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After fpNavigator->SetNavigatorState with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
496#endif
497
498 //-------
499 // Finale
500 //-------
501
502 // Update 'TrackLength' and remeber the Step length of the current Step
503 fpTrack->AddTrackLength(fpStep->GetStepLength());
505
506//#ifdef G4VERBOSE
507// // !!!!! Verbose
508// if(fpVerbose) fpVerbose->StepInfo();
509//#endif
510
511#ifdef G4VERBOSE
512 if(fpVerbose != nullptr) fpVerbose->PostStepVerbose(fpTrack);
513#endif
514
515// G4cout << " G4ITStepProcessor::DoStepping -- " <<fpTrack->GetTrackID() << " tps = " << fpTrack->GetGlobalTime() << G4endl;
516
517 // Send G4Step information to Hit/Dig if the volume is sensitive
518 /***
519 fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
520 StepControlFlag = fpStep->GetControlFlag();
521
522 if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
523 {
524 fpSensitive = fpStep->GetPreStepPoint()->
525 GetSensitiveDetector();
526 if( fpSensitive != 0 )
527 {
528 fpSensitive->Hit(fpStep);
529 }
530 }
531
532 User intervention process.
533 if( fpUserSteppingAction != 0 )
534 {
535 fpUserSteppingAction->UserSteppingAction(fpStep);
536 }
537 G4UserSteppingAction* regionalAction
538 = fpStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()
539 ->GetRegionalSteppingAction();
540 if( regionalAction ) regionalAction->UserSteppingAction(fpStep);
541 ***/
542 fpTrackingManager->AppendStep(fpTrack, fpStep);
543 // Stepping process finish. Return the value of the StepStatus.
544
545#ifdef DEBUG_MEM
546 MemStat mem_intermediaire = MemoryUsage();
547 mem_diff = mem_intermediaire-mem_first;
548 G4cout << "\t\t\t >> || MEM || End of DoStepping() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
549#endif
550
551 // return fpState->fStepStatus;
552}
@ fAtRestDoItProc
void AppendStep(G4Track *track, G4Step *step)
void SetSafety(const G4double aValue)
void SetStepStatus(const G4StepStatus aValue)
void SetStepLength(G4double value)
G4double GetStepLength() const
void SetStepLength(G4double value)
void AddTrackLength(const G4double aValue)
void IncrementCurrentStepNumber()
G4ITNavigatorState_Lock * GetNavigatorState() const
virtual void StepInfoForLeadingTrack()=0
virtual void PostStepDoItAllDone()=0
virtual void PreStepVerbose(G4Track *)=0
virtual void PostStepVerbose(G4Track *)=0
virtual void AtRestDoItInvoked()=0
virtual void NewStep()=0
virtual void AlongStepDoItAllDone()=0
G4ProcessVector * fpAtRestDoItVector

Referenced by Stepping().

◆ ExtractDoItData()

void G4ITStepProcessor::ExtractDoItData ( )

Definition at line 191 of file G4ITStepProcessor2.cc.

192{
193 if (fpTrack == nullptr)
194 {
196 return;
197 }
198
199 G4TrackStatus status = fpTrack->GetTrackStatus();
200
201 switch (status)
202 {
203 case fAlive:
204 case fStopButAlive:
205 case fSuspend:
207 default:
209 break;
210
211 case fStopAndKill:
214// G4TrackList::Pop(fpTrack);
215 fpTrackingManager->EndTracking(fpTrack);
216// fTrackContainer->PushToKill(fpTrack);
217 break;
218
221 if (fpSecondary != nullptr)
222 {
223 for (auto & i : *fpSecondary)
224 {
225 delete i;
226 }
227 fpSecondary->clear();
228 }
229// G4TrackList::Pop(fpTrack);
230 fpTrackingManager->EndTracking(fpTrack);
231// fTrackContainer->PushToKill(fpTrack);
232 break;
233 }
234
236}
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fPostponeToNextEvent
void RemoveReactionSet(G4Track *track)
static G4ITReactionSet * Instance()

Referenced by DoIt().

◆ ExtractILData()

void G4ITStepProcessor::ExtractILData ( )
protected

Definition at line 630 of file G4ITStepProcessor.cc.

631{
632 assert(fpTrack != 0);
633 if (fpTrack == nullptr)
634 {
636 return;
637 }
638
639 // assert(fpTrack->GetTrackStatus() != fStopAndKill);
640
641 if (fpTrack->GetTrackStatus() == fStopAndKill)
642 {
643// trackContainer->GetMainList()->pop(fpTrack);
644 fpTrackingManager->EndTracking(fpTrack);
646 return;
647 }
648
649 if (IsInf(fTimeStep))
650 {
651 // G4cout << "!!!!!!!!!!!! IS INF " << track->GetTrackID() << G4endl;
653 return;
654 }
655 if (fTimeStep < fILTimeStep - DBL_EPSILON)
656 {
657 // G4cout << "!!!!!!!!!!!! TEMPS DIFFERENTS "
658 // << track->GetTrackID() << G4endl;
659
660 fLeadingTracks.Reset();
661
662 fILTimeStep = GetInteractionTime();
663
664// G4cout << "Will set leading step to true for time :"
665// << SP -> GetInteractionTime() << " against fTimeStep : "
666// << G4BestUnit(fILTimeStep, "Time") << " the trackID is : " << track->GetTrackID()
667// << G4endl;
668
669// GetIT(fpTrack)->GetTrackingInfo()->SetLeadingStep(true);
670 fLeadingTracks.Push(fpTrack);
671 }
672 else if(fabs(fILTimeStep - fTimeStep) < DBL_EPSILON )
673 {
674
675 // G4cout << "!!!!!!!!!!!! MEME TEMPS " << track->GetTrackID() << G4endl;
676 // G4cout << "Will set leading step to true for time :"
677 // << SP -> GetInteractionTime() << " against fTimeStep : "
678 // << fTimeStep << " the trackID is : " << track->GetTrackID()<< G4endl;
679// GetIT(fpTrack)->GetTrackingInfo()->SetLeadingStep(true);
680 fLeadingTracks.Push(fpTrack);
681 }
682 // else
683 // {
684 // G4cout << "!!!! Bigger time : " << "currentTime : "<<fILTimeStep
685 // << " proposedTime : " << SP -> GetInteractionTime() << G4endl;
686 // }
687
689}
bool IsInf(T value)
void Push(G4Track *)
#define DBL_EPSILON
Definition templates.hh:66

Referenced by ComputeInteractionLength().

◆ FindTransportationStep()

void G4ITStepProcessor::FindTransportationStep ( )

Definition at line 801 of file G4ITStepProcessor2.cc.

802{
803 double physicalStep(0.);
804
805 fpTransportation = fpProcessInfo->fpTransportation;
806 // dynamic_cast<G4ITTransportation*>((fpProcessInfo->fpAlongStepGetPhysIntVector)[MAXofAlongStepLoops-1]);
807
808 if(fpTrack == nullptr)
809 {
810 G4ExceptionDescription exceptionDescription;
811 exceptionDescription << "No G4ITStepProcessor::fpTrack found";
812 G4Exception("G4ITStepProcessor::FindTransportationStep",
813 "ITStepProcessor0013",
815 exceptionDescription);
816 return;
817
818 }
819 if(fpITrack == nullptr)
820 {
821 G4ExceptionDescription exceptionDescription;
822 exceptionDescription << "No G4ITStepProcessor::fITrack";
823 G4Exception("G4ITStepProcessor::FindTransportationStep",
824 "ITStepProcessor0014",
826 exceptionDescription);
827 return;
828 }
829 if((fpITrack->GetTrack()) == nullptr)
830 {
831 G4ExceptionDescription exceptionDescription;
832 exceptionDescription << "No G4ITStepProcessor::fITrack->GetTrack()";
833 G4Exception("G4ITStepProcessor::FindTransportationStep",
834 "ITStepProcessor0015",
836 exceptionDescription);
837 return;
838 }
839
840 if(fpTransportation != nullptr)
841 {
842 fpTransportation->SetProcessState(fpTrackingInfo->GetProcessState(fpTransportation
843 ->GetProcessID()));
844 fpTransportation->ComputeStep(*fpTrack, *fpStep, fTimeStep, physicalStep);
845
846// fpTransportation->SetProcessState(0);
847 fpTransportation->ResetProcessState();
848 }
849
850 if(physicalStep >= DBL_MAX)
851 {
852// G4cout << "---- 2) Setting stop and kill for " << GetIT(fpTrack)->GetName() << G4endl;
853 fpTrack -> SetTrackStatus(fStopAndKill);
854 return;
855 }
856
857 fpState->fPhysicalStep = physicalStep;
858}
G4Track * GetTrack()
Definition G4IT.hh:218
G4ITTransportation * fpTransportation

Referenced by DoStepping().

◆ ForceReInitialization()

void G4ITStepProcessor::ForceReInitialization ( )

Definition at line 184 of file G4ITStepProcessor.cc.

185{
186 fInitialized = false;
188 Initialize();
189}
virtual void Initialize()

◆ GetAtRestDoItProcTriggered()

std::size_t G4ITStepProcessor::GetAtRestDoItProcTriggered ( ) const
inline

Definition at line 219 of file G4ITStepProcessor.hh.

220 {
221 return fAtRestDoItProcTriggered;
222 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetAtRestIL()

void G4ITStepProcessor::GetAtRestIL ( )
protected

Definition at line 535 of file G4ITStepProcessor.cc.

536{
537 // Select the rest process which has the shortest time before
538 // it is invoked. In rest processes, GPIL()
539 // returns the time before a process occurs.
540 G4double lifeTime(DBL_MAX), shortestLifeTime (DBL_MAX);
541
542 fAtRestDoItProcTriggered = 0;
543 shortestLifeTime = DBL_MAX;
544
545 unsigned int NofInactiveProc=0;
546
547 for( G4int ri=0; ri < (G4int)fpProcessInfo->MAXofAtRestLoops; ++ri )
548 {
549 fpCurrentProcess = dynamic_cast<G4VITProcess*>((*fpProcessInfo->fpAtRestGetPhysIntVector)[ri]);
550 if (fpCurrentProcess== nullptr)
551 {
552 (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
553 NofInactiveProc++;
554 continue;
555 } // NULL means the process is inactivated by a user on fly.
556
557 fCondition=NotForced;
558 fpCurrentProcess->SetProcessState(
559 fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
560
561 lifeTime = fpCurrentProcess->AtRestGPIL( *fpTrack, &fCondition );
562 fpCurrentProcess->ResetProcessState();
563
564 if(fCondition==Forced)
565 {
566 (fpState->fSelectedAtRestDoItVector)[ri] = Forced;
567 }
568 else
569 {
570 (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
571 if(lifeTime < shortestLifeTime )
572 {
573 shortestLifeTime = lifeTime;
574 fAtRestDoItProcTriggered = G4int(ri);
575 }
576 }
577 }
578
579 (fpState->fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
580
581// G4cout << " --> Selected at rest process : "
582// << (*fpProcessInfo->fpAtRestGetPhysIntVector)[fAtRestDoItProcTriggered]
583// ->GetProcessName()
584// << G4endl;
585
586 fTimeStep = shortestLifeTime;
587
588 // at least one process is necessary to destroy the particle
589 // exit with warning
590 if(NofInactiveProc==fpProcessInfo->MAXofAtRestLoops)
591 {
592 G4cerr << "ERROR - G4ITStepProcessor::InvokeAtRestDoItProcs()" << G4endl
593 << " No AtRestDoIt process is active!" << G4endl;
594 }
595}
G4SelectedAtRestDoItVector fSelectedAtRestDoItVector
G4double AtRestGPIL(const G4Track &track, G4ForceCondition *condition)
G4ProcessVector * fpAtRestGetPhysIntVector

Referenced by DoDefinePhysicalStepLength().

◆ GetCondition()

G4ForceCondition G4ITStepProcessor::GetCondition ( ) const
inline

Definition at line 284 of file G4ITStepProcessor.hh.

285 {
286 return fCondition;
287 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetCurrentProcess()

const G4VITProcess * G4ITStepProcessor::GetCurrentProcess ( ) const
inline

Definition at line 244 of file G4ITStepProcessor.hh.

245 {
246 return fpCurrentProcess;
247 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetCurrentProcessInfo()

const ProcessGeneralInfo * G4ITStepProcessor::GetCurrentProcessInfo ( ) const
inline

Definition at line 264 of file G4ITStepProcessor.hh.

265 {
266 return fpProcessInfo;
267 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetCurrentVolume()

const G4VPhysicalVolume * G4ITStepProcessor::GetCurrentVolume ( ) const
inline

Definition at line 279 of file G4ITStepProcessor.hh.

280 {
281 return fpCurrentVolume;
282 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetGPILSelection()

G4GPILSelection G4ITStepProcessor::GetGPILSelection ( ) const
inline

Definition at line 224 of file G4ITStepProcessor.hh.

225 {
226 return fGPILSelection;
227 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetILTimeStep()

G4double G4ITStepProcessor::GetILTimeStep ( )
inline

Definition at line 202 of file G4ITStepProcessor.hh.

203 {
204 return fILTimeStep;
205 }

◆ GetInteractionTime()

double G4ITStepProcessor::GetInteractionTime ( )
inline

Definition at line 490 of file G4ITStepProcessor.hh.

491{
492 return fTimeStep;
493}

Referenced by ExtractILData().

◆ GetN2ndariesAlongStepDoIt()

G4int G4ITStepProcessor::GetN2ndariesAlongStepDoIt ( ) const
inline

Definition at line 229 of file G4ITStepProcessor.hh.

230 {
231 return fN2ndariesAlongStepDoIt;
232 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetN2ndariesAtRestDoIt()

G4int G4ITStepProcessor::GetN2ndariesAtRestDoIt ( ) const
inline

Definition at line 234 of file G4ITStepProcessor.hh.

235 {
236 return fN2ndariesAtRestDoIt;
237 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetN2ndariesPostStepDoIt()

G4int G4ITStepProcessor::GetN2ndariesPostStepDoIt ( ) const
inline

Definition at line 239 of file G4ITStepProcessor.hh.

240 {
241 return fN2ndariesPostStepDoIt;
242 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetParticleChange()

const G4VParticleChange * G4ITStepProcessor::GetParticleChange ( ) const
inline

Definition at line 274 of file G4ITStepProcessor.hh.

275 {
276 return fpParticleChange;
277 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetPhysIntLength()

G4double G4ITStepProcessor::GetPhysIntLength ( ) const
inline

Definition at line 249 of file G4ITStepProcessor.hh.

250 {
251 return fPhysIntLength;
252 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetPostStepAtTimeDoItProcTriggered()

std::size_t G4ITStepProcessor::GetPostStepAtTimeDoItProcTriggered ( ) const
inline

Definition at line 254 of file G4ITStepProcessor.hh.

255 {
256 return fPostStepAtTimeDoItProcTriggered;
257 }

◆ GetPostStepDoItProcTriggered()

std::size_t G4ITStepProcessor::GetPostStepDoItProcTriggered ( ) const
inline

Definition at line 259 of file G4ITStepProcessor.hh.

260 {
261 return fPostStepDoItProcTriggered;
262 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetProcessInfo()

void G4ITStepProcessor::GetProcessInfo ( )
protected

Definition at line 482 of file G4ITStepProcessor.cc.

483{
484 G4ParticleDefinition* particle = fpTrack->GetDefinition();
485 auto it =
486 fProcessGeneralInfoMap.find(particle);
487
488 if(it == fProcessGeneralInfoMap.end())
489 {
491 fpTrack->GetDefinition()->GetProcessManager());
492 if(fpProcessInfo == nullptr)
493 {
494 G4ExceptionDescription exceptionDescription("...");
495 G4Exception("G4ITStepProcessor::GetProcessNumber",
496 "ITStepProcessor0008",
498 exceptionDescription);
499 return;
500 }
501 }
502 else
503 {
504 fpProcessInfo = it->second;
505 }
506}
void SetupGeneralProcessInfo(G4ParticleDefinition *, G4ProcessManager *)

Referenced by SetupMembers().

◆ GetProcessorState()

const G4ITStepProcessorState * G4ITStepProcessor::GetProcessorState ( ) const
inline

Definition at line 269 of file G4ITStepProcessor.hh.

270 {
271 return fpState;
272 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetSecondaries()

G4TrackVector * G4ITStepProcessor::GetSecondaries ( ) const
inline

Definition at line 178 of file G4ITStepProcessor.hh.

179 {
180 return fpSecondary;
181 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetStep() [1/2]

G4Step * G4ITStepProcessor::GetStep ( )
inline

Definition at line 165 of file G4ITStepProcessor.hh.

166 {
167 return fpStep;
168 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetStep() [2/2]

const G4Step * G4ITStepProcessor::GetStep ( ) const
inline

Definition at line 169 of file G4ITStepProcessor.hh.

170 {
171 return fpStep;
172 }

◆ GetTrack() [1/2]

G4Track * G4ITStepProcessor::GetTrack ( )
inline

Definition at line 161 of file G4ITStepProcessor.hh.

162 {
163 return fpTrack;
164 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetTrack() [2/2]

const G4Track * G4ITStepProcessor::GetTrack ( ) const
inline

Definition at line 433 of file G4ITStepProcessor.hh.

434{
435 return fpTrack;
436}

◆ GetTrackingManager()

G4ITTrackingManager * G4ITStepProcessor::GetTrackingManager ( )
inline

Definition at line 186 of file G4ITStepProcessor.hh.

187 {
188 return fpTrackingManager;
189 }

◆ InitDefineStep()

void G4ITStepProcessor::InitDefineStep ( )
protected

ADDED BACK

ADDED BACK

Definition at line 845 of file G4ITStepProcessor.cc.

846{
847
848 if(fpStep == nullptr)
849 {
850 // Create new Step and give it to the track
851 fpStep = new G4Step();
852 fpTrack->SetStep(fpStep);
853 fpSecondary = fpStep->NewSecondaryVector();
854
855 // Create new state and set it in the trackingInfo
856 fpState = new G4ITStepProcessorState();
858
859 SetupMembers();
861
862 fpTrackingManager->StartTracking(fpTrack);
863 }
864 else
865 {
866 SetupMembers();
867
868 fpState->fPreviousStepSize = fpTrack->GetStepLength();
869 /***
870 // Send G4Step information to Hit/Dig if the volume is sensitive
871 fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
872 StepControlFlag = fpStep->GetControlFlag();
873 if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
874 {
875 fpSensitive = fpStep->GetPreStepPoint()->
876 GetSensitiveDetector();
877
878 // if( fSensitive != 0 ) {
879 // fSensitive->Hit(fStep);
880 // }
881 }
882 ***/
883 // Store last PostStepPoint to PreStepPoint, and swap current and next
884 // volume information of G4Track. Reset total energy deposit in one Step.
885 fpStep->CopyPostToPreStepPoint();
886 fpStep->ResetTotalEnergyDeposit();
887
888 //JA Set the volume before it is used (in DefineStepLength() for User Limit)
889 fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
890 /*
891 G4cout << G4endl;
892 G4cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!" << G4endl;
893 G4cout << "PreStepPoint Volume : "
894 << fpCurrentVolume->GetName() << G4endl;
895 G4cout << "Track Touchable : "
896 << fpTrack->GetTouchableHandle()->GetVolume()->GetName() << G4endl;
897 G4cout << "Track NextTouchable : "
898 << fpTrack->GetNextTouchableHandle()->GetVolume()->GetName()
899 << G4endl;
900 */
901 // Reset the step's auxiliary points vector pointer
903
904 // Switch next touchable in track to current one
905 fpTrack->SetTouchableHandle(fpTrack->GetNextTouchableHandle());
906 fpState->fTouchableHandle = fpTrack->GetTouchableHandle();
907 fpTrack->SetNextTouchableHandle(fpState->fTouchableHandle);
908
909 //! ADDED BACK
910 /*
911 G4VPhysicalVolume* oldTopVolume =
912 fpTrack->GetTouchableHandle()->GetVolume();
913 fpNavigator->SetNavigatorState(
914 fpITrack->GetTrackingInfo()->GetNavigatorState());
915
916 G4VPhysicalVolume* newTopVolume = fpNavigator->ResetHierarchyAndLocate(
917 fpTrack->GetPosition(), fpTrack->GetMomentumDirection(),
918 *((G4TouchableHistory*) fpTrack->GetTouchableHandle()()));
919
920 // G4VPhysicalVolume* newTopVolume=
921 // fpNavigator->LocateGlobalPointAndSetup(fpTrack->GetPosition(),
922 // &fpTrack->GetMomentumDirection(),
923 // true, false);
924
925 // G4cout << "New Top Volume : " << newTopVolume->GetName() << G4endl;
926
927 if (newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId()
928 == 1)
929 {
930 fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
931 fpTrack->SetTouchableHandle(fpState->fTouchableHandle);
932 fpTrack->SetNextTouchableHandle(fpState->fTouchableHandle);
933 }
934
935 */
936 //! ADDED BACK
937 //==========================================================================
938 // Only reset navigator state + reset volume hierarchy (internal)
939 // No need to relocate
940 //==========================================================================
941 fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()
943 }
944}
void StartTracking(G4Track *)
G4VPhysicalVolume * GetPhysicalVolume() const
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *vec)
void ResetTotalEnergyDeposit()
void CopyPostToPreStepPoint()
G4StepPoint * GetPreStepPoint() const
G4TrackVector * NewSecondaryVector()
void SetStep(const G4Step *aValue)
const G4TouchableHandle & GetNextTouchableHandle() const
void SetNextTouchableHandle(const G4TouchableHandle &apValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4double GetStepLength() const
void SetStepProcessorState(G4ITStepProcessorState_Lock *)

Referenced by DoDefinePhysicalStepLength().

◆ Initialize()

void G4ITStepProcessor::Initialize ( )
virtual

Definition at line 193 of file G4ITStepProcessor.cc.

194{
196 if(fInitialized) return;
197 // ActiveOnlyITProcess();
198
200 ->GetNavigatorForTracking());
201
202 fPhysIntLength = DBL_MAX;
203 kCarTolerance = 0.5
205
206 if(fpVerbose == nullptr)
207 {
208 G4ITTrackingInteractivity* interactivity = fpTrackingManager->GetInteractivity();
209
210 if(interactivity != nullptr)
211 {
212 fpVerbose = interactivity->GetSteppingVerbose();
213 fpVerbose->SetStepProcessor(this);
214 }
215 }
216
217 fpTrackContainer = G4ITTrackHolder::Instance();
218
219 fInitialized = true;
220}
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
void SetNavigator(G4ITNavigator *value)
static G4ITTrackHolder * Instance()
G4VITSteppingVerbose * GetSteppingVerbose()
G4ITTrackingInteractivity * GetInteractivity()
static G4ITTransportationManager * GetTransportationManager()
void SetStepProcessor(const G4ITStepProcessor *stepProcessor)

Referenced by ForceReInitialization(), and G4Scheduler::Process().

◆ InvokeAlongStepDoItProcs()

void G4ITStepProcessor::InvokeAlongStepDoItProcs ( )
protected

Definition at line 627 of file G4ITStepProcessor2.cc.

628{
629
630#ifdef DEBUG_MEM
631 MemStat mem_first, mem_second, mem_diff;
632#endif
633
634#ifdef DEBUG_MEM
635 mem_first = MemoryUsage();
636#endif
637
638 // If the current Step is defined by a 'ExclusivelyForced'
639 // PostStepDoIt, then don't invoke any AlongStepDoIt
640 if(fpState->fStepStatus == fExclusivelyForcedProc)
641 {
642 return; // Take note 'return' at here !!!
643 }
644
645 // Invoke the all active continuous processes
646 for(std::size_t ci = 0; ci < fpProcessInfo->MAXofAlongStepLoops; ++ci)
647 {
648 fpCurrentProcess =
649 (G4VITProcess*) (*fpProcessInfo->fpAlongStepDoItVector)[(G4int)ci];
650 if(fpCurrentProcess == nullptr) continue;
651 // NULL means the process is inactivated by a user on fly.
652
653 fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
654 ->GetProcessID()));
655 fpParticleChange = fpCurrentProcess->AlongStepDoIt(*fpTrack, *fpStep);
656
657#ifdef DEBUG_MEM
658 MemStat mem_intermediaire = MemoryUsage();
659 mem_diff = mem_intermediaire-mem_first;
660 G4cout << "\t\t\t >> || MEM || After calling AlongStepDoIt for " << fpCurrentProcess->GetProcessName() << " and track "<< fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
661#endif
662
663// fpCurrentProcess->SetProcessState(0);
664 fpCurrentProcess->ResetProcessState();
665 // Update the PostStepPoint of Step according to ParticleChange
666
667 fpParticleChange->UpdateStepForAlongStep(fpStep);
668
669#ifdef G4VERBOSE
670 // !!!!! Verbose
671 if(fpVerbose != nullptr) fpVerbose->AlongStepDoItOneByOne();
672#endif
673
674 // Now Store the secondaries from ParticleChange to SecondaryList
675 DealWithSecondaries(fN2ndariesAlongStepDoIt);
676
677 // Set the track status according to what the process defined
678 // if kinetic energy >0, otherwise set fStopButAlive
679 fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
680
681 // clear ParticleChange
682 fpParticleChange->Clear();
683 }
684
685#ifdef DEBUG_MEM
686 MemStat mem_intermediaire = MemoryUsage();
687 mem_diff = mem_intermediaire-mem_first;
688 G4cout << "\t\t\t >> || MEM || After looping on processes with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
689#endif
690
691 fpStep->UpdateTrack();
692
693 G4TrackStatus fNewStatus = fpTrack->GetTrackStatus();
694
695 if(fNewStatus == fAlive && fpTrack->GetKineticEnergy() <= DBL_MIN)
696 {
697 // G4cout << "G4ITStepProcessor::InvokeAlongStepDoItProcs : Track will be killed" << G4endl;
698 if(fpProcessInfo->MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
699 else fNewStatus = fStopAndKill;
700 fpTrack->SetTrackStatus( fNewStatus );
701 }
702
703}
void DealWithSecondaries(G4int &)
void UpdateTrack()
virtual void AlongStepDoItOneByOne()=0
virtual G4Step * UpdateStepForAlongStep(G4Step *Step)
G4TrackStatus GetTrackStatus() const
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
const G4String & GetProcessName() const
G4ProcessVector * fpAlongStepDoItVector

Referenced by DoStepping().

◆ InvokeAtRestDoItProcs()

void G4ITStepProcessor::InvokeAtRestDoItProcs ( )
protected

Definition at line 560 of file G4ITStepProcessor2.cc.

561{
562 fpStep->SetStepLength(0.); //the particle has stopped
563 fpTrack->SetStepLength(0.);
564
565 G4SelectedAtRestDoItVector& selectedAtRestDoItVector =
567
568 // invoke selected process
569 for(std::size_t np = 0; np < fpProcessInfo->MAXofAtRestLoops; ++np)
570 {
571 //
572 // Note: DoItVector has inverse order against GetPhysIntVector
573 // and SelectedAtRestDoItVector.
574 //
575 if(selectedAtRestDoItVector[fpProcessInfo->MAXofAtRestLoops - np - 1] != InActivated)
576 {
577 fpCurrentProcess =
578 (G4VITProcess*) (*fpProcessInfo->fpAtRestDoItVector)[(G4int)np];
579
580// G4cout << " Invoke : "
581// << fpCurrentProcess->GetProcessName()
582// << G4endl;
583
584 // if(fpVerbose)
585 // {
586 // fpVerbose->AtRestDoItOneByOne();
587 // }
588
589 fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
590 ->GetProcessID()));
591 fpParticleChange = fpCurrentProcess->AtRestDoIt(*fpTrack, *fpStep);
592 fpCurrentProcess->ResetProcessState();
593
594 // Set the current process as a process which defined this Step length
595 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpCurrentProcess);
596
597 // Update Step
598 fpParticleChange->UpdateStepForAtRest(fpStep);
599
600 // Now Store the secondaries from ParticleChange to SecondaryList
601 DealWithSecondaries(fN2ndariesAtRestDoIt);
602
603 // Set the track status according to what the process defined
604 // if kinetic energy >0, otherwise set fStopButAlive
605 fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
606
607 // clear ParticleChange
608 fpParticleChange->Clear();
609
610 } //if(fSelectedAtRestDoItVector[np] != InActivated){
611 } //for(std::size_t np=0; np < MAXofAtRestLoops; ++np){
612 fpStep->UpdateTrack();
613
614 // Modification par rapport au transport standard :
615 // fStopAndKill doit etre propose par le modele
616 // sinon d autres processus AtRest seront appeles
617 // au pas suivant
618 // fpTrack->SetTrackStatus(fStopAndKill);
619}
std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0

Referenced by DoStepping().

◆ InvokePostStepDoItProcs()

void G4ITStepProcessor::InvokePostStepDoItProcs ( )
protected

Definition at line 711 of file G4ITStepProcessor2.cc.

712{
713 std::size_t _MAXofPostStepLoops = fpProcessInfo->MAXofPostStepLoops;
714 G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState
716 G4StepStatus& stepStatus = fpState->fStepStatus;
717
718 // Invoke the specified discrete processes
719 for(std::size_t np = 0; np < _MAXofPostStepLoops; ++np)
720 {
721 //
722 // Note: DoItVector has inverse order against GetPhysIntVector
723 // and SelectedPostStepDoItVector.
724 //
725 G4int Cond = selectedPostStepDoItVector[_MAXofPostStepLoops
726 - np - 1];
727 if(Cond != InActivated)
728 {
729 if(((Cond == NotForced) && (stepStatus == fPostStepDoItProc)) ||
730 ((Cond == Forced) && (stepStatus != fExclusivelyForcedProc))
731 ||
732 // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
733 ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc))
734 || ((Cond == StronglyForced)))
735 {
736
737 InvokePSDIP(np);
738 }
739 } //if(*fSelectedPostStepDoItVector(np)........
740
741 // Exit from PostStepLoop if the track has been killed,
742 // but extra treatment for processes with Strongly Forced flag
743 if(fpTrack->GetTrackStatus() == fStopAndKill)
744 {
745 for(std::size_t np1 = np + 1; np1 < _MAXofPostStepLoops; ++np1)
746 {
747 G4int Cond2 = selectedPostStepDoItVector[_MAXofPostStepLoops
748 - np1 - 1];
749 if(Cond2 == StronglyForced)
750 {
751 InvokePSDIP(np1);
752 }
753 }
754 break;
755 }
756 } //for(std::size_t np=0; np < MAXofPostStepLoops; ++np){
757}
std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
G4StepStatus
void InvokePSDIP(std::size_t)

Referenced by DoStepping().

◆ InvokePSDIP()

void G4ITStepProcessor::InvokePSDIP ( std::size_t np)
protected

Definition at line 761 of file G4ITStepProcessor2.cc.

762{
763 fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpPostStepDoItVector)[(G4int)np];
764
765 fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
766 ->GetProcessID()));
767 fpParticleChange = fpCurrentProcess->PostStepDoIt(*fpTrack, *fpStep);
768// fpCurrentProcess->SetProcessState(0);
769 fpCurrentProcess->ResetProcessState();
770
771 // Update PostStepPoint of Step according to ParticleChange
772 fpParticleChange->UpdateStepForPostStep(fpStep);
773
774#ifdef G4VERBOSE
775 // !!!!! Verbose
776 if(fpVerbose != nullptr) fpVerbose->PostStepDoItOneByOne();
777#endif
778
779 // Update G4Track according to ParticleChange after each PostStepDoIt
780 fpStep->UpdateTrack();
781
782 // Update safety after each invocation of PostStepDoIts
784
785 // Now Store the secondaries from ParticleChange to SecondaryList
786 DealWithSecondaries(fN2ndariesPostStepDoIt);
787
788 // Set the track status according to what the process defined
789 fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
790
791 // clear ParticleChange
792 fpParticleChange->Clear();
793}
virtual void PostStepDoItOneByOne()=0
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
G4ProcessVector * fpPostStepDoItVector

Referenced by InvokePostStepDoItProcs(), and InvokeTransportationProc().

◆ InvokeTransportationProc()

void G4ITStepProcessor::InvokeTransportationProc ( )
protected

Definition at line 862 of file G4ITStepProcessor2.cc.

863{
864 std::size_t _MAXofPostStepLoops = fpProcessInfo->MAXofPostStepLoops;
865 G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState
867 G4StepStatus& stepStatus = fpState->fStepStatus;
868
869 // Invoke the specified discrete processes
870 for(std::size_t np = 0; np < _MAXofPostStepLoops; ++np)
871 {
872 //
873 // Note: DoItVector has inverse order against GetPhysIntVector
874 // and SelectedPostStepDoItVector.
875 //
876 G4int Cond = selectedPostStepDoItVector[_MAXofPostStepLoops - np - 1];
877 if(Cond != InActivated)
878 {
879 if(((Cond == Forced) && (stepStatus != fExclusivelyForcedProc)) ||
880 // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
881 ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc))
882 || ((Cond == StronglyForced)))
883 {
884
885 InvokePSDIP(np);
886 }
887 } //if(Cond != InActivated)
888
889 // Exit from PostStepLoop if the track has been killed,
890 // but extra treatment for processes with Strongly Forced flag
891 if(fpTrack->GetTrackStatus() == fStopAndKill)
892 {
893 for(std::size_t np1 = np + 1; np1 < _MAXofPostStepLoops; ++np1)
894 {
895 G4int Cond2 = selectedPostStepDoItVector[_MAXofPostStepLoops - np1 - 1];
896 if(Cond2 == StronglyForced)
897 {
898 InvokePSDIP(np1);
899 }
900 }
901 break;
902 }
903 }
904}

Referenced by DoStepping().

◆ operator=()

G4ITStepProcessor & G4ITStepProcessor::operator= ( const G4ITStepProcessor & other)
protected

Definition at line 275 of file G4ITStepProcessor.cc.

276{
277 if(this == &rhs) return *this; // handle self assignment
278 //assignment operator
279 return *this;
280}

◆ PrepareLeadingTracks()

void G4ITStepProcessor::PrepareLeadingTracks ( )

Definition at line 268 of file G4ITStepProcessor.cc.

269{
270 fLeadingTracks.PrepareLeadingTracks();
271}

Referenced by G4Scheduler::Stepping().

◆ PushSecondaries()

void G4ITStepProcessor::PushSecondaries ( )
protected

Definition at line 240 of file G4ITStepProcessor2.cc.

241{
242 if ((fpSecondary == nullptr) || fpSecondary->empty())
243 {
244 // DEBUG
245 // G4cout << "NO SECONDARIES !!! " << G4endl;
246 return;
247 }
248
249 // DEBUG
250 // G4cout << "There are secondaries : "<< secondaries -> size() << G4endl ;
251
252 auto secondaries_i = fpSecondary->begin();
253
254 for (; secondaries_i != fpSecondary->end(); ++secondaries_i)
255 {
256 G4Track* secondary = *secondaries_i;
257 fpTrackContainer->_PushTrack(secondary);
258 }
259}
void _PushTrack(G4Track *track)

Referenced by ExtractDoItData().

◆ ResetLeadingTracks()

void G4ITStepProcessor::ResetLeadingTracks ( )

Definition at line 261 of file G4ITStepProcessor.cc.

262{
263 fLeadingTracks.Reset();
264}

Referenced by G4Scheduler::Stepping().

◆ ResetSecondaries()

void G4ITStepProcessor::ResetSecondaries ( )
protected

Definition at line 525 of file G4ITStepProcessor.cc.

526{
527 // Reset the secondary particles
528 fN2ndariesAtRestDoIt = 0;
529 fN2ndariesAlongStepDoIt = 0;
530 fN2ndariesPostStepDoIt = 0;
531}

Referenced by G4ITStepProcessor(), G4ITStepProcessor(), and SetupMembers().

◆ SetInitialStep()

void G4ITStepProcessor::SetInitialStep ( )
protected

Definition at line 701 of file G4ITStepProcessor.cc.

702{
703 // DEBUG
704 // G4cout << "SetInitialStep for : " << fpITrack-> GetName() << G4endl;
705 //________________________________________________________
706 // Initialize geometry
707
708 if(!fpTrack->GetTouchableHandle())
709 {
710 //==========================================================================
711 // Create navigator state and Locate particle in geometry
712 //==========================================================================
713 /*
714 fpNavigator->NewNavigatorStateAndLocate(fpTrack->GetPosition(),
715 fpTrack->GetMomentumDirection());
716
717 fpITrack->GetTrackingInfo()->
718 SetNavigatorState(fpNavigator->GetNavigatorState());
719 */
720 fpNavigator->NewNavigatorState();
721 fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator
722 ->GetNavigatorState());
723
724 G4ThreeVector direction = fpTrack->GetMomentumDirection();
725 fpNavigator->LocateGlobalPointAndSetup(fpTrack->GetPosition(),
726 &direction,
727 false,
728 false); // was false, false
729
730 fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
731
732 fpTrack->SetTouchableHandle(fpState->fTouchableHandle);
733 fpTrack->SetNextTouchableHandle(fpState->fTouchableHandle);
734 }
735 else
736 {
737 fpState->fTouchableHandle = fpTrack->GetTouchableHandle();
738 fpTrack->SetNextTouchableHandle(fpState->fTouchableHandle);
739
740 //==========================================================================
741 // Create OR set navigator state
742 //==========================================================================
743
744 if(fpITrack->GetTrackingInfo()->GetNavigatorState() != nullptr)
745 {
746 fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()
748 fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator
749 ->GetNavigatorState());
750 }
751 else
752 {
753 fpNavigator->NewNavigatorState(*((G4TouchableHistory*) fpState
754 ->fTouchableHandle()));
755 fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator
756 ->GetNavigatorState());
757 }
758
759 G4VPhysicalVolume* oldTopVolume =
760 fpTrack->GetTouchableHandle()->GetVolume();
761
762 //==========================================================================
763 // Locate particle in geometry
764 //==========================================================================
765
766// G4VPhysicalVolume* newTopVolume =
767// fpNavigator->LocateGlobalPointAndSetup(
768// fpTrack->GetPosition(),
769// &fpTrack->GetMomentumDirection(),
770// true, false);
771
772 G4VPhysicalVolume* newTopVolume =
773 fpNavigator->ResetHierarchyAndLocate(fpTrack->GetPosition(),
774 fpTrack->GetMomentumDirection(),
775 *((G4TouchableHistory*) fpTrack
776 ->GetTouchableHandle()()));
777
778 if(newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId()
779 == 1)
780 {
781 fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
782 fpTrack->SetTouchableHandle(fpState->fTouchableHandle);
783 fpTrack->SetNextTouchableHandle(fpState->fTouchableHandle);
784 }
785 }
786
787 fpCurrentVolume = fpState->fTouchableHandle->GetVolume();
788
789 //________________________________________________________
790 // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
791 // set the track state to 'Alive'.
792 if((fpTrack->GetTrackStatus() == fSuspend) || (fpTrack->GetTrackStatus()
794 {
795 fpTrack->SetTrackStatus(fAlive);
796 }
797
798 //HoangTRAN: it's better to check the status here
799 if(fpTrack->GetTrackStatus() == fStopAndKill) return;
800
801 // If the primary track has 'zero' kinetic energy, set the track
802 // state to 'StopButAlive'.
803 if(fpTrack->GetKineticEnergy() <= 0.0)
804 {
806 }
807 //________________________________________________________
808 // Set vertex information of G4Track at here
809 if(fpTrack->GetCurrentStepNumber() == 0)
810 {
811 fpTrack->SetVertexPosition(fpTrack->GetPosition());
813 fpTrack->SetVertexKineticEnergy(fpTrack->GetKineticEnergy());
815 }
816 //________________________________________________________
817 // If track is already outside the world boundary, kill it
818 if(fpCurrentVolume == nullptr)
819 {
820 // If the track is a primary, stop processing
821 if(fpTrack->GetParentID() == 0)
822 {
823 G4cerr << "ERROR - G4ITStepProcessor::SetInitialStep()" << G4endl<< " Primary particle starting at - "
824 << fpTrack->GetPosition()
825 << " - is outside of the world volume." << G4endl;
826 G4Exception("G4ITStepProcessor::SetInitialStep()", "ITStepProcessor0011",
827 FatalException, "Primary vertex outside of the world!");
828 }
829
830 fpTrack->SetTrackStatus( fStopAndKill );
831 G4cout << "WARNING - G4ITStepProcessor::SetInitialStep()" << G4endl
832 << " Initial track position is outside world! - "
833 << fpTrack->GetPosition() << G4endl;
834 }
835 else
836 {
837 // Initial set up for attribues of 'Step'
838 fpStep->InitializeStep( fpTrack );
839 }
840
841 fpState->fStepStatus = fUndefined;
842}
@ fUndefined
G4TouchableHandle fTouchableHandle
void InitializeStep(G4Track *aValue)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
G4VPhysicalVolume * GetVolume() const
void SetVertexPosition(const G4ThreeVector &aValue)
void SetVertexMomentumDirection(const G4ThreeVector &aValue)
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
void SetVertexKineticEnergy(const G4double aValue)
G4int GetParentID() const
void SetLogicalVolumeAtVertex(const G4LogicalVolume *)
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetRegularStructureId() const =0

Referenced by InitDefineStep().

◆ SetNavigator()

void G4ITStepProcessor::SetNavigator ( G4ITNavigator * value)
inlineprotected

Definition at line 449 of file G4ITStepProcessor.hh.

450{
451 fpNavigator = value;
452}

Referenced by Initialize().

◆ SetPreviousStepTime()

void G4ITStepProcessor::SetPreviousStepTime ( G4double previousTimeStep)
inline

Definition at line 426 of file G4ITStepProcessor.hh.

427{
428 fPreviousTimeStep = previousTimeStep;
429}

Referenced by ComputeInteractionLength().

◆ SetStep()

void G4ITStepProcessor::SetStep ( G4Step * val)
inline

Definition at line 173 of file G4ITStepProcessor.hh.

174 {
175 fpStep = val;
176 }

◆ SetTrack()

void G4ITStepProcessor::SetTrack ( G4Track * track)
protected

Definition at line 449 of file G4ITStepProcessor.cc.

450{
451 fpTrack = track;
452 if(fpTrack != nullptr)
453 {
454 fpITrack = GetIT(fpTrack);
455 fpStep = const_cast<G4Step*>(fpTrack->GetStep());
456
457 if(fpITrack != nullptr)
458 {
459 fpTrackingInfo = fpITrack->GetTrackingInfo();
460 }
461 else
462 {
463 fpTrackingInfo = nullptr;
464 G4cerr << "Track ID : " << fpTrack->GetTrackID() << G4endl;
465
467 errMsg << "No IT pointer was attached to the track you try to process.";
468 G4Exception("G4ITStepProcessor::SetTrack",
469 "ITStepProcessor0007",
471 errMsg);
472 }
473 }
474 else
475 {
476 fpITrack = nullptr;
477 fpStep = nullptr;
478 }
479}
const G4Step * GetStep() const

Referenced by DefinePhysicalStepLength(), and Stepping().

◆ SetTrackingManager()

void G4ITStepProcessor::SetTrackingManager ( G4ITTrackingManager * trackMan)
inline

Definition at line 182 of file G4ITStepProcessor.hh.

183 {
184 fpTrackingManager = trackMan;
185 }

Referenced by G4Scheduler::Initialize().

◆ SetupGeneralProcessInfo()

void G4ITStepProcessor::SetupGeneralProcessInfo ( G4ParticleDefinition * particle,
G4ProcessManager * pm )
protected

Definition at line 338 of file G4ITStepProcessor.cc.

340{
341
342#ifdef debug
343 G4cout<<"G4ITStepProcessor::GetProcessNumber: is called track"<<G4endl;
344#endif
345 if(pm == nullptr)
346 {
347 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl<< " ProcessManager is NULL for particle = "
348 << particle->GetParticleName() << ", PDG_code = "
349 << particle->GetPDGEncoding() << G4endl;
350 G4Exception("G4SteppingManager::GetProcessNumber()", "ITStepProcessor0002",
351 FatalException, "Process Manager is not found.");
352 return;
353 }
354
355 auto it =
356 fProcessGeneralInfoMap.find(particle);
357 if(it != fProcessGeneralInfoMap.end())
358 {
359 G4Exception("G4SteppingManager::SetupGeneralProcessInfo()",
360 "ITStepProcessor0003",
361 FatalException, "Process info already registered.");
362 return;
363 }
364
365 // here used as temporary
366 fpProcessInfo = new ProcessGeneralInfo();
367
368 // AtRestDoits
369 fpProcessInfo->MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
371 fpProcessInfo->fpAtRestGetPhysIntVector =
373#ifdef debug
374 G4cout << "G4ITStepProcessor::GetProcessNumber: #ofAtRest="
375 << fpProcessInfo->MAXofAtRestLoops << G4endl;
376#endif
377
378 // AlongStepDoits
379 fpProcessInfo->MAXofAlongStepLoops =
381 fpProcessInfo->fpAlongStepDoItVector =
383 fpProcessInfo->fpAlongStepGetPhysIntVector =
385#ifdef debug
386 G4cout << "G4ITStepProcessor::GetProcessNumber:#ofAlongStp="
387 << fpProcessInfo->MAXofAlongStepLoops << G4endl;
388#endif
389
390 // PostStepDoits
391 fpProcessInfo->MAXofPostStepLoops =
394 fpProcessInfo->fpPostStepGetPhysIntVector =
396#ifdef debug
397 G4cout << "G4ITStepProcessor::GetProcessNumber: #ofPostStep="
398 << fpProcessInfo->MAXofPostStepLoops << G4endl;
399#endif
400
401 if (SizeOfSelectedDoItVector<fpProcessInfo->MAXofAtRestLoops ||
402 SizeOfSelectedDoItVector<fpProcessInfo->MAXofAlongStepLoops ||
403 SizeOfSelectedDoItVector<fpProcessInfo->MAXofPostStepLoops )
404 {
405 G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl
406 << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
407 << " ; is smaller then one of MAXofAtRestLoops= "
408 << fpProcessInfo->MAXofAtRestLoops << G4endl
409 << " or MAXofAlongStepLoops= " << fpProcessInfo->MAXofAlongStepLoops
410 << " or MAXofPostStepLoops= " << fpProcessInfo->MAXofPostStepLoops << G4endl;
411 G4Exception("G4ITStepProcessor::GetProcessNumber()",
412 "ITStepProcessor0004", FatalException,
413 "The array size is smaller than the actual No of processes.");
414 }
415
416 if((fpProcessInfo->fpAtRestDoItVector == nullptr) &&
417 (fpProcessInfo->fpAlongStepDoItVector == nullptr) &&
418 (fpProcessInfo->fpPostStepDoItVector == nullptr))
419 {
420 G4ExceptionDescription exceptionDescription;
421 exceptionDescription << "No DoIt process found ";
422 G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessor0005",
423 FatalErrorInArgument,exceptionDescription);
424 return;
425 }
426
427 if((fpProcessInfo->fpAlongStepGetPhysIntVector != nullptr)
428 && fpProcessInfo->MAXofAlongStepLoops>0)
429 {
430 fpProcessInfo->fpTransportation = dynamic_cast<G4ITTransportation*>
431 ((*fpProcessInfo->fpAlongStepGetPhysIntVector)
432 [G4int(fpProcessInfo->MAXofAlongStepLoops-1)]);
433
434 if(fpProcessInfo->fpTransportation == nullptr)
435 {
436 G4ExceptionDescription exceptionDescription;
437 exceptionDescription << "No transportation process found ";
438 G4Exception("G4ITStepProcessor::SetupGeneralProcessInfo",
439 "ITStepProcessor0006",
440 FatalErrorInArgument,exceptionDescription);
441 }
442 }
443 fProcessGeneralInfoMap[particle] = fpProcessInfo;
444 // fpProcessInfo = 0;
445}
@ typeGPIL
@ typeDoIt
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const

Referenced by GetProcessInfo().

◆ SetupMembers()

void G4ITStepProcessor::SetupMembers ( )
protected

Definition at line 510 of file G4ITStepProcessor.cc.

511{
512 fpSecondary = fpStep->GetfSecondary();
513 fpPreStepPoint = fpStep->GetPreStepPoint();
514 fpPostStepPoint = fpStep->GetPostStepPoint();
515
516 fpState = (G4ITStepProcessorState*) fpITrack->GetTrackingInfo()
518
521}
G4TrackVector * GetfSecondary()
G4ITStepProcessorState_Lock * GetStepProcessorState()

Referenced by DoStepping(), and InitDefineStep().

◆ Stepping()

void G4ITStepProcessor::Stepping ( G4Track * track,
const double & timeStep )

Definition at line 263 of file G4ITStepProcessor2.cc.

264{
265
266#ifdef DEBUG_MEM
267 MemStat mem_first, mem_second, mem_diff;
268#endif
269
270#ifdef DEBUG_MEM
271 mem_first = MemoryUsage();
272#endif
273
275
276#ifdef DEBUG_MEM
277 MemStat mem_intermediaire = MemoryUsage();
278 mem_diff = mem_intermediaire-mem_first;
279 G4cout << "\t\t\t >> || MEM || After CleanProcessor " << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
280#endif
281
282 if(track == nullptr) return; // maybe put an exception here
283 fTimeStep = timeStep;
284 SetTrack(track);
285 DoStepping();
286}

Referenced by DoIt().

Friends And Related Symbol Documentation

◆ G4Scheduler

friend class G4Scheduler
friend

Definition at line 154 of file G4ITStepProcessor.hh.


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