File indexing completed on 2025-01-18 09:59:26
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 #ifndef G4VProcess_hh
0039 #define G4VProcess_hh 1
0040
0041 #include <cmath>
0042
0043 #include "globals.hh"
0044 #include "G4ios.hh"
0045 #include "Randomize.hh"
0046
0047 #include "G4PhysicsTable.hh"
0048 #include "G4VParticleChange.hh"
0049 #include "G4ForceCondition.hh"
0050 #include "G4GPILSelection.hh"
0051 #include "G4ParticleChange.hh"
0052 #include "G4ProcessType.hh"
0053
0054 class G4ParticleDefinition;
0055 class G4DynamicParticle;
0056 class G4Track;
0057 class G4Step;
0058 class G4ProcessTable;
0059
0060 class G4VProcess
0061 {
0062
0063 public:
0064
0065 G4VProcess(const G4String& aName = "NoName",
0066 G4ProcessType aType = fNotDefined);
0067
0068
0069 G4VProcess(const G4VProcess& right);
0070
0071
0072
0073 virtual ~G4VProcess();
0074
0075
0076 G4VProcess& operator=(const G4VProcess&) = delete;
0077
0078 G4bool operator==(const G4VProcess& right) const;
0079 G4bool operator!=(const G4VProcess& right) const;
0080
0081
0082
0083
0084
0085
0086 virtual G4VParticleChange* PostStepDoIt(
0087 const G4Track& track,
0088 const G4Step& stepData
0089 ) = 0;
0090
0091 virtual G4VParticleChange* AlongStepDoIt(
0092 const G4Track& track,
0093 const G4Step& stepData
0094 ) = 0;
0095 virtual G4VParticleChange* AtRestDoIt(
0096 const G4Track& track,
0097 const G4Step& stepData
0098 ) = 0;
0099
0100
0101
0102
0103
0104
0105
0106
0107
0108
0109
0110
0111
0112
0113 virtual G4double AlongStepGetPhysicalInteractionLength(
0114 const G4Track& track,
0115 G4double previousStepSize,
0116 G4double currentMinimumStep,
0117 G4double& proposedSafety,
0118 G4GPILSelection* selection) = 0;
0119
0120 virtual G4double AtRestGetPhysicalInteractionLength(
0121 const G4Track& track,
0122 G4ForceCondition* condition ) = 0;
0123
0124 virtual G4double PostStepGetPhysicalInteractionLength(
0125 const G4Track& track,
0126 G4double previousStepSize,
0127 G4ForceCondition* condition ) = 0;
0128
0129
0130
0131
0132
0133
0134
0135
0136
0137
0138
0139
0140
0141
0142
0143
0144
0145
0146
0147
0148
0149
0150
0151
0152
0153
0154 inline G4double GetCurrentInteractionLength() const;
0155
0156
0157
0158
0159 inline void SetPILfactor(G4double value);
0160 inline G4double GetPILfactor() const;
0161
0162
0163
0164
0165
0166
0167
0168
0169 inline G4double AlongStepGPIL( const G4Track& track,
0170 G4double previousStepSize,
0171 G4double currentMinimumStep,
0172 G4double& proposedSafety,
0173 G4GPILSelection* selection );
0174
0175 inline G4double AtRestGPIL( const G4Track& track,
0176 G4ForceCondition* condition );
0177
0178 inline G4double PostStepGPIL( const G4Track& track,
0179 G4double previousStepSize,
0180 G4ForceCondition* condition );
0181
0182 virtual G4bool IsApplicable(const G4ParticleDefinition&) { return true; }
0183
0184
0185
0186
0187 virtual void BuildPhysicsTable(const G4ParticleDefinition&) {}
0188
0189
0190
0191
0192
0193
0194 virtual void PreparePhysicsTable(const G4ParticleDefinition&) {}
0195
0196
0197
0198
0199
0200
0201
0202
0203
0204
0205
0206 virtual G4bool StorePhysicsTable(const G4ParticleDefinition* ,
0207 const G4String&, G4bool) { return true; }
0208
0209
0210
0211 virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition* ,
0212 const G4String&, G4bool) { return false; }
0213
0214
0215
0216
0217
0218
0219 const G4String& GetPhysicsTableFileName(const G4ParticleDefinition* ,
0220 const G4String& directory,
0221 const G4String& tableName,
0222 G4bool ascii = false);
0223
0224
0225 inline const G4String& GetProcessName() const;
0226
0227
0228 inline G4ProcessType GetProcessType() const;
0229
0230
0231 inline void SetProcessType(G4ProcessType);
0232
0233
0234 inline G4int GetProcessSubType() const;
0235
0236
0237 inline void SetProcessSubType(G4int);
0238
0239
0240 static const G4String& GetProcessTypeName(G4ProcessType);
0241
0242
0243 virtual const G4VProcess* GetCreatorProcess() const;
0244
0245
0246
0247 virtual void StartTracking(G4Track*);
0248 virtual void EndTracking();
0249
0250
0251 virtual void SetProcessManager(const G4ProcessManager*);
0252
0253
0254 virtual const G4ProcessManager* GetProcessManager();
0255
0256
0257 virtual void ResetNumberOfInteractionLengthLeft();
0258
0259
0260 inline G4double GetNumberOfInteractionLengthLeft() const;
0261
0262
0263 inline G4double GetTotalNumberOfInteractionLengthTraversed() const;
0264
0265
0266
0267 inline G4bool isAtRestDoItIsEnabled() const;
0268 inline G4bool isAlongStepDoItIsEnabled() const;
0269 inline G4bool isPostStepDoItIsEnabled() const;
0270
0271
0272
0273
0274 virtual void DumpInfo() const;
0275
0276
0277 virtual void ProcessDescription(std::ostream& outfile) const;
0278
0279
0280 inline void SetVerboseLevel(G4int value);
0281 inline G4int GetVerboseLevel() const;
0282
0283
0284
0285
0286
0287 virtual void SetMasterProcess(G4VProcess* masterP);
0288
0289 inline const G4VProcess* GetMasterProcess() const;
0290
0291
0292
0293
0294
0295
0296 virtual void BuildWorkerPhysicsTable(const G4ParticleDefinition& part);
0297
0298
0299
0300
0301
0302
0303
0304 virtual void PrepareWorkerPhysicsTable(const G4ParticleDefinition&);
0305
0306
0307
0308
0309
0310
0311
0312 protected:
0313
0314 inline void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize);
0315
0316
0317
0318 inline void ClearNumberOfInteractionLengthLeft();
0319
0320
0321 protected:
0322
0323 const G4ProcessManager* aProcessManager = nullptr;
0324
0325 G4VParticleChange* pParticleChange = nullptr;
0326
0327
0328
0329
0330
0331 G4ParticleChange aParticleChange;
0332
0333
0334
0335 G4double theNumberOfInteractionLengthLeft = -1.0;
0336
0337
0338
0339 G4double currentInteractionLength = -1.0;
0340
0341
0342 G4double theInitialNumberOfInteractionLength = -1.0;
0343
0344
0345 G4String theProcessName;
0346
0347
0348 G4String thePhysicsTableFileName;
0349
0350 G4ProcessType theProcessType = fNotDefined;
0351
0352
0353 G4int theProcessSubType = -1;
0354
0355
0356 G4double thePILfactor = 1.0;
0357
0358
0359
0360 G4int verboseLevel = 0;
0361
0362
0363 G4bool enableAtRestDoIt = true;
0364 G4bool enableAlongStepDoIt = true;
0365 G4bool enablePostStepDoIt = true;
0366
0367 private:
0368
0369 G4VProcess();
0370
0371
0372 private:
0373
0374 G4VProcess* masterProcessShadow = nullptr;
0375
0376
0377
0378 G4ProcessTable* fProcessTable = nullptr;
0379 };
0380
0381
0382
0383
0384
0385 inline
0386 const G4String& G4VProcess::GetProcessName() const
0387 {
0388 return theProcessName;
0389 }
0390
0391 inline
0392 G4ProcessType G4VProcess::GetProcessType() const
0393 {
0394 return theProcessType;
0395 }
0396
0397 inline
0398 void G4VProcess::SetProcessType(G4ProcessType aType)
0399 {
0400 theProcessType = aType;
0401 }
0402
0403 inline
0404 G4int G4VProcess::GetProcessSubType() const
0405 {
0406 return theProcessSubType;
0407 }
0408
0409 inline
0410 void G4VProcess::SetProcessSubType(G4int value)
0411 {
0412 theProcessSubType = value;
0413 }
0414
0415 inline
0416 void G4VProcess::SetVerboseLevel(G4int value)
0417 {
0418 verboseLevel = value;
0419 }
0420
0421 inline
0422 G4int G4VProcess::GetVerboseLevel() const
0423 {
0424 return verboseLevel;
0425 }
0426
0427 inline
0428 void G4VProcess::ClearNumberOfInteractionLengthLeft()
0429 {
0430 theInitialNumberOfInteractionLength = -1.0;
0431 theNumberOfInteractionLengthLeft = -1.0;
0432 }
0433
0434 inline
0435 G4double G4VProcess::GetNumberOfInteractionLengthLeft() const
0436 {
0437 return theNumberOfInteractionLengthLeft;
0438 }
0439
0440 inline
0441 G4double G4VProcess::GetTotalNumberOfInteractionLengthTraversed() const
0442 {
0443 return theInitialNumberOfInteractionLength - theNumberOfInteractionLengthLeft;
0444 }
0445
0446 inline
0447 G4double G4VProcess::GetCurrentInteractionLength() const
0448 {
0449 return currentInteractionLength;
0450 }
0451
0452 inline
0453 void G4VProcess::SetPILfactor(G4double value)
0454 {
0455 if (value>0.) { thePILfactor = value; }
0456 }
0457
0458 inline
0459 G4double G4VProcess::GetPILfactor() const
0460 {
0461 return thePILfactor;
0462 }
0463
0464 inline
0465 G4double G4VProcess::AlongStepGPIL( const G4Track& track,
0466 G4double previousStepSize,
0467 G4double currentMinimumStep,
0468 G4double& proposedSafety,
0469 G4GPILSelection* selection )
0470 {
0471 return AlongStepGetPhysicalInteractionLength(track, previousStepSize,
0472 currentMinimumStep, proposedSafety, selection);
0473 }
0474
0475 inline
0476 G4double G4VProcess::AtRestGPIL( const G4Track& track,
0477 G4ForceCondition* condition )
0478 {
0479 return thePILfactor * AtRestGetPhysicalInteractionLength(track, condition);
0480 }
0481
0482 inline
0483 G4double G4VProcess::PostStepGPIL( const G4Track& track,
0484 G4double previousStepSize,
0485 G4ForceCondition* condition )
0486 {
0487 return thePILfactor *
0488 PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
0489 }
0490
0491 inline
0492 void G4VProcess::SetProcessManager(const G4ProcessManager* procMan)
0493 {
0494 aProcessManager = procMan;
0495 }
0496
0497 inline
0498 const G4ProcessManager* G4VProcess::GetProcessManager()
0499 {
0500 return aProcessManager;
0501 }
0502
0503 inline
0504 G4bool G4VProcess::isAtRestDoItIsEnabled() const
0505 {
0506 return enableAtRestDoIt;
0507 }
0508
0509 inline
0510 G4bool G4VProcess::isAlongStepDoItIsEnabled() const
0511 {
0512 return enableAlongStepDoIt;
0513 }
0514
0515 inline
0516 G4bool G4VProcess::isPostStepDoItIsEnabled() const
0517 {
0518 return enablePostStepDoIt;
0519 }
0520
0521 inline
0522 const G4VProcess* G4VProcess::GetMasterProcess() const
0523 {
0524 return masterProcessShadow;
0525 }
0526
0527 inline
0528 void G4VProcess::SubtractNumberOfInteractionLengthLeft( G4double prevStepSize )
0529 {
0530 if (currentInteractionLength>0.0)
0531 {
0532 theNumberOfInteractionLengthLeft -= prevStepSize/currentInteractionLength;
0533 if(theNumberOfInteractionLengthLeft<0.)
0534 {
0535 theNumberOfInteractionLengthLeft=CLHEP::perMillion;
0536 }
0537 }
0538 else
0539 {
0540 #ifdef G4VERBOSE
0541 if (verboseLevel>0)
0542 {
0543 G4cerr << "G4VProcess::SubtractNumberOfInteractionLengthLeft()";
0544 G4cerr << " [" << theProcessName << "]" <<G4endl;
0545 G4cerr << " currentInteractionLength = "
0546 << currentInteractionLength << " [mm]";
0547 G4cerr << " previousStepSize = " << prevStepSize << " [mm]";
0548 G4cerr << G4endl;
0549 }
0550 #endif
0551 G4String msg = "Negative currentInteractionLength for ";
0552 msg += theProcessName;
0553 G4Exception("G4VProcess::SubtractNumberOfInteractionLengthLeft()",
0554 "ProcMan201", EventMustBeAborted, msg);
0555 }
0556 }
0557
0558 #endif