Geant4 11.1.1
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 153 of file G4ITStepProcessor.hh.

Constructor & Destructor Documentation

◆ G4ITStepProcessor() [1/2]

G4ITStepProcessor::G4ITStepProcessor ( )

Definition at line 75 of file G4ITStepProcessor.cc.

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

◆ ~G4ITStepProcessor()

G4ITStepProcessor::~G4ITStepProcessor ( )
virtual

Definition at line 227 of file G4ITStepProcessor.cc.

228{
229 if(fpStep)
230 {
231 fpStep->DeleteSecondaryVector();
232 delete fpStep;
233 }
234
235 if(fpSecondary) delete fpSecondary;
237 //G4ITTransportationManager::DeleteInstance();
238
239 // if(fpUserSteppingAction) delete fpUserSteppingAction;
240}
void DeleteSecondaryVector()

◆ G4ITStepProcessor() [2/2]

G4ITStepProcessor::G4ITStepProcessor ( const G4ITStepProcessor other)
protected

Definition at line 243 of file G4ITStepProcessor.cc.

244{
245 fpVerbose = rhs.fpVerbose;
246 fStoreTrajectory = rhs.fStoreTrajectory;
247
248 // fpUserSteppingAction = 0 ;
249 fpTrackingManager = 0;
250 fpNavigator = 0;
251 fInitialized = false;
252
253 kCarTolerance = rhs.kCarTolerance;
254 fInitialized = false;
255 fPreviousTimeStep = DBL_MAX;
256
259 fpTrackContainer = 0;
260 fILTimeStep = DBL_MAX;
261}

Member Function Documentation

◆ ActiveOnlyITProcess() [1/2]

void G4ITStepProcessor::ActiveOnlyITProcess ( )
protected

Definition at line 288 of file G4ITStepProcessor.cc.

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

Referenced by ActiveOnlyITProcess().

◆ ActiveOnlyITProcess() [2/2]

void G4ITStepProcessor::ActiveOnlyITProcess ( G4ProcessManager processManager)
protected

Definition at line 322 of file G4ITStepProcessor.cc.

323{
324 // Method not used for the time being
325 G4ProcessVector* processVector = processManager->GetProcessList();
326
327 G4VITProcess* itProcess = 0;
328 for(G4int i = 0; i < (G4int)processVector->size(); ++i)
329 {
330 G4VProcess* base_process = (*processVector)[i];
331 itProcess = dynamic_cast<G4VITProcess*>(base_process);
332
333 if(!itProcess)
334 {
335 processManager->SetProcessActivation(base_process, false);
336 }
337 }
338}
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 913 of file G4ITStepProcessor2.cc.

914{
915 G4bool tBelowCutEnergyAndSafety = false;
916 G4int tPtclIdx = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
917 if(tPtclIdx < 0)
918 {
919 return;
920 }
921 G4ProductionCutsTable* tCutsTbl =
923 G4int tCoupleIdx = tCutsTbl->GetCoupleIndex(fpPreStepPoint
924 ->GetMaterialCutsCouple());
925 G4double tProdThreshold =
926 (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
927 if(aSecondary->GetKineticEnergy() < tProdThreshold)
928 {
929 tBelowCutEnergyAndSafety = true;
930 if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
931 {
932 G4double currentRange
934 aSecondary->GetKineticEnergy(),
935 fpPreStepPoint->GetMaterialCutsCouple());
936 tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
937 }
938 }
939
940 if(tBelowCutEnergyAndSafety)
941 {
942 if(!(aSecondary->IsGoodForTracking()))
943 {
944 // Add kinetic energy to the total energy deposit
945 fpStep->AddTotalEnergyDeposit(aSecondary->GetKineticEnergy());
946 aSecondary->SetKineticEnergy(0.0);
947 }
948 }
949}
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 441 of file G4ITStepProcessor.hh.

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

Referenced by ApplyProductionCut(), and InvokePSDIP().

◆ CleanProcessor()

void G4ITStepProcessor::CleanProcessor ( )
inline

Definition at line 457 of file G4ITStepProcessor.hh.

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

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

◆ ClearProcessInfo()

void G4ITStepProcessor::ClearProcessInfo ( )
protected

Definition at line 169 of file G4ITStepProcessor.cc.

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

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

◆ ComputeInteractionLength()

G4double G4ITStepProcessor::ComputeInteractionLength ( double  previousTimeStep)

Definition at line 602 of file G4ITStepProcessor.cc.

