Geant4  10.02.p03
G4ITStepProcessor Class Reference

#include <G4ITStepProcessor.hh>

Collaboration diagram for G4ITStepProcessor:

Public Member Functions

 G4ITStepProcessor ()
 
virtual ~G4ITStepProcessor ()
 
void SetPreviousStepTime (G4double)
 
G4Track * GetTrack ()
 
G4Step * GetStep ()
 
const G4Step * GetStep () const
 
void SetStep (G4Step *val)
 
G4TrackVector * GetSecondaries () 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 G4Track * GetTrack () const
 
void CleanProcessor ()
 
size_t GetAtRestDoItProcTriggered () const
 
G4GPILSelection GetGPILSelection () const
 
G4int GetN2ndariesAlongStepDoIt () const
 
G4int GetN2ndariesAtRestDoIt () const
 
G4int GetN2ndariesPostStepDoIt () const
 
const G4VITProcessGetCurrentProcess () const
 
G4double GetPhysIntLength () const
 
size_t GetPostStepAtTimeDoItProcTriggered () const
 
size_t GetPostStepDoItProcTriggered () const
 
const ProcessGeneralInfoGetCurrentProcessInfo () const
 
const G4ITStepProcessorStateGetProcessorState () const
 
const G4VParticleChange * GetParticleChange () 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 (size_t)
 
void InvokeTransportationProc ()
 
void SetNavigator (G4ITNavigator *value)
 
G4double CalculateSafety ()
 
void ApplyProductionCut (G4Track *)
 
 G4ITStepProcessor (const G4ITStepProcessor &other)
 
G4ITStepProcessoroperator= (const G4ITStepProcessor &other)
 

Private Attributes

G4bool fInitialized
 
G4ITTrackingManagerfpTrackingManager
 
G4double kCarTolerance
 
G4ITNavigator * fpNavigator
 
G4int fStoreTrajectory
 
G4VITSteppingVerbosefpVerbose
 
G4ITTrackHolderfpTrackContainer
 
G4ITLeadingTracks fLeadingTracks
 
G4double fTimeStep
 
G4double fILTimeStep
 
G4double fPreviousTimeStep
 
G4TrackVector * fpSecondary
 
G4VParticleChange * fpParticleChange
 
G4VITProcessfpCurrentProcess
 
G4int fN2ndariesAtRestDoIt
 
G4int fN2ndariesAlongStepDoIt
 
G4int fN2ndariesPostStepDoIt
 
size_t fAtRestDoItProcTriggered
 
size_t fPostStepDoItProcTriggered
 
size_t fPostStepAtTimeDoItProcTriggered
 
G4ForceCondition fCondition
 
G4GPILSelection fGPILSelection
 
G4double fPhysIntLength
 
G4VPhysicalVolumefpCurrentVolume
 
std::map< const G4ParticleDefinition *, ProcessGeneralInfo * > fProcessGeneralInfoMap
 
ProcessGeneralInfofpProcessInfo
 
G4ITTransportationfpTransportation
 
G4Track * fpTrack
 
G4ITfpITrack
 
G4TrackingInformationfpTrackingInfo
 
G4ITStepProcessorStatefpState
 
G4Step * fpStep
 
G4StepPoint * fpPreStepPoint
 
G4StepPoint * fpPostStepPoint
 

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 154 of file G4ITStepProcessor.hh.

Constructor & Destructor Documentation

◆ G4ITStepProcessor() [1/2]

G4ITStepProcessor::G4ITStepProcessor ( )

Definition at line 77 of file G4ITStepProcessor.cc.

78 {
79  fpVerbose = 0;
80  // fpUserSteppingAction = 0 ;
81  fStoreTrajectory = 0;
83  fpNavigator = 0;
84  kCarTolerance = -1.;
85  fInitialized = false;
88  fpTrackContainer = 0;
89 
92 }
G4ITTrackingManager * fpTrackingManager
G4ITTrackHolder * fpTrackContainer
G4ITNavigator * fpNavigator
G4VITSteppingVerbose * fpVerbose
#define DBL_MAX
Definition: templates.hh:83

◆ ~G4ITStepProcessor()

G4ITStepProcessor::~G4ITStepProcessor ( )
virtual

Definition at line 229 of file G4ITStepProcessor.cc.

230 {
231  if(fpStep)
232  {
233  fpStep->DeleteSecondaryVector();
234  delete fpStep;
235  }
236 
237  if(fpSecondary) delete fpSecondary;
239  //G4ITTransportationManager::DeleteInstance();
240 
241  // if(fpUserSteppingAction) delete fpUserSteppingAction;
242 }
G4TrackVector * fpSecondary

◆ G4ITStepProcessor() [2/2]

G4ITStepProcessor::G4ITStepProcessor ( const G4ITStepProcessor other)
protected

Definition at line 245 of file G4ITStepProcessor.cc.

246 {
247  fpVerbose = rhs.fpVerbose;
248  fStoreTrajectory = rhs.fStoreTrajectory;
249 
250  // fpUserSteppingAction = 0 ;
251  fpTrackingManager = 0;
252  fpNavigator = 0;
253  fInitialized = false;
254 
255  kCarTolerance = rhs.kCarTolerance;
256  fInitialized = false;
258 
259  CleanProcessor();
261  fpTrackContainer = 0;
263 }
G4ITTrackingManager * fpTrackingManager
G4ITTrackHolder * fpTrackContainer
G4ITNavigator * fpNavigator
G4VITSteppingVerbose * fpVerbose
#define DBL_MAX
Definition: templates.hh:83

