File indexing completed on 2025-01-18 09:58:12
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
0047
0048
0049
0050
0051
0052
0053
0054 #ifndef G4EmCalculator_h
0055 #define G4EmCalculator_h 1
0056
0057 #include <vector>
0058 #include "globals.hh"
0059 #include "G4VAtomDeexcitation.hh"
0060
0061 class G4LossTableManager;
0062 class G4NistManager;
0063 class G4Material;
0064 class G4MaterialCutsCouple;
0065 class G4ParticleDefinition;
0066 class G4PhysicsTable;
0067 class G4VEmModel;
0068 class G4VEnergyLossProcess;
0069 class G4VEmProcess;
0070 class G4VMultipleScattering;
0071 class G4VProcess;
0072 class G4ionEffectiveCharge;
0073 class G4Region;
0074 class G4Element;
0075 class G4EmCorrections;
0076 class G4EmParameters;
0077 class G4IonTable;
0078
0079 class G4EmCalculator
0080 {
0081
0082 public:
0083
0084 G4EmCalculator();
0085
0086 ~G4EmCalculator();
0087
0088
0089
0090
0091
0092
0093 G4double GetDEDX(G4double kinEnergy, const G4ParticleDefinition*,
0094 const G4Material*,
0095 const G4Region* r = nullptr);
0096 inline G4double GetDEDX(G4double kinEnergy, const G4String& part,
0097 const G4String& mat,
0098 const G4String& regname = "world");
0099
0100 G4double GetRangeFromRestricteDEDX(G4double kinEnergy,
0101 const G4ParticleDefinition*,
0102 const G4Material*,
0103 const G4Region* r = nullptr);
0104 inline G4double GetRangeFromRestricteDEDX(G4double kinEnergy,
0105 const G4String& part,
0106 const G4String& mat,
0107 const G4String& regname = "world");
0108
0109 G4double GetCSDARange(G4double kinEnergy, const G4ParticleDefinition*,
0110 const G4Material*,
0111 const G4Region* r = nullptr);
0112 inline G4double GetCSDARange(G4double kinEnergy, const G4String& part,
0113 const G4String& mat,
0114 const G4String& regname = "world");
0115
0116 G4double GetRange(G4double kinEnergy, const G4ParticleDefinition*,
0117 const G4Material*,
0118 const G4Region* r = nullptr);
0119 inline G4double GetRange(G4double kinEnergy, const G4String& part,
0120 const G4String& mat,
0121 const G4String& regname = "world");
0122
0123 G4double GetKinEnergy(G4double range, const G4ParticleDefinition*,
0124 const G4Material*,
0125 const G4Region* r = nullptr);
0126 inline G4double GetKinEnergy(G4double range, const G4String& part,
0127 const G4String& mat,
0128 const G4String& regname = "world");
0129
0130 G4double GetCrossSectionPerVolume(
0131 G4double kinEnergy, const G4ParticleDefinition*,
0132 const G4String& processName, const G4Material*,
0133 const G4Region* r = nullptr);
0134 inline G4double GetCrossSectionPerVolume(
0135 G4double kinEnergy, const G4String& part, const G4String& proc,
0136 const G4String& mat, const G4String& regname = "world");
0137
0138 G4double GetShellIonisationCrossSectionPerAtom(
0139 const G4String& part, G4int Z,
0140 G4AtomicShellEnumerator shell,
0141 G4double kinEnergy);
0142
0143 G4double GetMeanFreePath(G4double kinEnergy, const G4ParticleDefinition*,
0144 const G4String& processName, const G4Material*,
0145 const G4Region* r = nullptr);
0146 inline G4double GetMeanFreePath(G4double kinEnergy, const G4String& part,
0147 const G4String& proc, const G4String& mat,
0148 const G4String& regname = "world");
0149
0150 void PrintDEDXTable(const G4ParticleDefinition*);
0151
0152 void PrintRangeTable(const G4ParticleDefinition*);
0153
0154 void PrintInverseRangeTable(const G4ParticleDefinition*);
0155
0156
0157
0158
0159
0160
0161 G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition*,
0162 const G4String& processName, const G4Material*,
0163 G4double cut = DBL_MAX);
0164 inline G4double ComputeDEDX(G4double kinEnergy, const G4String& part,
0165 const G4String& proc,
0166 const G4String& mat, G4double cut = DBL_MAX);
0167
0168 G4double ComputeElectronicDEDX(G4double kinEnergy,
0169 const G4ParticleDefinition*,
0170 const G4Material* mat, G4double cut = DBL_MAX);
0171 inline G4double ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
0172 const G4String& mat, G4double cut = DBL_MAX);
0173
0174 G4double ComputeDEDXForCutInRange(G4double kinEnergy,
0175 const G4ParticleDefinition*,
0176 const G4Material* mat, G4double rangecut = DBL_MAX);
0177 inline G4double ComputeDEDXForCutInRange(G4double kinEnergy, const G4String& part,
0178 const G4String& mat,
0179 G4double rangecut = DBL_MAX);
0180
0181 G4double ComputeNuclearDEDX(G4double kinEnergy, const G4ParticleDefinition*,
0182 const G4Material*);
0183 inline G4double ComputeNuclearDEDX(G4double kinEnergy, const G4String& part,
0184 const G4String& mat);
0185
0186 G4double ComputeTotalDEDX(G4double kinEnergy, const G4ParticleDefinition*,
0187 const G4Material*, G4double cut = DBL_MAX);
0188 inline G4double ComputeTotalDEDX(G4double kinEnergy, const G4String& part,
0189 const G4String& mat, G4double cut = DBL_MAX);
0190
0191 G4double ComputeCrossSectionPerVolume(
0192 G4double kinEnergy, const G4ParticleDefinition*,
0193 const G4String& processName, const G4Material*,
0194 G4double cut = 0.0);
0195 inline G4double ComputeCrossSectionPerVolume(
0196 G4double kinEnergy, const G4String& part,
0197 const G4String& proc,
0198 const G4String& mat, G4double cut = 0.0);
0199
0200 G4double ComputeCrossSectionPerAtom(
0201 G4double kinEnergy, const G4ParticleDefinition*,
0202 const G4String& processName, G4double Z, G4double A,
0203 G4double cut = 0.0);
0204 inline G4double ComputeCrossSectionPerAtom(
0205 G4double kinEnergy, const G4String& part,
0206 const G4String& processName, const G4Element*,
0207 G4double cut = 0.0);
0208
0209 G4double ComputeCrossSectionPerShell(
0210 G4double kinEnergy, const G4ParticleDefinition*,
0211 const G4String& processName, G4int Z, G4int shellIdx,
0212 G4double cut = 0.0);
0213 inline G4double ComputeCrossSectionPerShell(
0214 G4double kinEnergy, const G4String& part,
0215 const G4String& processName, const G4Element*,
0216 G4int shellIdx,
0217 G4double cut = 0.0);
0218
0219 G4double ComputeGammaAttenuationLength(G4double kinEnergy,
0220 const G4Material*);
0221
0222 G4double ComputeShellIonisationCrossSectionPerAtom(
0223 const G4String& part, G4int Z,
0224 G4AtomicShellEnumerator shell,
0225 G4double kinEnergy,
0226 const G4Material* mat = nullptr);
0227
0228 G4double ComputeMeanFreePath(
0229 G4double kinEnergy, const G4ParticleDefinition*,
0230 const G4String& processName, const G4Material*,
0231 G4double cut = 0.0);
0232 inline G4double ComputeMeanFreePath(
0233 G4double kinEnergy, const G4String&, const G4String&,
0234 const G4String& processName, G4double cut = 0.0);
0235
0236 G4double ComputeEnergyCutFromRangeCut(
0237 G4double range, const G4ParticleDefinition*,
0238 const G4Material*);
0239 inline G4double ComputeEnergyCutFromRangeCut(
0240 G4double range, const G4String&,
0241 const G4String&);
0242
0243
0244
0245
0246
0247 const G4ParticleDefinition* FindParticle(const G4String&);
0248
0249 const G4ParticleDefinition* FindIon(G4int Z, G4int A);
0250
0251 const G4Material* FindMaterial(const G4String&);
0252
0253 const G4Region* FindRegion(const G4String&);
0254
0255 const G4MaterialCutsCouple* FindCouple(const G4Material*,
0256 const G4Region* r = nullptr);
0257
0258 G4VProcess* FindProcess(const G4ParticleDefinition* part,
0259 const G4String& processName);
0260
0261 void SetupMaterial(const G4Material*);
0262
0263 void SetupMaterial(const G4String&);
0264
0265 void SetVerbose(G4int val);
0266
0267 inline void SetApplySmoothing(G4int val);
0268
0269
0270 G4EmCalculator & operator=(const G4EmCalculator &right) = delete;
0271 G4EmCalculator(const G4EmCalculator&) = delete;
0272
0273 private:
0274
0275 G4bool UpdateParticle(const G4ParticleDefinition*, G4double kinEnergy);
0276
0277 G4bool UpdateCouple(const G4Material*, G4double cut);
0278
0279 void FindLambdaTable(const G4ParticleDefinition*,
0280 const G4String& processName,
0281 G4double kinEnergy, G4int& proctype);
0282
0283 G4bool FindEmModel(const G4ParticleDefinition*,
0284 const G4String& processName,
0285 G4double kinEnergy);
0286
0287 G4VEnergyLossProcess* FindEnLossProcess(const G4ParticleDefinition*,
0288 const G4String& processName);
0289
0290 G4VEmProcess* FindDiscreteProcess(const G4ParticleDefinition*,
0291 const G4String& processName);
0292
0293 G4VMultipleScattering* FindMscProcess(const G4ParticleDefinition*,
0294 const G4String& processName);
0295
0296 G4bool ActiveForParticle(const G4ParticleDefinition* part,
0297 G4VProcess* proc);
0298
0299 void CheckMaterial(G4int Z);
0300
0301 G4EmParameters* theParameters;
0302 G4LossTableManager* manager;
0303 G4NistManager* nist;
0304 G4IonTable* ionTable;
0305 G4EmCorrections* corr;
0306
0307
0308 const G4MaterialCutsCouple* currentCouple = nullptr;
0309 const G4Material* currentMaterial = nullptr;
0310 const G4Material* cutMaterial = nullptr;
0311 const G4ParticleDefinition* currentParticle = nullptr;
0312 const G4ParticleDefinition* lambdaParticle = nullptr;
0313 const G4ParticleDefinition* baseParticle = nullptr;
0314 const G4PhysicsTable* currentLambda = nullptr;
0315
0316 G4VEmModel* currentModel = nullptr;
0317 G4VEmModel* loweModel = nullptr;
0318 G4VEnergyLossProcess* currentProcess = nullptr;
0319 G4VProcess* curProcess = nullptr;
0320 G4DynamicParticle* dynParticle = nullptr;
0321
0322 const G4ParticleDefinition* theGenericIon;
0323 G4ionEffectiveCharge* ionEffCharge;
0324
0325 G4double currentCut = DBL_MAX;
0326 G4double chargeSquare = 1.0;
0327 G4double massRatio = 1.0;
0328 G4double mass = 0;
0329 G4double cutenergy[3];
0330
0331 G4int currentCoupleIndex = 0;
0332 G4int nLocalMaterials = 0;
0333 G4int verbose = 0;
0334
0335 G4bool isIon = false;
0336 G4bool isApplicable = false;
0337 G4bool applySmoothing = true;
0338
0339 std::vector<const G4Material*> localMaterials;
0340 std::vector<const G4MaterialCutsCouple*> localCouples;
0341 std::vector<G4double> localCuts;
0342
0343 G4String currentName = "";
0344 G4String lambdaName = "";
0345 G4String currentParticleName = "";
0346 G4String currentMaterialName = "";
0347 G4String currentProcessName = "";
0348 };
0349
0350
0351
0352
0353 inline
0354 G4double G4EmCalculator::GetDEDX(G4double kinEnergy, const G4String& particle,
0355 const G4String& material, const G4String& reg)
0356 {
0357 return GetDEDX(kinEnergy,FindParticle(particle),
0358 FindMaterial(material),FindRegion(reg));
0359 }
0360
0361
0362
0363 inline
0364 G4double G4EmCalculator::GetRangeFromRestricteDEDX(G4double kinEnergy,
0365 const G4String& particle,
0366 const G4String& material,
0367 const G4String& reg)
0368 {
0369 return GetRangeFromRestricteDEDX(kinEnergy,FindParticle(particle),
0370 FindMaterial(material),FindRegion(reg));
0371 }
0372
0373
0374
0375 inline
0376 G4double G4EmCalculator::GetCSDARange(G4double kinEnergy,
0377 const G4String& particle,
0378 const G4String& material,
0379 const G4String& reg)
0380 {
0381 return GetCSDARange(kinEnergy,FindParticle(particle),
0382 FindMaterial(material),FindRegion(reg));
0383 }
0384
0385
0386
0387 inline
0388 G4double G4EmCalculator::GetRange(G4double kinEnergy,
0389 const G4String& particle,
0390 const G4String& material,
0391 const G4String& reg)
0392 {
0393 return GetRange(kinEnergy,FindParticle(particle),
0394 FindMaterial(material),FindRegion(reg));
0395 }
0396
0397
0398
0399 inline
0400 G4double G4EmCalculator::GetKinEnergy(G4double range, const G4String& particle,
0401 const G4String& material, const G4String& reg)
0402 {
0403 return GetKinEnergy(range,FindParticle(particle),
0404 FindMaterial(material),FindRegion(reg));
0405 }
0406
0407
0408
0409 inline
0410 G4double G4EmCalculator::GetCrossSectionPerVolume(G4double kinEnergy,
0411 const G4String& particle,
0412 const G4String& processName,
0413 const G4String& material,
0414 const G4String& reg)
0415 {
0416 return GetCrossSectionPerVolume(kinEnergy,FindParticle(particle),processName,
0417 FindMaterial(material),FindRegion(reg));
0418 }
0419
0420
0421
0422 inline
0423 G4double G4EmCalculator::GetMeanFreePath(G4double kinEnergy,
0424 const G4String& particle,
0425 const G4String& processName,
0426 const G4String& material,
0427 const G4String& reg)
0428 {
0429 return GetMeanFreePath(kinEnergy,FindParticle(particle),processName,
0430 FindMaterial(material),FindRegion(reg));
0431 }
0432
0433
0434
0435 inline G4double
0436 G4EmCalculator::ComputeElectronicDEDX(G4double kinEnergy, const G4String& part,
0437 const G4String& mat, G4double cut)
0438 {
0439 return
0440 ComputeElectronicDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
0441 }
0442
0443
0444
0445 inline G4double
0446 G4EmCalculator::ComputeDEDXForCutInRange(G4double kinEnergy,
0447 const G4String& part,
0448 const G4String& mat,
0449 G4double rangecut)
0450 {
0451 return ComputeDEDXForCutInRange(kinEnergy,FindParticle(part),
0452 FindMaterial(mat), rangecut);
0453 }
0454
0455
0456
0457 inline
0458 G4double G4EmCalculator::ComputeTotalDEDX(G4double kinEnergy,
0459 const G4String& part,
0460 const G4String& mat,
0461 G4double cut)
0462 {
0463 return ComputeTotalDEDX(kinEnergy,FindParticle(part),FindMaterial(mat),cut);
0464 }
0465
0466
0467
0468 inline
0469 G4double G4EmCalculator::ComputeDEDX(G4double kinEnergy,
0470 const G4String& particle,
0471 const G4String& processName,
0472 const G4String& material,
0473 G4double cut)
0474 {
0475 return ComputeDEDX(kinEnergy,FindParticle(particle),processName,
0476 FindMaterial(material),cut);
0477 }
0478
0479
0480
0481 inline
0482 G4double G4EmCalculator::ComputeNuclearDEDX(G4double kinEnergy,
0483 const G4String& particle,
0484 const G4String& material)
0485 {
0486 return ComputeNuclearDEDX(kinEnergy,FindParticle(particle),
0487 FindMaterial(material));
0488 }
0489
0490
0491
0492 inline
0493 G4double G4EmCalculator::ComputeCrossSectionPerVolume(
0494 G4double kinEnergy,
0495 const G4String& particle,
0496 const G4String& processName,
0497 const G4String& material,
0498 G4double cut)
0499 {
0500 return ComputeCrossSectionPerVolume(kinEnergy,FindParticle(particle),
0501 processName,
0502 FindMaterial(material),cut);
0503 }
0504
0505
0506
0507 inline
0508 G4double G4EmCalculator::ComputeCrossSectionPerAtom(G4double kinEnergy,
0509 const G4String& particle,
0510 const G4String& processName,
0511 const G4Element* elm,
0512 G4double cut)
0513 {
0514 return ComputeCrossSectionPerAtom(kinEnergy,FindParticle(particle),
0515 processName,
0516 elm->GetZ(),elm->GetN(),cut);
0517 }
0518
0519
0520
0521 inline G4double G4EmCalculator::ComputeCrossSectionPerShell(
0522 G4double kinEnergy, const G4String& part,
0523 const G4String& processName, const G4Element* elm,
0524 G4int shellIdx, G4double cut)
0525 {
0526 return ComputeCrossSectionPerShell(kinEnergy, FindParticle(part),
0527 processName, elm->GetZasInt(),
0528 shellIdx, cut);
0529 }
0530
0531
0532
0533 inline
0534 G4double G4EmCalculator::ComputeEnergyCutFromRangeCut(
0535 G4double range,
0536 const G4String& particle,
0537 const G4String& material)
0538 {
0539 return ComputeEnergyCutFromRangeCut(range,FindParticle(particle),
0540 FindMaterial(material));
0541 }
0542
0543
0544
0545 inline
0546 G4double G4EmCalculator::ComputeMeanFreePath(G4double kinEnergy,
0547 const G4String& particle,
0548 const G4String& processName,
0549 const G4String& material,
0550 G4double cut)
0551 {
0552 return ComputeMeanFreePath(kinEnergy,FindParticle(particle),processName,
0553 FindMaterial(material),cut);
0554 }
0555
0556
0557
0558 inline void G4EmCalculator::SetApplySmoothing(G4int val)
0559 {
0560 applySmoothing = val;
0561 }
0562
0563
0564
0565 #endif