603{
604 G4TrackManyList* mainList = fpTrackContainer->GetMainList();
605 G4TrackManyList::iterator it = mainList ->begin();
606 G4TrackManyList::iterator end = mainList ->end();
607
608 SetPreviousStepTime(previousTimeStep);
609
610 fILTimeStep = DBL_MAX;
611
612 for (; it != end; )
613 {
614 G4Track * track = *it;
615
616#ifdef DEBUG
617 G4cout << "*CIL* " << GetIT(track)->GetName()
618 << " ID: " << track->GetTrackID()
619 << " at time : " << track->GetGlobalTime()
620 << G4endl;
621#endif
622
623 ++it;
625
627 }
628
629 return fILTimeStep;
630}
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 64 of file G4ITStepProcessor2.cc.

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

698{
699 SetTrack(track);
701}
void SetTrack(G4Track *)

Referenced by ComputeInteractionLength().

◆ DoDefinePhysicalStepLength()

void G4ITStepProcessor::DoDefinePhysicalStepLength ( )
protected

Definition at line 955 of file G4ITStepProcessor.cc.

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

Referenced by DefinePhysicalStepLength().

◆ DoIt()

void G4ITStepProcessor::DoIt ( double  timeStep)

Definition at line 111 of file G4ITStepProcessor2.cc.

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

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

Referenced by Stepping().

◆ ExtractDoItData()

void G4ITStepProcessor::ExtractDoItData ( )

Definition at line 192 of file G4ITStepProcessor2.cc.

193{
194 if (!fpTrack)
195 {
197 return;
198 }
199
200 G4TrackStatus status = fpTrack->GetTrackStatus();
201
202 switch (status)
203 {
204 case fAlive:
205 case fStopButAlive:
206 case fSuspend:
208 default:
210 break;
211
212 case fStopAndKill:
215// G4TrackList::Pop(fpTrack);
216 fpTrackingManager->EndTracking(fpTrack);
217// fTrackContainer->PushToKill(fpTrack);
218 break;
219
222 if (fpSecondary)
223 {
224 for (std::size_t i = 0; i < fpSecondary->size(); ++i)
225 {
226 delete (*fpSecondary)[i];
227 }
228 fpSecondary->clear();
229 }
230// G4TrackList::Pop(fpTrack);
231 fpTrackingManager->EndTracking(fpTrack);
232// fTrackContainer->PushToKill(fpTrack);
233 break;
234 }
235
237}
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fPostponeToNextEvent
void RemoveReactionSet(G4Track *track)
static G4ITReactionSet * Instance()

Referenced by DoIt().

◆ ExtractILData()

void G4ITStepProcessor::ExtractILData ( )
protected

Definition at line 634 of file G4ITStepProcessor.cc.

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

Referenced by ComputeInteractionLength().

◆ FindTransportationStep()

void G4ITStepProcessor::FindTransportationStep ( )

Definition at line 802 of file G4ITStepProcessor2.cc.

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

Referenced by DoStepping().

◆ ForceReInitialization()

void G4ITStepProcessor::ForceReInitialization ( )

Definition at line 188 of file G4ITStepProcessor.cc.

189{
190 fInitialized = false;
192 Initialize();
193}
virtual void Initialize()

◆ GetAtRestDoItProcTriggered()

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

Definition at line 220 of file G4ITStepProcessor.hh.