Member Function Documentation

◆ ActiveOnlyITProcess() [1/2]

void G4ITStepProcessor::ActiveOnlyITProcess ( )
protected

Definition at line 290 of file G4ITStepProcessor.cc.

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

◆ ActiveOnlyITProcess() [2/2]

void G4ITStepProcessor::ActiveOnlyITProcess ( G4ProcessManager processManager)
protected

Definition at line 324 of file G4ITStepProcessor.cc.

325 {
326  // Method not used for the time being
327  G4ProcessVector* processVector = processManager->GetProcessList();
328 
329  G4VITProcess* itProcess = 0;
330  for(int i = 0; i < processVector->size(); i++)
331  {
332  G4VProcess* base_process = (*processVector)[i];
333  itProcess = dynamic_cast<G4VITProcess*>(base_process);
334 
335  if(!itProcess)
336  {
337  processManager->SetProcessActivation(base_process, false);
338  }
339  }
340 }
G4ProcessVector * GetProcessList() const
G4VProcess * SetProcessActivation(G4VProcess *aProcess, G4bool fActive)
G4int size() const
Here is the call graph for this function:

◆ ApplyProductionCut()

void G4ITStepProcessor::ApplyProductionCut ( G4Track *  aSecondary)
protected

Definition at line 914 of file G4ITStepProcessor2.cc.

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

◆ CalculateSafety()

G4double G4ITStepProcessor::CalculateSafety ( )
inlineprotected

Definition at line 442 of file G4ITStepProcessor.hh.

443 {
445  - fpPostStepPoint->GetPosition()).mag(),
446  kCarTolerance);
447 }
G4ThreeVector fEndpointSafOrigin
G4ITStepProcessorState * fpState
G4StepPoint * fpPostStepPoint

◆ CleanProcessor()

void G4ITStepProcessor::CleanProcessor ( )
inline

Definition at line 458 of file G4ITStepProcessor.hh.

459 {
460  fTimeStep = DBL_MAX;
462 
463  fpState = 0;
464  fpTrack = 0;
465  fpTrackingInfo = 0;
466  fpITrack = 0;
467  fpStep = 0;
468  fpPreStepPoint = 0;
469  fpPostStepPoint = 0;
470 
471  fpParticleChange = 0;
472 
473  fpCurrentVolume = 0;
474  // fpSensitive = 0;
475 
476  fpSecondary = 0;
477 
478  fpTransportation = 0;
479 
480  fpCurrentProcess= 0;
481  fpProcessInfo = 0;
482 
486  fGPILSelection = NotCandidateForSelection;
487  fCondition = NotForced;
488 }
G4TrackVector * fpSecondary
G4GPILSelection fGPILSelection
G4VITProcess * fpCurrentProcess
G4ForceCondition fCondition
G4TrackingInformation * fpTrackingInfo
G4StepPoint * fpPreStepPoint
G4VPhysicalVolume * fpCurrentVolume
G4ITStepProcessorState * fpState
G4StepPoint * fpPostStepPoint
#define INT_MAX
Definition: templates.hh:111
G4VParticleChange * fpParticleChange
G4ITTransportation * fpTransportation
ProcessGeneralInfo * fpProcessInfo
#define DBL_MAX
Definition: templates.hh:83
size_t fPostStepAtTimeDoItProcTriggered

◆ ClearProcessInfo()

void G4ITStepProcessor::ClearProcessInfo ( )
protected

Definition at line 171 of file G4ITStepProcessor.cc.

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

◆ ComputeInteractionLength()

G4double G4ITStepProcessor::ComputeInteractionLength ( double  previousTimeStep)

Definition at line 604 of file G4ITStepProcessor.cc.

605 {
607  G4TrackManyList::iterator it = mainList ->begin();
608  G4TrackManyList::iterator end = mainList ->end();
609 
610  SetPreviousStepTime(previousTimeStep);
611 
613 
614  for (; it != end; )
615  {
616  G4Track * track = *it;
617 
618 #ifdef DEBUG
619  G4cout << "*CIL* " << GetIT(track)->GetName()
620  << " ID: " << track->GetTrackID()
621  << " at time : " << track->GetGlobalTime()
622  << G4endl;
623 #endif
624 
625  ++it;
627 
628  ExtractILData();
629  }
630 
631  return fILTimeStep;
632 }
virtual const G4String & GetName() const =0
void SetPreviousStepTime(G4double)
G4ITTrackHolder * fpTrackContainer
void DefinePhysicalStepLength(G4Track *)
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:49
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4TrackList * GetMainList(Key)
#define DBL_MAX
Definition: templates.hh:83
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DealWithSecondaries()

void G4ITStepProcessor::DealWithSecondaries ( G4int counter)
protected

Definition at line 65 of file G4ITStepProcessor2.cc.

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

◆ DefinePhysicalStepLength()

void G4ITStepProcessor::DefinePhysicalStepLength ( G4Track *  track)

Definition at line 699 of file G4ITStepProcessor.cc.

700 {
701  SetTrack(track);
703 }
void SetTrack(G4Track *)

◆ DoDefinePhysicalStepLength()

void G4ITStepProcessor::DoDefinePhysicalStepLength ( )
protected

Definition at line 956 of file G4ITStepProcessor.cc.

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

◆ DoIt()

void G4ITStepProcessor::DoIt ( double  timeStep)

Definition at line 112 of file G4ITStepProcessor2.cc.

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

