File indexing completed on 2025-01-18 09:58:33
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046 #ifndef G4ITSTEPPROCESSOR_H
0047 #define G4ITSTEPPROCESSOR_H
0048
0049 #include "G4ios.hh" // Include from 'system'
0050 #include "globals.hh" // Include from 'global'
0051 #include "Randomize.hh" // Include from 'global'
0052
0053 #include "G4LogicalVolume.hh" // Include from 'geometry'
0054 #include "G4VPhysicalVolume.hh" // Include from 'geometry'
0055 #include "G4ProcessManager.hh" // Include from 'piim'
0056
0057 #include "G4Track.hh" // Include from 'track'
0058 #include "G4TrackVector.hh" // Include from 'track'
0059 #include "G4TrackStatus.hh" // Include from 'track'
0060 #include "G4StepStatus.hh" // Include from 'track'
0061 #include "G4Step.hh" // Include from 'track'
0062 #include "G4StepPoint.hh" // Include from 'track'
0063 #include "G4TouchableHandle.hh" // Include from 'geometry'
0064
0065 #include "G4ITStepProcessorState_Lock.hh"
0066 #include "G4ITLeadingTracks.hh"
0067
0068 #include <vector>
0069
0070 class G4ITNavigator;
0071 class G4ParticleDefinition;
0072 class G4ITTrackingManager;
0073 class G4IT;
0074 class G4TrackingInformation;
0075 class G4ITTransportation;
0076 class G4VITProcess;
0077 class G4VITSteppingVerbose;
0078 class G4ITTrackHolder;
0079 using G4SelectedAtRestDoItVector = std::vector<int, std::allocator<int>>;
0080 using G4SelectedAlongStepDoItVector = std::vector<int, std::allocator<int>>;
0081 using G4SelectedPostStepDoItVector = std::vector<int, std::allocator<int>>;
0082
0083
0084
0085
0086
0087
0088 struct ProcessGeneralInfo
0089 {
0090 G4ProcessVector* fpAtRestDoItVector;
0091 G4ProcessVector* fpAlongStepDoItVector;
0092 G4ProcessVector* fpPostStepDoItVector;
0093
0094 G4ProcessVector* fpAtRestGetPhysIntVector;
0095 G4ProcessVector* fpAlongStepGetPhysIntVector;
0096 G4ProcessVector* fpPostStepGetPhysIntVector;
0097
0098
0099
0100
0101
0102 std::size_t MAXofAtRestLoops;
0103 std::size_t MAXofAlongStepLoops;
0104 std::size_t MAXofPostStepLoops;
0105
0106
0107
0108
0109 G4ITTransportation* fpTransportation;
0110 };
0111
0112
0113
0114
0115
0116 class G4ITStepProcessorState : public G4ITStepProcessorState_Lock
0117 {
0118 public:
0119 G4ITStepProcessorState();
0120 ~G4ITStepProcessorState() override;
0121
0122 G4ITStepProcessorState(const G4ITStepProcessorState&);
0123 G4ITStepProcessorState& operator=(const G4ITStepProcessorState&);
0124
0125
0126 G4SelectedAtRestDoItVector fSelectedAtRestDoItVector;
0127 G4SelectedPostStepDoItVector fSelectedPostStepDoItVector;
0128
0129 G4double fPhysicalStep;
0130 G4double fPreviousStepSize;
0131 G4double fSafety;
0132
0133 G4StepStatus fStepStatus;
0134
0135
0136 G4double fProposedSafety;
0137
0138 G4ThreeVector fEndpointSafOrigin;
0139 G4double fEndpointSafety;
0140
0141
0142
0143 G4TouchableHandle fTouchableHandle;
0144 };
0145
0146
0147
0148
0149
0150
0151
0152 class G4ITStepProcessor
0153 {
0154 friend class G4Scheduler;
0155 public:
0156 G4ITStepProcessor();
0157 virtual ~G4ITStepProcessor();
0158
0159 inline void SetPreviousStepTime(G4double);
0160
0161 inline G4Track* GetTrack()
0162 {
0163 return fpTrack;
0164 }
0165 inline G4Step* GetStep()
0166 {
0167 return fpStep;
0168 }
0169 inline const G4Step* GetStep() const
0170 {
0171 return fpStep;
0172 }
0173 inline void SetStep(G4Step* val)
0174 {
0175 fpStep = val;
0176 }
0177
0178 inline G4TrackVector* GetSecondaries() const
0179 {
0180 return fpSecondary;
0181 }
0182 inline void SetTrackingManager(G4ITTrackingManager* trackMan)
0183 {
0184 fpTrackingManager = trackMan;
0185 }
0186 inline G4ITTrackingManager* GetTrackingManager()
0187 {
0188 return fpTrackingManager;
0189 }
0190
0191
0192
0193 virtual void Initialize();
0194 void ForceReInitialization();
0195
0196 void ResetLeadingTracks();
0197 void PrepareLeadingTracks();
0198
0199
0200 G4double ComputeInteractionLength(double previousTimeStep);
0201 void DefinePhysicalStepLength(G4Track*);
0202 G4double GetILTimeStep()
0203 {
0204 return fILTimeStep;
0205 }
0206
0207
0208
0209 void DoIt(double timeStep);
0210 void ExtractDoItData();
0211 void Stepping(G4Track*, const double&);
0212 void FindTransportationStep();
0213
0214
0215 inline double GetInteractionTime();
0216 inline const G4Track* GetTrack() const;
0217 inline void CleanProcessor();
0218
0219 std::size_t GetAtRestDoItProcTriggered() const
0220 {
0221 return fAtRestDoItProcTriggered;
0222 }
0223
0224 G4GPILSelection GetGPILSelection() const
0225 {
0226 return fGPILSelection;
0227 }
0228
0229 G4int GetN2ndariesAlongStepDoIt() const
0230 {
0231 return fN2ndariesAlongStepDoIt;
0232 }
0233
0234 G4int GetN2ndariesAtRestDoIt() const
0235 {
0236 return fN2ndariesAtRestDoIt;
0237 }
0238
0239 G4int GetN2ndariesPostStepDoIt() const
0240 {
0241 return fN2ndariesPostStepDoIt;
0242 }
0243
0244 const G4VITProcess* GetCurrentProcess() const
0245 {
0246 return fpCurrentProcess;
0247 }
0248
0249 G4double GetPhysIntLength() const
0250 {
0251 return fPhysIntLength;
0252 }
0253
0254 std::size_t GetPostStepAtTimeDoItProcTriggered() const
0255 {
0256 return fPostStepAtTimeDoItProcTriggered;
0257 }
0258
0259 std::size_t GetPostStepDoItProcTriggered() const
0260 {
0261 return fPostStepDoItProcTriggered;
0262 }
0263
0264 const ProcessGeneralInfo* GetCurrentProcessInfo() const
0265 {
0266 return fpProcessInfo;
0267 }
0268
0269 const G4ITStepProcessorState* GetProcessorState() const
0270 {
0271 return fpState;
0272 }
0273
0274 const G4VParticleChange* GetParticleChange() const
0275 {
0276 return fpParticleChange;
0277 }
0278
0279 const G4VPhysicalVolume* GetCurrentVolume() const
0280 {
0281 return fpCurrentVolume;
0282 }
0283
0284 G4ForceCondition GetCondition() const
0285 {
0286 return fCondition;
0287 }
0288
0289 protected:
0290
0291 void ExtractILData();
0292
0293 void SetupGeneralProcessInfo(G4ParticleDefinition*, G4ProcessManager*);
0294 void ClearProcessInfo();
0295 void SetTrack(G4Track*);
0296
0297 void GetProcessInfo();
0298
0299 void SetupMembers();
0300 void ResetSecondaries();
0301 void InitDefineStep();
0302
0303 void SetInitialStep();
0304
0305 void GetAtRestIL();
0306 void DoDefinePhysicalStepLength();
0307 void DoStepping();
0308 void PushSecondaries();
0309
0310
0311
0312
0313
0314 void ActiveOnlyITProcess();
0315 void ActiveOnlyITProcess(G4ProcessManager*);
0316
0317 void DealWithSecondaries(G4int&);
0318 void InvokeAtRestDoItProcs();
0319 void InvokeAlongStepDoItProcs();
0320 void InvokePostStepDoItProcs();
0321 void InvokePSDIP(std::size_t);
0322 void InvokeTransportationProc();
0323 void SetNavigator(G4ITNavigator *value);
0324 G4double CalculateSafety();
0325
0326
0327 void ApplyProductionCut(G4Track*);
0328
0329 G4ITStepProcessor(const G4ITStepProcessor& other);
0330 G4ITStepProcessor& operator=(const G4ITStepProcessor& other);
0331
0332 private:
0333
0334
0335
0336
0337
0338 G4bool fInitialized;
0339
0340 G4ITTrackingManager* fpTrackingManager;
0341
0342 G4double kCarTolerance;
0343
0344
0345 G4ITNavigator* fpNavigator;
0346 G4int fStoreTrajectory;
0347 G4VITSteppingVerbose* fpVerbose;
0348
0349 G4ITTrackHolder* fpTrackContainer;
0350 G4ITLeadingTracks fLeadingTracks;
0351
0352
0353
0354
0355
0356
0357 G4double fTimeStep;
0358 G4double fILTimeStep;
0359
0360 G4double fPreviousTimeStep;
0361 G4TrackVector* fpSecondary;
0362 G4VParticleChange* fpParticleChange;
0363
0364 G4VITProcess* fpCurrentProcess;
0365
0366
0367
0368
0369 G4int fN2ndariesAtRestDoIt;
0370 G4int fN2ndariesAlongStepDoIt;
0371 G4int fN2ndariesPostStepDoIt;
0372
0373
0374
0375
0376 std::size_t fAtRestDoItProcTriggered;
0377 std::size_t fPostStepDoItProcTriggered;
0378 std::size_t fPostStepAtTimeDoItProcTriggered;
0379
0380
0381 G4ForceCondition fCondition;
0382 G4GPILSelection fGPILSelection;
0383
0384
0385
0386
0387
0388 G4double fPhysIntLength;
0389
0390
0391
0392
0393
0394
0395 G4VPhysicalVolume* fpCurrentVolume;
0396
0397
0398
0399
0400
0401
0402
0403
0404 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*> fProcessGeneralInfoMap;
0405 ProcessGeneralInfo* fpProcessInfo;
0406 G4ITTransportation* fpTransportation;
0407
0408
0409
0410
0411
0412
0413 G4Track* fpTrack;
0414 G4IT* fpITrack;
0415 G4TrackingInformation* fpTrackingInfo;
0416
0417 G4ITStepProcessorState* fpState;
0418 G4Step* fpStep;
0419
0420 G4StepPoint* fpPreStepPoint;
0421 G4StepPoint* fpPostStepPoint;
0422 };
0423
0424
0425
0426 inline void G4ITStepProcessor::SetPreviousStepTime(G4double previousTimeStep)
0427 {
0428 fPreviousTimeStep = previousTimeStep;
0429 }
0430
0431
0432
0433 inline const G4Track* G4ITStepProcessor::GetTrack() const
0434 {
0435 return fpTrack;
0436 }
0437
0438
0439
0440 inline G4double G4ITStepProcessor::CalculateSafety()
0441 {
0442 return std::max(fpState->fEndpointSafety - (fpState->fEndpointSafOrigin
0443 - fpPostStepPoint->GetPosition()).mag(),
0444 kCarTolerance);
0445 }
0446
0447
0448
0449 inline void G4ITStepProcessor::SetNavigator(G4ITNavigator *value)
0450 {
0451 fpNavigator = value;
0452 }
0453
0454
0455
0456 inline void G4ITStepProcessor::CleanProcessor()
0457 {
0458 fTimeStep = DBL_MAX;
0459 fPhysIntLength = DBL_MAX;
0460
0461 fpState = nullptr;
0462 fpTrack = nullptr;
0463 fpTrackingInfo = nullptr;
0464 fpITrack = nullptr;
0465 fpStep = nullptr;
0466 fpPreStepPoint = nullptr;
0467 fpPostStepPoint = nullptr;
0468
0469 fpParticleChange = nullptr;
0470
0471 fpCurrentVolume = nullptr;
0472
0473
0474 fpSecondary = nullptr;
0475
0476 fpTransportation = nullptr;
0477
0478 fpCurrentProcess= nullptr;
0479 fpProcessInfo = nullptr;
0480
0481 fAtRestDoItProcTriggered = INT_MAX;
0482 fPostStepDoItProcTriggered = INT_MAX;
0483 fPostStepAtTimeDoItProcTriggered = INT_MAX;
0484 fGPILSelection = NotCandidateForSelection;
0485 fCondition = NotForced;
0486 }
0487
0488
0489
0490 inline double G4ITStepProcessor::GetInteractionTime()
0491 {
0492 return fTimeStep;
0493 }
0494
0495 #endif