221 {
222 return fAtRestDoItProcTriggered;
223 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetAtRestIL()

void G4ITStepProcessor::GetAtRestIL ( )
protected

Definition at line 539 of file G4ITStepProcessor.cc.

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

Referenced by DoDefinePhysicalStepLength().

◆ GetCondition()

G4ForceCondition G4ITStepProcessor::GetCondition ( ) const
inline

Definition at line 285 of file G4ITStepProcessor.hh.

286 {
287 return fCondition;
288 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetCurrentProcess()

const G4VITProcess * G4ITStepProcessor::GetCurrentProcess ( ) const
inline

Definition at line 245 of file G4ITStepProcessor.hh.

246 {
247 return fpCurrentProcess;
248 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetCurrentProcessInfo()

const ProcessGeneralInfo * G4ITStepProcessor::GetCurrentProcessInfo ( ) const
inline

Definition at line 265 of file G4ITStepProcessor.hh.

266 {
267 return fpProcessInfo;
268 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetCurrentVolume()

const G4VPhysicalVolume * G4ITStepProcessor::GetCurrentVolume ( ) const
inline

Definition at line 280 of file G4ITStepProcessor.hh.

281 {
282 return fpCurrentVolume;
283 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetGPILSelection()

G4GPILSelection G4ITStepProcessor::GetGPILSelection ( ) const
inline

Definition at line 225 of file G4ITStepProcessor.hh.

226 {
227 return fGPILSelection;
228 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetILTimeStep()

G4double G4ITStepProcessor::GetILTimeStep ( )
inline

Definition at line 203 of file G4ITStepProcessor.hh.

204 {
205 return fILTimeStep;
206 }

◆ GetInteractionTime()

double G4ITStepProcessor::GetInteractionTime ( )
inline

Definition at line 491 of file G4ITStepProcessor.hh.

492{
493 return fTimeStep;
494}

Referenced by ExtractILData().

◆ GetN2ndariesAlongStepDoIt()

G4int G4ITStepProcessor::GetN2ndariesAlongStepDoIt ( ) const
inline

Definition at line 230 of file G4ITStepProcessor.hh.

231 {
232 return fN2ndariesAlongStepDoIt;
233 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetN2ndariesAtRestDoIt()

G4int G4ITStepProcessor::GetN2ndariesAtRestDoIt ( ) const
inline

Definition at line 235 of file G4ITStepProcessor.hh.

236 {
237 return fN2ndariesAtRestDoIt;
238 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetN2ndariesPostStepDoIt()

G4int G4ITStepProcessor::GetN2ndariesPostStepDoIt ( ) const
inline

Definition at line 240 of file G4ITStepProcessor.hh.

241 {
242 return fN2ndariesPostStepDoIt;
243 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetParticleChange()

const G4VParticleChange * G4ITStepProcessor::GetParticleChange ( ) const
inline

Definition at line 275 of file G4ITStepProcessor.hh.

276 {
277 return fpParticleChange;
278 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetPhysIntLength()

G4double G4ITStepProcessor::GetPhysIntLength ( ) const
inline

Definition at line 250 of file G4ITStepProcessor.hh.

251 {
252 return fPhysIntLength;
253 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetPostStepAtTimeDoItProcTriggered()

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

Definition at line 255 of file G4ITStepProcessor.hh.

256 {
257 return fPostStepAtTimeDoItProcTriggered;
258 }

◆ GetPostStepDoItProcTriggered()

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

Definition at line 260 of file G4ITStepProcessor.hh.

261 {
262 return fPostStepDoItProcTriggered;
263 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetProcessInfo()

void G4ITStepProcessor::GetProcessInfo ( )
protected

Definition at line 486 of file G4ITStepProcessor.cc.

487{
488 G4ParticleDefinition* particle = fpTrack->GetDefinition();
489 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it =
490 fProcessGeneralInfoMap.find(particle);
491
492 if(it == fProcessGeneralInfoMap.end())
493 {
495 fpTrack->GetDefinition()->GetProcessManager());
496 if(fpProcessInfo == 0)
497 {
498 G4ExceptionDescription exceptionDescription("...");
499 G4Exception("G4ITStepProcessor::GetProcessNumber",
500 "ITStepProcessor0008",
502 exceptionDescription);
503 return;
504 }
505 }
506 else
507 {
508 fpProcessInfo = it->second;
509 }
510}
void SetupGeneralProcessInfo(G4ParticleDefinition *, G4ProcessManager *)

Referenced by SetupMembers().

◆ GetProcessorState()

const G4ITStepProcessorState * G4ITStepProcessor::GetProcessorState ( ) const
inline

Definition at line 270 of file G4ITStepProcessor.hh.

271 {
272 return fpState;
273 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetSecondaries()

G4TrackVector * G4ITStepProcessor::GetSecondaries ( ) const
inline

Definition at line 179 of file G4ITStepProcessor.hh.

180 {
181 return fpSecondary;
182 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetStep() [1/2]

G4Step * G4ITStepProcessor::GetStep ( )
inline

Definition at line 166 of file G4ITStepProcessor.hh.

167 {
168 return fpStep;
169 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetStep() [2/2]

const G4Step * G4ITStepProcessor::GetStep ( ) const
inline

Definition at line 170 of file G4ITStepProcessor.hh.

171 {
172 return fpStep;
173 }

◆ GetTrack() [1/2]

G4Track * G4ITStepProcessor::GetTrack ( )
inline

Definition at line 162 of file G4ITStepProcessor.hh.

163 {
164 return fpTrack;
165 }

Referenced by G4VITSteppingVerbose::CopyState().

◆ GetTrack() [2/2]

const G4Track * G4ITStepProcessor::GetTrack ( ) const
inline

Definition at line 434 of file G4ITStepProcessor.hh.

435{
436 return fpTrack;
437}

◆ GetTrackingManager()

G4ITTrackingManager * G4ITStepProcessor::GetTrackingManager ( )
inline

Definition at line 187 of file G4ITStepProcessor.hh.

188 {
189 return fpTrackingManager;
190 }

◆ InitDefineStep()

void G4ITStepProcessor::InitDefineStep ( )
protected

ADDED BACK

ADDED BACK

Definition at line 849 of file G4ITStepProcessor.cc.

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

198{
200 if(fInitialized) return;
201 // ActiveOnlyITProcess();
202
204 ->GetNavigatorForTracking());
205
206 fPhysIntLength = DBL_MAX;
207 kCarTolerance = 0.5
209
210 if(fpVerbose == 0)
211 {
212 G4ITTrackingInteractivity* interactivity = fpTrackingManager->GetInteractivity();
213
214 if(interactivity)
215 {
216 fpVerbose = interactivity->GetSteppingVerbose();
217 fpVerbose->SetStepProcessor(this);
218 }
219 }
220
221 fpTrackContainer = G4ITTrackHolder::Instance();
222
223 fInitialized = true;
224}
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 628 of file G4ITStepProcessor2.cc.

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

Referenced by DoStepping().

◆ InvokeAtRestDoItProcs()

void G4ITStepProcessor::InvokeAtRestDoItProcs ( )
protected

Definition at line 561 of file G4ITStepProcessor2.cc.

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

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

Referenced by DoStepping().

◆ InvokePSDIP()

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

Definition at line 762 of file G4ITStepProcessor2.cc.

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

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

Referenced by DoStepping().

◆ operator=()

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

Definition at line 279 of file G4ITStepProcessor.cc.

280{
281 if(this == &rhs) return *this; // handle self assignment
282 //assignment operator
283 return *this;
284}

◆ PrepareLeadingTracks()

void G4ITStepProcessor::PrepareLeadingTracks ( )

Definition at line 272 of file G4ITStepProcessor.cc.

273{
274 fLeadingTracks.PrepareLeadingTracks();
275}

Referenced by G4Scheduler::Stepping().

◆ PushSecondaries()

void G4ITStepProcessor::PushSecondaries ( )
protected

Definition at line 241 of file G4ITStepProcessor2.cc.

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

Referenced by ExtractDoItData().

◆ ResetLeadingTracks()

void G4ITStepProcessor::ResetLeadingTracks ( )

Definition at line 265 of file G4ITStepProcessor.cc.

266{
267 fLeadingTracks.Reset();
268}

Referenced by G4Scheduler::Stepping().

◆ ResetSecondaries()

void G4ITStepProcessor::ResetSecondaries ( )
protected

Definition at line 529 of file G4ITStepProcessor.cc.

530{
531 // Reset the secondary particles
532 fN2ndariesAtRestDoIt = 0;
533 fN2ndariesAlongStepDoIt = 0;
534 fN2ndariesPostStepDoIt = 0;
535}

Referenced by G4ITStepProcessor(), and SetupMembers().

◆ SetInitialStep()

void G4ITStepProcessor::SetInitialStep ( )
protected

Definition at line 705 of file G4ITStepProcessor.cc.

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

Referenced by InitDefineStep().

◆ SetNavigator()

void G4ITStepProcessor::SetNavigator ( G4ITNavigator *  value)
inlineprotected

Definition at line 450 of file G4ITStepProcessor.hh.

451{
452 fpNavigator = value;
453}

Referenced by Initialize().

◆ SetPreviousStepTime()

void G4ITStepProcessor::SetPreviousStepTime ( G4double  previousTimeStep)
inline

Definition at line 427 of file G4ITStepProcessor.hh.

428{
429 fPreviousTimeStep = previousTimeStep;
430}

Referenced by ComputeInteractionLength().

◆ SetStep()

void G4ITStepProcessor::SetStep ( G4Step val)
inline

Definition at line 174 of file G4ITStepProcessor.hh.

175 {
176 fpStep = val;
177 }

◆ SetTrack()

void G4ITStepProcessor::SetTrack ( G4Track track)
protected

Definition at line 453 of file G4ITStepProcessor.cc.

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

Referenced by DefinePhysicalStepLength(), and Stepping().

◆ SetTrackingManager()

void G4ITStepProcessor::SetTrackingManager ( G4ITTrackingManager trackMan)
inline

Definition at line 183 of file G4ITStepProcessor.hh.

184 {
185 fpTrackingManager = trackMan;
186 }

Referenced by G4Scheduler::Initialize().

◆ SetupGeneralProcessInfo()

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

Definition at line 342 of file G4ITStepProcessor.cc.

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

515{
516 fpSecondary = fpStep->GetfSecondary();
517 fpPreStepPoint = fpStep->GetPreStepPoint();
518 fpPostStepPoint = fpStep->GetPostStepPoint();
519
520 fpState = (G4ITStepProcessorState*) fpITrack->GetTrackingInfo()
522
525}
G4TrackVector * GetfSecondary()
G4ITStepProcessorState_Lock * GetStepProcessorState()

Referenced by DoStepping(), and InitDefineStep().

◆ Stepping()

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

Definition at line 264 of file G4ITStepProcessor2.cc.

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

Referenced by DoIt().

Friends And Related Function Documentation

◆ G4Scheduler

friend class G4Scheduler
friend

Definition at line 155 of file G4ITStepProcessor.hh.


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