◆ DoStepping()

void G4ITStepProcessor::DoStepping ( )
protected

Definition at line 294 of file G4ITStepProcessor2.cc.

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

◆ ExtractDoItData()

void G4ITStepProcessor::ExtractDoItData ( )

Definition at line 193 of file G4ITStepProcessor2.cc.

194 {
195  if (!fpTrack)
196  {
197  CleanProcessor();
198  return;
199  }
200 
201  G4TrackStatus status = fpTrack->GetTrackStatus();
202 
203  switch (status)
204  {
205  case fAlive:
206  case fStopButAlive:
207  case fSuspend:
208  case fPostponeToNextEvent:
209  default:
210  PushSecondaries();
211  break;
212 
213  case fStopAndKill:
215  PushSecondaries();
216 // G4TrackList::Pop(fpTrack);
218 // fTrackContainer->PushToKill(fpTrack);
219  break;
220 
221  case fKillTrackAndSecondaries:
223  if (fpSecondary)
224  {
225  for (size_t i = 0; i < fpSecondary->size(); ++i)
226  {
227  delete (*fpSecondary)[i];
228  }
229  fpSecondary->clear();
230  }
231 // G4TrackList::Pop(fpTrack);
233 // fTrackContainer->PushToKill(fpTrack);
234  break;
235  }
236 
237  CleanProcessor();
238 }
G4TrackVector * fpSecondary
G4ITTrackingManager * fpTrackingManager
void RemoveReactionSet(G4Track *track)
static G4ITReactionSet * Instance()
void EndTracking(G4Track *)
Here is the call graph for this function:

◆ ExtractILData()

void G4ITStepProcessor::ExtractILData ( )
protected

Definition at line 636 of file G4ITStepProcessor.cc.

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

◆ FindTransportationStep()

void G4ITStepProcessor::FindTransportationStep ( )

Definition at line 803 of file G4ITStepProcessor2.cc.

