File indexing completed on 2025-01-18 09:59:24
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 G4VITProcess_H
0047 #define G4VITProcess_H
0048
0049 #include <G4VProcess.hh>
0050 #include "AddClone_def.hh"
0051 #include "G4ReferenceCast.hh"
0052 #include "G4memory.hh"
0053 #include <typeinfo>
0054
0055 class G4IT;
0056 class G4TrackingInformation;
0057
0058 struct G4ProcessState_Lock
0059 {
0060 inline virtual ~G4ProcessState_Lock()
0061 {
0062 ;
0063 }
0064 };
0065
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076 #define InitProcessState(destinationType,source) \
0077 reference_cast<destinationType>(source)
0078
0079 #define DowncastProcessState(destinationType) \
0080 G4dynamic_pointer_cast<destinationType>(G4VITProcess::fpState)
0081
0082 #define UpcastProcessState(destinationType) \
0083 G4dynamic_pointer_cast<destinationType>(G4VITProcess::fpState)
0084
0085 #define DowncastState(destinationType,source) \
0086 G4dynamic_pointer_cast<destinationType>(source)
0087
0088 #define UpcastState(destinationType,source) \
0089 G4dynamic_pointer_cast<destinationType>(source)
0090
0091
0092
0093
0094
0095
0096
0097
0098 class G4VITProcess : public G4VProcess
0099 {
0100 public:
0101
0102
0103 G4VITProcess(const G4String& name, G4ProcessType type = fNotDefined);
0104
0105 ~G4VITProcess() override;
0106 G4VITProcess(const G4VITProcess& other);
0107 G4VITProcess& operator=(const G4VITProcess& other);
0108
0109
0110 G4bool operator==(const G4VITProcess &right) const;
0111 G4bool operator!=(const G4VITProcess &right) const;
0112
0113 G4IT_TO_BE_CLONED(G4VITProcess)
0114
0115 size_t GetProcessID() const
0116 {
0117 return fProcessID;
0118 }
0119
0120 G4shared_ptr<G4ProcessState_Lock> GetProcessState()
0121 {
0122 return UpcastProcessState(G4ProcessState_Lock);
0123 }
0124
0125 void SetProcessState(G4shared_ptr<G4ProcessState_Lock> aProcInfo)
0126 {
0127 fpState = DowncastState(G4ProcessState, aProcInfo);
0128 }
0129
0130 void ResetProcessState()
0131 {
0132 fpState.reset();
0133 }
0134
0135
0136
0137
0138 void StartTracking(G4Track*) override;
0139
0140 void BuildPhysicsTable(const G4ParticleDefinition&) override
0141 {
0142 }
0143
0144 inline G4double GetInteractionTimeLeft();
0145
0146
0147
0148
0149 void ResetNumberOfInteractionLengthLeft() override;
0150
0151 inline G4bool ProposesTimeStep() const;
0152
0153 inline static const size_t& GetMaxProcessIndex();
0154
0155 protected:
0156
0157
0158 void RetrieveProcessInfo();
0159 void CreateInfo();
0160
0161
0162
0163
0164
0165 struct G4ProcessState : public G4ProcessState_Lock
0166 {
0167 public:
0168 G4ProcessState();
0169 ~G4ProcessState() override;
0170
0171 virtual G4String GetType()
0172 {
0173 return "G4ProcessState";
0174 }
0175
0176 G4double theNumberOfInteractionLengthLeft;
0177
0178
0179
0180 G4double theInteractionTimeLeft;
0181
0182
0183 G4double currentInteractionLength;
0184
0185
0186 template<typename T>
0187 T* GetState()
0188 {
0189 return dynamic_cast<T*>(this);
0190 }
0191 };
0192
0193 template<typename T>
0194 class G4ProcessStateBase : public G4ProcessState
0195 {
0196 public:
0197 G4ProcessStateBase() :
0198 G4ProcessState()
0199 {
0200 }
0201 ~G4ProcessStateBase() override
0202 = default;
0203
0204 G4String GetType() override
0205 {
0206 return typeid(T).name();
0207 }
0208 };
0209
0210 template<typename T>
0211 T* GetState()
0212 {
0213 return fpState->GetState<T>();
0214 }
0215
0216 G4shared_ptr<G4ProcessState> fpState;
0217
0218 void virtual SubtractNumberOfInteractionLengthLeft(G4double previousStepSize);
0219
0220 inline virtual void ClearInteractionTimeLeft();
0221
0222 inline virtual void ClearNumberOfInteractionLengthLeft();
0223
0224
0225
0226
0227
0228 void SetInstantiateProcessState(G4bool flag)
0229 {
0230 fInstantiateProcessState = flag;
0231 }
0232
0233 G4bool InstantiateProcessState()
0234 {
0235 return fInstantiateProcessState;
0236 }
0237
0238 G4bool fProposesTimeStep;
0239
0240 private:
0241
0242 size_t fProcessID;
0243
0244
0245
0246
0247 staticsize_t *fNbProcess;
0248
0249 G4bool fInstantiateProcessState;
0250
0251
0252 G4double* theNumberOfInteractionLengthLeft;
0253 G4double* currentInteractionLength;
0254 G4double* theInteractionTimeLeft;
0255 };
0256
0257 inline void G4VITProcess::ClearInteractionTimeLeft()
0258 {
0259 fpState->theInteractionTimeLeft = -1.0;
0260 }
0261
0262 inline void G4VITProcess::ClearNumberOfInteractionLengthLeft()
0263 {
0264 fpState->theNumberOfInteractionLengthLeft = -1.0;
0265 }
0266
0267 inline void G4VITProcess::ResetNumberOfInteractionLengthLeft()
0268 {
0269 fpState->theNumberOfInteractionLengthLeft = -std::log( G4UniformRand());
0270 }
0271
0272 inline G4double G4VITProcess::GetInteractionTimeLeft()
0273 {
0274 if (fpState) return fpState->theInteractionTimeLeft;
0275
0276 return -1;
0277 }
0278
0279 inline G4bool G4VITProcess::ProposesTimeStep() const
0280 {
0281 return fProposesTimeStep;
0282 }
0283
0284 inline const size_t& G4VITProcess::GetMaxProcessIndex()
0285 {
0286 if (fNbProcess == nullptr) fNbProcess = new size_t(0);
0287 return *fNbProcess;
0288 }
0289
0290 inline
0291 void G4VITProcess::SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
0292 {
0293 if (fpState->currentInteractionLength > 0.0)
0294 {
0295 fpState->theNumberOfInteractionLengthLeft -= previousStepSize
0296 / fpState->currentInteractionLength;
0297 if (fpState->theNumberOfInteractionLengthLeft < 0.)
0298 {
0299 fpState->theNumberOfInteractionLengthLeft = CLHEP::perMillion;
0300 }
0301
0302 }
0303 else
0304 {
0305 #ifdef G4VERBOSE
0306 if (verboseLevel > 0)
0307 {
0308 G4cerr << "G4VITProcess::SubtractNumberOfInteractionLengthLeft()";
0309 G4cerr << " [" << theProcessName << "]" << G4endl;
0310 G4cerr << " currentInteractionLength = "
0311 << fpState->currentInteractionLength << " [mm]";
0312 G4cerr << " previousStepSize = " << previousStepSize << " [mm]";
0313 G4cerr << G4endl;
0314 }
0315 #endif
0316 G4String msg = "Negative currentInteractionLength for ";
0317 msg += theProcessName;
0318 G4Exception("G4VITProcess::SubtractNumberOfInteractionLengthLeft()",
0319 "ProcMan201",EventMustBeAborted,
0320 msg);
0321 }
0322 }
0323 #endif