804 {
805  double physicalStep(0.);
806 
808  // dynamic_cast<G4ITTransportation*>((fpProcessInfo->fpAlongStepGetPhysIntVector)[MAXofAlongStepLoops-1]);
809 
810  if(!fpTrack)
811  {
812  G4ExceptionDescription exceptionDescription;
813  exceptionDescription << "No G4ITStepProcessor::fpTrack found";
814  G4Exception("G4ITStepProcessor::FindTransportationStep",
815  "ITStepProcessor0013",
817  exceptionDescription);
818  return;
819 
820  }
821  if(!fpITrack)
822  {
823  G4ExceptionDescription exceptionDescription;
824  exceptionDescription << "No G4ITStepProcessor::fITrack";
825  G4Exception("G4ITStepProcessor::FindTransportationStep",
826  "ITStepProcessor0014",
828  exceptionDescription);
829  return;
830  }
831  if(!(fpITrack->GetTrack()))
832  {
833  G4ExceptionDescription exceptionDescription;
834  exceptionDescription << "No G4ITStepProcessor::fITrack->GetTrack()";
835  G4Exception("G4ITStepProcessor::FindTransportationStep",
836  "ITStepProcessor0015",
838  exceptionDescription);
839  return;
840  }
841 
842  if(fpTransportation)
843  {
845  ->GetProcessID()));
847 
848 // fpTransportation->SetProcessState(0);
850  }
851 
852  if(physicalStep >= DBL_MAX)
853  {
854 // G4cout << "---- 2) Setting stop and kill for " << GetIT(fpTrack)->GetName() << G4endl;
855  fpTrack -> SetTrackStatus(fStopAndKill);
856  return;
857  }
858 
859  fpState->fPhysicalStep = physicalStep;
860 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4ITTransportation * fpTransportation
size_t GetProcessID() const
G4TrackingInformation * fpTrackingInfo
virtual void ComputeStep(const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
G4ITStepProcessorState * fpState
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetProcessState(G4shared_ptr< G4ProcessState_Lock > aProcInfo)
G4shared_ptr< G4ProcessState_Lock > GetProcessState(size_t index)
G4ITTransportation * fpTransportation
G4Track * GetTrack()
Definition: G4IT.hh:215
ProcessGeneralInfo * fpProcessInfo
#define DBL_MAX
Definition: templates.hh:83
void ResetProcessState()
Here is the call graph for this function:

◆ ForceReInitialization()

void G4ITStepProcessor::ForceReInitialization ( )

Definition at line 190 of file G4ITStepProcessor.cc.

191 {
192  fInitialized = false;
194  Initialize();
195 }
virtual void Initialize()
Here is the call graph for this function:

◆ GetAtRestDoItProcTriggered()

size_t G4ITStepProcessor::GetAtRestDoItProcTriggered ( ) const
inline

Definition at line 221 of file G4ITStepProcessor.hh.

222  {
224  }
Here is the caller graph for this function:

◆ GetAtRestIL()

void G4ITStepProcessor::GetAtRestIL ( )
protected

Definition at line 541 of file G4ITStepProcessor.cc.

542 {
543  // Select the rest process which has the shortest time before
544  // it is invoked. In rest processes, GPIL()
545  // returns the time before a process occurs.
546  G4double lifeTime(DBL_MAX), shortestLifeTime (DBL_MAX);
547 
549  shortestLifeTime = DBL_MAX;
550 
551  unsigned int NofInactiveProc=0;
552 
553  for( size_t ri=0; ri < fpProcessInfo->MAXofAtRestLoops; ri++ )
554  {
556  if (fpCurrentProcess== 0)
557  {
558  (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
559  NofInactiveProc++;
560  continue;
561  } // NULL means the process is inactivated by a user on fly.
562 
563  fCondition=NotForced;
566 
567  lifeTime = fpCurrentProcess->AtRestGPIL( *fpTrack, &fCondition );
569 
570  if(fCondition==Forced)
571  {
572  (fpState->fSelectedAtRestDoItVector)[ri] = Forced;
573  }
574  else
575  {
576  (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
577  if(lifeTime < shortestLifeTime )
578  {
579  shortestLifeTime = lifeTime;
581  }
582  }
583  }
584 
586 
587 // G4cout << " --> Selected at rest process : "
588 // << (*fpProcessInfo->fpAtRestGetPhysIntVector)[fAtRestDoItProcTriggered]
589 // ->GetProcessName()
590 // << G4endl;
591 
592  fTimeStep = shortestLifeTime;
593 
594  // at least one process is necessary to destroy the particle
595  // exit with warning
596  if(NofInactiveProc==fpProcessInfo->MAXofAtRestLoops)
597  {
598  G4cerr << "ERROR - G4ITStepProcessor::InvokeAtRestDoItProcs()" << G4endl
599  << " No AtRestDoIt process is active!" << G4endl;
600  }
601 }
G4ProcessVector * fpAtRestGetPhysIntVector
size_t GetProcessID() const
G4VITProcess * fpCurrentProcess
int G4int
Definition: G4Types.hh:78
G4SelectedAtRestDoItVector fSelectedAtRestDoItVector
G4ForceCondition fCondition
G4double AtRestGPIL(const G4Track &track, G4ForceCondition *condition)
Definition: G4VProcess.hh:490
G4TrackingInformation * fpTrackingInfo
G4ITStepProcessorState * fpState
void SetProcessState(G4shared_ptr< G4ProcessState_Lock > aProcInfo)
G4shared_ptr< G4ProcessState_Lock > GetProcessState(size_t index)
ProcessGeneralInfo * fpProcessInfo
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void ResetProcessState()
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:

◆ GetCondition()

G4ForceCondition G4ITStepProcessor::GetCondition ( ) const
inline

Definition at line 286 of file G4ITStepProcessor.hh.

287  {
288  return fCondition;
289  }
G4ForceCondition fCondition
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCurrentProcess()

const G4VITProcess* G4ITStepProcessor::GetCurrentProcess ( ) const
inline

Definition at line 246 of file G4ITStepProcessor.hh.

247  {
248  return fpCurrentProcess;
249  }
G4VITProcess * fpCurrentProcess
Here is the caller graph for this function:

◆ GetCurrentProcessInfo()

const ProcessGeneralInfo* G4ITStepProcessor::GetCurrentProcessInfo ( ) const
inline

Definition at line 266 of file G4ITStepProcessor.hh.

267  {
268  return fpProcessInfo;
269  }
ProcessGeneralInfo * fpProcessInfo
Here is the caller graph for this function:

◆ GetCurrentVolume()

const G4VPhysicalVolume* G4ITStepProcessor::GetCurrentVolume ( ) const
inline

Definition at line 281 of file G4ITStepProcessor.hh.

282  {
283  return fpCurrentVolume;
284  }
G4VPhysicalVolume * fpCurrentVolume
Here is the caller graph for this function:

◆ GetGPILSelection()

G4GPILSelection G4ITStepProcessor::GetGPILSelection ( ) const
inline

Definition at line 226 of file G4ITStepProcessor.hh.

227  {
228  return fGPILSelection;
229  }
G4GPILSelection fGPILSelection
Here is the caller graph for this function:

◆ GetILTimeStep()

G4double G4ITStepProcessor::GetILTimeStep ( )
inline

Definition at line 204 of file G4ITStepProcessor.hh.

205  {
206  return fILTimeStep;
207  }
Here is the call graph for this function:

◆ GetInteractionTime()

double G4ITStepProcessor::GetInteractionTime ( )
inline

Definition at line 492 of file G4ITStepProcessor.hh.

493 {
494  return fTimeStep;
495 }

◆ GetN2ndariesAlongStepDoIt()

G4int G4ITStepProcessor::GetN2ndariesAlongStepDoIt ( ) const
inline

Definition at line 231 of file G4ITStepProcessor.hh.

232  {
234  }
Here is the caller graph for this function:

◆ GetN2ndariesAtRestDoIt()

G4int G4ITStepProcessor::GetN2ndariesAtRestDoIt ( ) const
inline

Definition at line 236 of file G4ITStepProcessor.hh.

237  {
238  return fN2ndariesAtRestDoIt;
239  }
Here is the caller graph for this function:

◆ GetN2ndariesPostStepDoIt()

G4int G4ITStepProcessor::GetN2ndariesPostStepDoIt ( ) const
inline

Definition at line 241 of file G4ITStepProcessor.hh.

242  {
243  return fN2ndariesPostStepDoIt;
244  }
Here is the caller graph for this function:

◆ GetParticleChange()

const G4VParticleChange* G4ITStepProcessor::GetParticleChange ( ) const
inline

Definition at line 276 of file G4ITStepProcessor.hh.

277  {
278  return fpParticleChange;
279  }
G4VParticleChange * fpParticleChange
Here is the caller graph for this function:

◆ GetPhysIntLength()

G4double G4ITStepProcessor::GetPhysIntLength ( ) const
inline

Definition at line 251 of file G4ITStepProcessor.hh.

252  {
253  return fPhysIntLength;
254  }
Here is the caller graph for this function:

◆ GetPostStepAtTimeDoItProcTriggered()

size_t G4ITStepProcessor::GetPostStepAtTimeDoItProcTriggered ( ) const
inline

Definition at line 256 of file G4ITStepProcessor.hh.

257  {
259  }
size_t fPostStepAtTimeDoItProcTriggered

◆ GetPostStepDoItProcTriggered()

size_t G4ITStepProcessor::GetPostStepDoItProcTriggered ( ) const
inline

Definition at line 261 of file G4ITStepProcessor.hh.

262  {
264  }
Here is the caller graph for this function:

◆ GetProcessInfo()

void G4ITStepProcessor::GetProcessInfo ( )
protected

Definition at line 488 of file G4ITStepProcessor.cc.

489 {
490  G4ParticleDefinition* particle = fpTrack->GetDefinition();
491  std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it =
492  fProcessGeneralInfoMap.find(particle);
493 
494  if(it == fProcessGeneralInfoMap.end())
495  {
496  SetupGeneralProcessInfo(particle,
497  fpTrack->GetDefinition()->GetProcessManager());
498  if(fpProcessInfo == 0)
499  {
500  G4ExceptionDescription exceptionDescription("...");
501  G4Exception("G4ITStepProcessor::GetProcessNumber",
502  "ITStepProcessor0008",
504  exceptionDescription);
505  return;
506  }
507  }
508  else
509  {
510  fpProcessInfo = it->second;
511  }
512 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
std::map< const G4ParticleDefinition *, ProcessGeneralInfo * > fProcessGeneralInfoMap
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
ProcessGeneralInfo * fpProcessInfo
void SetupGeneralProcessInfo(G4ParticleDefinition *, G4ProcessManager *)
Here is the call graph for this function:

◆ GetProcessorState()

const G4ITStepProcessorState* G4ITStepProcessor::GetProcessorState ( ) const
inline

Definition at line 271 of file G4ITStepProcessor.hh.

272  {
273  return fpState;
274  }
G4ITStepProcessorState * fpState
Here is the caller graph for this function:

◆ GetSecondaries()

G4TrackVector* G4ITStepProcessor::GetSecondaries ( ) const
inline

Definition at line 180 of file G4ITStepProcessor.hh.

181  {
182  return fpSecondary;
183  }
G4TrackVector * fpSecondary
Here is the caller graph for this function:

◆ GetStep() [1/2]

G4Step* G4ITStepProcessor::GetStep ( )
inline

Definition at line 167 of file G4ITStepProcessor.hh.

168  {
169  return fpStep;
170  }
Here is the caller graph for this function:

◆ GetStep() [2/2]

const G4Step* G4ITStepProcessor::GetStep ( ) const
inline

Definition at line 171 of file G4ITStepProcessor.hh.

172  {
173  return fpStep;
174  }

◆ GetTrack() [1/2]

G4Track* G4ITStepProcessor::GetTrack ( )
inline

Definition at line 163 of file G4ITStepProcessor.hh.

164  {
165  return fpTrack;
166  }
Here is the caller graph for this function:

◆ GetTrack() [2/2]

const G4Track * G4ITStepProcessor::GetTrack ( ) const
inline

Definition at line 435 of file G4ITStepProcessor.hh.

436 {
437  return fpTrack;
438 }

◆ GetTrackingManager()

G4ITTrackingManager* G4ITStepProcessor::GetTrackingManager ( )
inline

Definition at line 188 of file G4ITStepProcessor.hh.

189  {
190  return fpTrackingManager;
191  }
G4ITTrackingManager * fpTrackingManager
Here is the call graph for this function:

◆ InitDefineStep()

void G4ITStepProcessor::InitDefineStep ( )
protected

ADDED BACK

ADDED BACK

Definition at line 850 of file G4ITStepProcessor.cc.

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

◆ Initialize()

void G4ITStepProcessor::Initialize ( )
virtual

Definition at line 199 of file G4ITStepProcessor.cc.

200 {
201  CleanProcessor();
202  if(fInitialized) return;
203  // ActiveOnlyITProcess();
204 
206  ->GetNavigatorForTracking());
207 
209  kCarTolerance = 0.5
211 
212  if(fpVerbose == 0)
213  {
215 
216  if(interactivity)
217  {
218  fpVerbose = interactivity->GetSteppingVerbose();
220  }
221  }
222 
224 
225  fInitialized = true;
226 }
G4ITTrackingInteractivity * GetInteractivity()
G4ITTrackingManager * fpTrackingManager
void SetNavigator(G4ITNavigator *value)
G4double GetSurfaceTolerance() const
G4ITTrackHolder * fpTrackContainer
static G4ITTransportationManager * GetTransportationManager()
G4VITSteppingVerbose * GetSteppingVerbose()
static G4ITTrackHolder * Instance()
G4VITSteppingVerbose * fpVerbose
void SetStepProcessor(const G4ITStepProcessor *stepProcessor)
#define DBL_MAX
Definition: templates.hh:83
static G4GeometryTolerance * GetInstance()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InvokeAlongStepDoItProcs()

void G4ITStepProcessor::InvokeAlongStepDoItProcs ( )
protected

Definition at line 629 of file G4ITStepProcessor2.cc.

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

◆ InvokeAtRestDoItProcs()

void G4ITStepProcessor::InvokeAtRestDoItProcs ( )
protected

Definition at line 562 of file G4ITStepProcessor2.cc.

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

◆ InvokePostStepDoItProcs()

void G4ITStepProcessor::InvokePostStepDoItProcs ( )
protected

Definition at line 713 of file G4ITStepProcessor2.cc.

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

◆ InvokePSDIP()

void G4ITStepProcessor::InvokePSDIP ( size_t  np)
protected

Definition at line 763 of file G4ITStepProcessor2.cc.

764 {
766 
768  ->GetProcessID()));
770 // fpCurrentProcess->SetProcessState(0);
772 
773  // Update PostStepPoint of Step according to ParticleChange
774  fpParticleChange->UpdateStepForPostStep(fpStep);
775 
776 #ifdef G4VERBOSE
777  // !!!!! Verbose
779 #endif
780 
781  // Update G4Track according to ParticleChange after each PostStepDoIt
782  fpStep->UpdateTrack();
783 
784  // Update safety after each invocation of PostStepDoIts
785  fpStep->GetPostStepPoint()->SetSafety(CalculateSafety());
786 
787  // Now Store the secondaries from ParticleChange to SecondaryList
789 
790  // Set the track status according to what the process defined
791  fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
792 
793  // clear ParticleChange
794  fpParticleChange->Clear();
795 }
virtual void PostStepDoItOneByOne()=0
G4ProcessVector * fpPostStepDoItVector
size_t GetProcessID() const
G4VITProcess * fpCurrentProcess
G4TrackingInformation * fpTrackingInfo
void DealWithSecondaries(G4int &)
G4VParticleChange * fpParticleChange
void SetProcessState(G4shared_ptr< G4ProcessState_Lock > aProcInfo)
G4shared_ptr< G4ProcessState_Lock > GetProcessState(size_t index)
ProcessGeneralInfo * fpProcessInfo
G4VITSteppingVerbose * fpVerbose
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
void ResetProcessState()

◆ InvokeTransportationProc()

void G4ITStepProcessor::InvokeTransportationProc ( )
protected

Definition at line 864 of file G4ITStepProcessor2.cc.

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

◆ operator=()

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

Definition at line 281 of file G4ITStepProcessor.cc.

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

◆ PrepareLeadingTracks()

void G4ITStepProcessor::PrepareLeadingTracks ( )

Definition at line 274 of file G4ITStepProcessor.cc.

275 {
277 }
G4ITLeadingTracks fLeadingTracks
Here is the caller graph for this function:

◆ PushSecondaries()

void G4ITStepProcessor::PushSecondaries ( )
protected

Definition at line 242 of file G4ITStepProcessor2.cc.

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

◆ ResetLeadingTracks()

void G4ITStepProcessor::ResetLeadingTracks ( )

Definition at line 267 of file G4ITStepProcessor.cc.

268 {
270 }
G4ITLeadingTracks fLeadingTracks
Here is the caller graph for this function:

◆ ResetSecondaries()

void G4ITStepProcessor::ResetSecondaries ( )
protected

Definition at line 531 of file G4ITStepProcessor.cc.

532 {
533  // Reset the secondary particles
537 }

◆ SetInitialStep()

void G4ITStepProcessor::SetInitialStep ( )
protected

Definition at line 707 of file G4ITStepProcessor.cc.

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

◆ SetNavigator()

void G4ITStepProcessor::SetNavigator ( G4ITNavigator *  value)
inlineprotected

Definition at line 451 of file G4ITStepProcessor.hh.

452 {
453  fpNavigator = value;
454 }
G4ITNavigator * fpNavigator

◆ SetPreviousStepTime()

void G4ITStepProcessor::SetPreviousStepTime ( G4double  previousTimeStep)
inline

Definition at line 428 of file G4ITStepProcessor.hh.

429 {
430  fPreviousTimeStep = previousTimeStep;
431 }

◆ SetStep()

void G4ITStepProcessor::SetStep ( G4Step *  val)
inline

Definition at line 175 of file G4ITStepProcessor.hh.

176  {
177  fpStep = val;
178  }

◆ SetTrack()

void G4ITStepProcessor::SetTrack ( G4Track *  track)
protected

Definition at line 455 of file G4ITStepProcessor.cc.

456 {
457  fpTrack = track;
458  if(fpTrack)
459  {
461  fpStep = const_cast<G4Step*>(fpTrack->GetStep());
462 
463  if(fpITrack)
464  {
466  }
467  else
468  {
469  fpTrackingInfo = 0;
470  G4cerr << "Track ID : " << fpTrack->GetTrackID() << G4endl;
471 
472  G4ExceptionDescription errMsg;
473  errMsg << "No IT pointer was attached to the track you try to process.";
474  G4Exception("G4ITStepProcessor::SetTrack",
475  "ITStepProcessor0007",
477  errMsg);
478  }
479  }
480  else
481  {
482  fpITrack = 0;
483  fpStep = 0;
484  }
485 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:49
G4TrackingInformation * fpTrackingInfo
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:140
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:

◆ SetTrackingManager()

void G4ITStepProcessor::SetTrackingManager ( G4ITTrackingManager trackMan)
inline

Definition at line 184 of file G4ITStepProcessor.hh.

185  {
186  fpTrackingManager = trackMan;
187  }
G4ITTrackingManager * fpTrackingManager
Here is the caller graph for this function:

◆ SetupGeneralProcessInfo()

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

Definition at line 344 of file G4ITStepProcessor.cc.

346 {
347 
348 #ifdef debug
349  G4cout<<"G4ITStepProcessor::GetProcessNumber: is called track"<<G4endl;
350 #endif
351  if(!pm)
352  {
353  G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl<< " ProcessManager is NULL for particle = "
354  << particle->GetParticleName() << ", PDG_code = "
355  << particle->GetPDGEncoding() << G4endl;
356  G4Exception("G4SteppingManager::GetProcessNumber()", "ITStepProcessor0002",
357  FatalException, "Process Manager is not found.");
358  return;
359  }
360 
361  std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it =
362  fProcessGeneralInfoMap.find(particle);
363  if(it != fProcessGeneralInfoMap.end())
364  {
365  G4Exception("G4SteppingManager::SetupGeneralProcessInfo()",
366  "ITStepProcessor0003",
367  FatalException, "Process info already registered.");
368  return;
369  }
370 
371  // here used as temporary
373 
374  // AtRestDoits
379 #ifdef debug
380  G4cout << "G4ITStepProcessor::GetProcessNumber: #ofAtRest="
382 #endif
383 
384  // AlongStepDoits
391 #ifdef debug
392  G4cout << "G4ITStepProcessor::GetProcessNumber:#ofAlongStp="
394 #endif
395 
396  // PostStepDoits
402 #ifdef debug
403  G4cout << "G4ITStepProcessor::GetProcessNumber: #ofPostStep="
405 #endif
406 
407  if (SizeOfSelectedDoItVector<fpProcessInfo->MAXofAtRestLoops ||
408  SizeOfSelectedDoItVector<fpProcessInfo->MAXofAlongStepLoops ||
409  SizeOfSelectedDoItVector<fpProcessInfo->MAXofPostStepLoops )
410  {
411  G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl
412  << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
413  << " ; is smaller then one of MAXofAtRestLoops= "
414  << fpProcessInfo->MAXofAtRestLoops << G4endl
415  << " or MAXofAlongStepLoops= " << fpProcessInfo->MAXofAlongStepLoops
416  << " or MAXofPostStepLoops= " << fpProcessInfo->MAXofPostStepLoops << G4endl;
417  G4Exception("G4ITStepProcessor::GetProcessNumber()",
418  "ITStepProcessor0004", FatalException,
419  "The array size is smaller than the actual No of processes.");
420  }
421 
425  {
426  G4ExceptionDescription exceptionDescription;
427  exceptionDescription << "No DoIt process found ";
428  G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessor0005",
429  FatalErrorInArgument,exceptionDescription);
430  return;
431  }
432 
435  {
439 
441  {
442  G4ExceptionDescription exceptionDescription;
443  exceptionDescription << "No transportation process found ";
444  G4Exception("G4ITStepProcessor::SetupGeneralProcessInfo",
445  "ITStepProcessor0006",
446  FatalErrorInArgument,exceptionDescription);
447  }
448  }
450  // fpProcessInfo = 0;
451 }
G4ProcessVector * fpAtRestGetPhysIntVector
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4ITTransportation * fpTransportation
std::map< const G4ParticleDefinition *, ProcessGeneralInfo * > fProcessGeneralInfoMap
G4ProcessVector * fpAlongStepDoItVector
G4ProcessVector * fpPostStepDoItVector
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4int entries() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4ProcessVector * fpAlongStepGetPhysIntVector
G4ProcessVector * fpAtRestDoItVector
G4ProcessVector * fpPostStepGetPhysIntVector
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
ProcessGeneralInfo * fpProcessInfo
#define G4endl
Definition: G4ios.hh:61
static const size_t SizeOfSelectedDoItVector
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:

◆ SetupMembers()

void G4ITStepProcessor::SetupMembers ( )
protected

Definition at line 516 of file G4ITStepProcessor.cc.

517 {
518  fpSecondary = fpStep->GetfSecondary();
519  fpPreStepPoint = fpStep->GetPreStepPoint();
520  fpPostStepPoint = fpStep->GetPostStepPoint();
521 
524 
525  GetProcessInfo();
527 }
G4TrackVector * fpSecondary
G4StepPoint * fpPreStepPoint
G4ITStepProcessorState_Lock * GetStepProcessorState()
G4ITStepProcessorState * fpState
G4StepPoint * fpPostStepPoint
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:140

◆ Stepping()

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

Definition at line 265 of file G4ITStepProcessor2.cc.

266 {
267 
268 #ifdef DEBUG_MEM
269  MemStat mem_first, mem_second, mem_diff;
270 #endif
271 
272 #ifdef DEBUG_MEM
273  mem_first = MemoryUsage();
274 #endif
275 
276  CleanProcessor();
277 
278 #ifdef DEBUG_MEM
279  MemStat mem_intermediaire = MemoryUsage();
280  mem_diff = mem_intermediaire-mem_first;
281  G4cout << "\t\t\t >> || MEM || After CleanProcessor " << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
282 #endif
283 
284  if(track == 0) return; // maybe put an exception here
285  fTimeStep = timeStep;
286  SetTrack(track);
287  DoStepping();
288 }
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
void SetTrack(G4Track *)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

Friends And Related Function Documentation

◆ G4Scheduler

friend class G4Scheduler
friend

Definition at line 156 of file G4ITStepProcessor.hh.

Member Data Documentation

◆ fAtRestDoItProcTriggered

size_t G4ITStepProcessor::fAtRestDoItProcTriggered
private

Definition at line 378 of file G4ITStepProcessor.hh.

◆ fCondition

G4ForceCondition G4ITStepProcessor::fCondition
private

Definition at line 383 of file G4ITStepProcessor.hh.

◆ fGPILSelection

G4GPILSelection G4ITStepProcessor::fGPILSelection
private

Definition at line 384 of file G4ITStepProcessor.hh.

◆ fILTimeStep

G4double G4ITStepProcessor::fILTimeStep
private

Definition at line 360 of file G4ITStepProcessor.hh.

◆ fInitialized

G4bool G4ITStepProcessor::fInitialized
private

Definition at line 340 of file G4ITStepProcessor.hh.

◆ fLeadingTracks

G4ITLeadingTracks G4ITStepProcessor::fLeadingTracks
private

Definition at line 352 of file G4ITStepProcessor.hh.

◆ fN2ndariesAlongStepDoIt

G4int G4ITStepProcessor::fN2ndariesAlongStepDoIt
private

Definition at line 372 of file G4ITStepProcessor.hh.

◆ fN2ndariesAtRestDoIt

G4int G4ITStepProcessor::fN2ndariesAtRestDoIt
private

Definition at line 371 of file G4ITStepProcessor.hh.

◆ fN2ndariesPostStepDoIt

G4int G4ITStepProcessor::fN2ndariesPostStepDoIt
private

Definition at line 373 of file G4ITStepProcessor.hh.

◆ fpCurrentProcess

G4VITProcess* G4ITStepProcessor::fpCurrentProcess
private

Definition at line 366 of file G4ITStepProcessor.hh.

◆ fpCurrentVolume

G4VPhysicalVolume* G4ITStepProcessor::fpCurrentVolume
private

Definition at line 397 of file G4ITStepProcessor.hh.

◆ fPhysIntLength

G4double G4ITStepProcessor::fPhysIntLength
private

Definition at line 390 of file G4ITStepProcessor.hh.

◆ fpITrack

G4IT* G4ITStepProcessor::fpITrack
private

Definition at line 416 of file G4ITStepProcessor.hh.

◆ fpNavigator

G4ITNavigator* G4ITStepProcessor::fpNavigator
private

Definition at line 347 of file G4ITStepProcessor.hh.

◆ fPostStepAtTimeDoItProcTriggered

size_t G4ITStepProcessor::fPostStepAtTimeDoItProcTriggered
private

Definition at line 380 of file G4ITStepProcessor.hh.

◆ fPostStepDoItProcTriggered

size_t G4ITStepProcessor::fPostStepDoItProcTriggered
private

Definition at line 379 of file G4ITStepProcessor.hh.

◆ fpParticleChange

G4VParticleChange* G4ITStepProcessor::fpParticleChange
private

Definition at line 364 of file G4ITStepProcessor.hh.

◆ fpPostStepPoint

G4StepPoint* G4ITStepProcessor::fpPostStepPoint
private

Definition at line 423 of file G4ITStepProcessor.hh.

◆ fpPreStepPoint

G4StepPoint* G4ITStepProcessor::fpPreStepPoint
private

Definition at line 422 of file G4ITStepProcessor.hh.

◆ fpProcessInfo

ProcessGeneralInfo* G4ITStepProcessor::fpProcessInfo
private

Definition at line 407 of file G4ITStepProcessor.hh.

◆ fPreviousTimeStep

G4double G4ITStepProcessor::fPreviousTimeStep
private

Definition at line 362 of file G4ITStepProcessor.hh.

◆ fProcessGeneralInfoMap

std::map<const G4ParticleDefinition*, ProcessGeneralInfo*> G4ITStepProcessor::fProcessGeneralInfoMap
private

Definition at line 406 of file G4ITStepProcessor.hh.

◆ fpSecondary

G4TrackVector* G4ITStepProcessor::fpSecondary
private

Definition at line 363 of file G4ITStepProcessor.hh.

◆ fpState

G4ITStepProcessorState* G4ITStepProcessor::fpState
private

Definition at line 419 of file G4ITStepProcessor.hh.

◆ fpStep

G4Step* G4ITStepProcessor::fpStep
private

Definition at line 420 of file G4ITStepProcessor.hh.

◆ fpTrack

G4Track* G4ITStepProcessor::fpTrack
private

Definition at line 415 of file G4ITStepProcessor.hh.

◆ fpTrackContainer

G4ITTrackHolder* G4ITStepProcessor::fpTrackContainer
private

Definition at line 351 of file G4ITStepProcessor.hh.

◆ fpTrackingInfo

G4TrackingInformation* G4ITStepProcessor::fpTrackingInfo
private

Definition at line 417 of file G4ITStepProcessor.hh.

◆ fpTrackingManager

G4ITTrackingManager* G4ITStepProcessor::fpTrackingManager
private

Definition at line 342 of file G4ITStepProcessor.hh.

◆ fpTransportation

G4ITTransportation* G4ITStepProcessor::fpTransportation
private

Definition at line 408 of file G4ITStepProcessor.hh.

◆ fpVerbose

G4VITSteppingVerbose* G4ITStepProcessor::fpVerbose
private

Definition at line 349 of file G4ITStepProcessor.hh.

◆ fStoreTrajectory

G4int G4ITStepProcessor::fStoreTrajectory
private

Definition at line 348 of file G4ITStepProcessor.hh.

◆ fTimeStep

G4double G4ITStepProcessor::fTimeStep
private

Definition at line 359 of file G4ITStepProcessor.hh.

◆ kCarTolerance

G4double G4ITStepProcessor::kCarTolerance
private

Definition at line 344 of file G4ITStepProcessor.hh.


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