File indexing completed on 2025-09-15 08:59:39
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 #ifndef G4VEnergyLossProcess_h
0052 #define G4VEnergyLossProcess_h 1
0053
0054 #include "G4VContinuousDiscreteProcess.hh"
0055 #include "globals.hh"
0056 #include "G4Material.hh"
0057 #include "G4MaterialCutsCouple.hh"
0058 #include "G4Track.hh"
0059 #include "G4EmModelManager.hh"
0060 #include "G4ParticleChangeForLoss.hh"
0061 #include "G4EmTableType.hh"
0062 #include "G4EmSecondaryParticleType.hh"
0063 #include "G4PhysicsTable.hh"
0064 #include "G4PhysicsVector.hh"
0065
0066 class G4Step;
0067 class G4ParticleDefinition;
0068 class G4EmParameters;
0069 class G4VEmModel;
0070 class G4VEmFluctuationModel;
0071 class G4DataVector;
0072 class G4Region;
0073 class G4SafetyHelper;
0074 class G4VAtomDeexcitation;
0075 class G4VSubCutProducer;
0076 class G4EmBiasingManager;
0077 class G4LossTableManager;
0078 class G4EmDataHandler;
0079
0080
0081
0082 class G4VEnergyLossProcess : public G4VContinuousDiscreteProcess
0083 {
0084 public:
0085
0086 G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
0087 G4ProcessType type = fElectromagnetic);
0088
0089 ~G4VEnergyLossProcess() override;
0090
0091
0092
0093
0094
0095 protected:
0096
0097
0098 virtual void StreamProcessInfo(std::ostream&) const {};
0099
0100 virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition*,
0101 const G4ParticleDefinition*) = 0;
0102
0103 public:
0104
0105
0106 virtual G4double MinPrimaryEnergy(const G4ParticleDefinition*,
0107 const G4Material*, G4double cut);
0108
0109
0110 void ProcessDescription(std::ostream& outFile) const override;
0111
0112
0113 void PreparePhysicsTable(const G4ParticleDefinition&) override;
0114
0115
0116 void BuildPhysicsTable(const G4ParticleDefinition&) override;
0117
0118
0119 G4PhysicsTable* BuildDEDXTable(G4EmTableType tType = fRestricted);
0120
0121
0122 G4PhysicsTable* BuildLambdaTable(G4EmTableType tType = fRestricted);
0123
0124
0125 void StartTracking(G4Track*) override;
0126
0127
0128 G4double AlongStepGetPhysicalInteractionLength(
0129 const G4Track&,
0130 G4double previousStepSize,
0131 G4double currentMinimumStep,
0132 G4double& currentSafety,
0133 G4GPILSelection* selection) override;
0134
0135
0136 G4double PostStepGetPhysicalInteractionLength(
0137 const G4Track& track,
0138 G4double previousStepSize,
0139 G4ForceCondition* condition) override;
0140
0141
0142 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
0143
0144
0145 G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
0146
0147
0148
0149 G4bool StorePhysicsTable(const G4ParticleDefinition*,
0150 const G4String& directory,
0151 G4bool ascii = false) override;
0152
0153
0154
0155
0156 G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
0157 const G4String& directory,
0158 G4bool ascii) override;
0159
0160 private:
0161
0162
0163 void StreamInfo(std::ostream& out, const G4ParticleDefinition& part,
0164 G4bool rst=false) const;
0165
0166
0167
0168
0169
0170
0171 public:
0172
0173
0174 G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple,
0175 const G4DynamicParticle* dp,
0176 G4double length);
0177
0178
0179 G4double CrossSectionPerVolume(G4double kineticEnergy,
0180 const G4MaterialCutsCouple* couple);
0181 G4double CrossSectionPerVolume(G4double kineticEnergy,
0182 const G4MaterialCutsCouple* couple,
0183 G4double logKineticEnergy);
0184
0185
0186 G4double MeanFreePath(const G4Track& track);
0187
0188
0189 G4double ContinuousStepLimit(const G4Track& track,
0190 G4double previousStepSize,
0191 G4double currentMinimumStep,
0192 G4double& currentSafety);
0193
0194 protected:
0195
0196
0197 G4double GetMeanFreePath(const G4Track& track,
0198 G4double previousStepSize,
0199 G4ForceCondition* condition) override;
0200
0201
0202 G4double GetContinuousStepLimit(const G4Track& track,
0203 G4double previousStepSize,
0204 G4double currentMinimumStep,
0205 G4double& currentSafety) override;
0206
0207
0208 G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*,
0209 G4double cut);
0210
0211 inline std::size_t CurrentMaterialCutsCoupleIndex() const;
0212
0213
0214
0215
0216
0217
0218 inline void SelectModel(G4double kinEnergy);
0219
0220 public:
0221
0222
0223 inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy,
0224 std::size_t& idxCouple) const;
0225
0226
0227
0228
0229 void AddEmModel(G4int, G4VEmModel*,
0230 G4VEmFluctuationModel* fluc = nullptr,
0231 const G4Region* region = nullptr);
0232
0233
0234
0235 void SetEmModel(G4VEmModel*, G4int index=0);
0236
0237
0238 inline std::size_t NumberOfModels() const;
0239
0240
0241 inline G4VEmModel* EmModel(std::size_t index=0) const;
0242
0243
0244 inline G4VEmModel* GetModelByIndex(std::size_t idx = 0, G4bool ver = false) const;
0245
0246
0247 inline void SetFluctModel(G4VEmFluctuationModel*);
0248
0249
0250 inline G4VEmFluctuationModel* FluctModel() const;
0251
0252
0253
0254
0255
0256 protected:
0257 inline void SetParticle(const G4ParticleDefinition* p);
0258 inline void SetSecondaryParticle(const G4ParticleDefinition* p);
0259
0260 public:
0261 inline void SetBaseParticle(const G4ParticleDefinition* p);
0262 inline const G4ParticleDefinition* Particle() const;
0263 inline const G4ParticleDefinition* BaseParticle() const;
0264 inline const G4ParticleDefinition* SecondaryParticle() const;
0265
0266
0267 G4VEnergyLossProcess(G4VEnergyLossProcess &) = delete;
0268 G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right) = delete;
0269
0270
0271
0272
0273
0274
0275 void ActivateSubCutoff(const G4Region* region);
0276
0277
0278 void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
0279
0280 void ActivateForcedInteraction(G4double length,
0281 const G4String& region,
0282 G4bool flag = true);
0283
0284 void ActivateSecondaryBiasing(const G4String& region, G4double factor,
0285 G4double energyLimit);
0286
0287 inline void SetLossFluctuations(G4bool val);
0288
0289 inline void SetSpline(G4bool val);
0290 inline void SetCrossSectionType(G4CrossSectionType val);
0291 inline G4CrossSectionType CrossSectionType() const;
0292
0293
0294 void SetIonisation(G4bool val);
0295 inline G4bool IsIonisationProcess() const;
0296
0297
0298 void SetLinearLossLimit(G4double val);
0299 void SetStepFunction(G4double v1, G4double v2);
0300 void SetLowestEnergyLimit(G4double);
0301
0302 inline G4int NumberOfSubCutoffRegions() const;
0303
0304
0305
0306
0307
0308 void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
0309 void SetCSDARangeTable(G4PhysicsTable* pRange);
0310 void SetRangeTableForLoss(G4PhysicsTable* p);
0311 void SetInverseRangeTable(G4PhysicsTable* p);
0312 void SetLambdaTable(G4PhysicsTable* p);
0313
0314 void SetTwoPeaksXS(std::vector<G4TwoPeaksXS*>*);
0315 void SetEnergyOfCrossSectionMax(std::vector<G4double>*);
0316
0317
0318
0319
0320
0321
0322 void SetDEDXBinning(G4int nbins);
0323
0324
0325 void SetMinKinEnergy(G4double e);
0326 inline G4double MinKinEnergy() const;
0327
0328
0329 void SetMaxKinEnergy(G4double e);
0330 inline G4double MaxKinEnergy() const;
0331
0332
0333 inline G4double CrossSectionBiasingFactor() const;
0334
0335
0336 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*);
0337 inline G4double GetCSDADEDX(G4double kineticEnergy,
0338 const G4MaterialCutsCouple*);
0339 inline G4double GetDEDX(G4double kineticEnergy, const G4MaterialCutsCouple*,
0340 G4double logKineticEnergy);
0341 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*);
0342 inline G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple*,
0343 G4double logKineticEnergy);
0344 inline G4double GetCSDARange(G4double kineticEnergy,
0345 const G4MaterialCutsCouple*);
0346 inline G4double GetKineticEnergy(G4double range,
0347 const G4MaterialCutsCouple*);
0348 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*);
0349 inline G4double GetLambda(G4double kineticEnergy,const G4MaterialCutsCouple*,
0350 G4double logKineticEnergy);
0351
0352 inline G4bool TablesAreBuilt() const;
0353
0354
0355 inline G4PhysicsTable* DEDXTable() const;
0356 inline G4PhysicsTable* DEDXunRestrictedTable() const;
0357 inline G4PhysicsTable* IonisationTable() const;
0358 inline G4PhysicsTable* CSDARangeTable() const;
0359 inline G4PhysicsTable* RangeTableForLoss() const;
0360 inline G4PhysicsTable* InverseRangeTable() const;
0361 inline G4PhysicsTable* LambdaTable() const;
0362 inline std::vector<G4TwoPeaksXS*>* TwoPeaksXS() const;
0363 inline std::vector<G4double>* EnergyOfCrossSectionMax() const;
0364
0365 inline G4bool UseBaseMaterial() const;
0366
0367
0368
0369
0370
0371
0372 const G4Element* GetCurrentElement() const;
0373
0374
0375 void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
0376
0377 private:
0378
0379 void FillSecondariesAlongStep(G4double weight);
0380
0381 void PrintWarning(const G4String&, G4double val) const;
0382
0383
0384 inline void DefineMaterial(const G4MaterialCutsCouple* couple);
0385
0386
0387
0388
0389 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE);
0390 inline G4double GetDEDXForScaledEnergy(G4double scaledKinE,
0391 G4double logScaledKinE);
0392 inline G4double GetIonisationForScaledEnergy(G4double scaledKinE);
0393 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinE);
0394 inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinE,
0395 G4double logScaledKinE);
0396
0397 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinE);
0398 inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinE,
0399 G4double logScaledKinE);
0400
0401 inline G4double ScaledKinEnergyForLoss(G4double range);
0402 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE);
0403 inline G4double GetLambdaForScaledEnergy(G4double scaledKinE,
0404 G4double logScaledKinE);
0405
0406 inline G4double LogScaledEkin(const G4Track& aTrack);
0407
0408 void ComputeLambdaForScaledEnergy(G4double scaledKinE,
0409 const G4Track& aTrack);
0410
0411 G4bool IsRegionForCubcutProcessor(const G4Track& aTrack);
0412
0413 protected:
0414
0415 G4ParticleChangeForLoss fParticleChange;
0416 const G4Material* currentMaterial = nullptr;
0417 const G4MaterialCutsCouple* currentCouple = nullptr;
0418
0419 private:
0420
0421 G4LossTableManager* lManager;
0422 G4EmModelManager* modelManager;
0423 G4VEmModel* currentModel = nullptr;
0424 G4EmBiasingManager* biasManager = nullptr;
0425 G4SafetyHelper* safetyHelper;
0426 G4EmParameters* theParameters;
0427 G4VEmFluctuationModel* fluctModel = nullptr;
0428 G4VAtomDeexcitation* atomDeexcitation = nullptr;
0429 G4VSubCutProducer* subcutProducer = nullptr;
0430
0431 const G4ParticleDefinition* particle = nullptr;
0432 const G4ParticleDefinition* baseParticle = nullptr;
0433 const G4ParticleDefinition* secondaryParticle = nullptr;
0434 G4EmDataHandler* theData = nullptr;
0435
0436 G4PhysicsTable* theDEDXTable = nullptr;
0437 G4PhysicsTable* theDEDXunRestrictedTable = nullptr;
0438 G4PhysicsTable* theIonisationTable = nullptr;
0439 G4PhysicsTable* theRangeTableForLoss = nullptr;
0440 G4PhysicsTable* theCSDARangeTable = nullptr;
0441 G4PhysicsTable* theInverseRangeTable = nullptr;
0442 G4PhysicsTable* theLambdaTable = nullptr;
0443
0444 std::vector<const G4Region*>* scoffRegions = nullptr;
0445 std::vector<G4VEmModel*>* emModels = nullptr;
0446 const std::vector<G4int>* theDensityIdx = nullptr;
0447 const std::vector<G4double>* theDensityFactor = nullptr;
0448 const G4DataVector* theCuts = nullptr;
0449
0450 std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
0451 std::vector<G4TwoPeaksXS*>* fXSpeaks = nullptr;
0452
0453 G4double lowestKinEnergy;
0454 G4double minKinEnergy;
0455 G4double maxKinEnergy;
0456 G4double maxKinEnergyCSDA;
0457
0458 G4double linLossLimit = 0.01;
0459 G4double dRoverRange = 0.2;
0460 G4double finalRange;
0461 G4double lambdaFactor = 0.8;
0462 G4double invLambdaFactor;
0463 G4double biasFactor = 1.0;
0464
0465 G4double massRatio = 1.0;
0466 G4double logMassRatio = 0.0;
0467 G4double fFactor = 1.0;
0468 G4double reduceFactor = 1.0;
0469 G4double chargeSqRatio = 1.0;
0470 G4double fRange = 0.0;
0471 G4double fRangeEnergy = 0.0;
0472
0473 protected:
0474
0475 G4double preStepLambda = 0.0;
0476 G4double preStepKinEnergy = 0.0;
0477 G4double preStepScaledEnergy = 0.0;
0478 G4double mfpKinEnergy = 0.0;
0479
0480 std::size_t currentCoupleIndex = 0;
0481
0482 private:
0483
0484 G4int nBins;
0485 G4int nBinsCSDA;
0486 G4int numberOfModels = 0;
0487 G4int nSCoffRegions = 0;
0488 G4int secID = _DeltaElectron;
0489 G4int tripletID = _TripletElectron;
0490 G4int biasID = _DeltaEBelowCut;
0491 G4int epixeID = _ePIXE;
0492 G4int gpixeID = _GammaPIXE;
0493 G4int mainSecondaries = 1;
0494
0495 std::size_t basedCoupleIndex = 0;
0496 std::size_t coupleIdxRange = 0;
0497 std::size_t idxDEDX = 0;
0498 std::size_t idxDEDXunRestricted = 0;
0499 std::size_t idxIonisation = 0;
0500 std::size_t idxRange = 0;
0501 std::size_t idxCSDA = 0;
0502 std::size_t idxSecRange = 0;
0503 std::size_t idxInverseRange = 0;
0504 std::size_t idxLambda = 0;
0505
0506 G4GPILSelection aGPILSelection;
0507 G4CrossSectionType fXSType = fEmOnePeak;
0508
0509 G4bool lossFluctuationFlag = true;
0510 G4bool useCutAsFinalRange = false;
0511 G4bool tablesAreBuilt = false;
0512 G4bool spline = true;
0513 G4bool isIon = false;
0514 G4bool isIonisation = false;
0515 G4bool useDeexcitation = false;
0516 G4bool biasFlag = false;
0517 G4bool weightFlag = false;
0518 G4bool isMaster = false;
0519 G4bool baseMat = false;
0520 G4bool actLinLossLimit = false;
0521 G4bool actLossFluc = false;
0522 G4bool actBinning = false;
0523 G4bool actMinKinEnergy = false;
0524 G4bool actMaxKinEnergy = false;
0525
0526 std::vector<G4DynamicParticle*> secParticles;
0527 std::vector<G4Track*> scTracks;
0528 };
0529
0530
0531
0532 inline std::size_t G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex() const
0533 {
0534 return currentCoupleIndex;
0535 }
0536
0537
0538
0539 inline void G4VEnergyLossProcess::SelectModel(G4double kinEnergy)
0540 {
0541 currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
0542 currentModel->SetCurrentCouple(currentCouple);
0543 }
0544
0545
0546
0547 inline G4VEmModel* G4VEnergyLossProcess::SelectModelForMaterial(
0548 G4double kinEnergy, std::size_t& idx) const
0549 {
0550 return modelManager->SelectModel(kinEnergy, idx);
0551 }
0552
0553
0554
0555 inline void
0556 G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
0557 {
0558 if(couple != currentCouple) {
0559 currentCouple = couple;
0560 currentMaterial = couple->GetMaterial();
0561 basedCoupleIndex = currentCoupleIndex = couple->GetIndex();
0562 fFactor = chargeSqRatio*biasFactor;
0563 mfpKinEnergy = DBL_MAX;
0564 idxLambda = 0;
0565 if(baseMat) {
0566 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
0567 fFactor *= (*theDensityFactor)[currentCoupleIndex];
0568 }
0569 reduceFactor = 1.0/(fFactor*massRatio);
0570 }
0571 }
0572
0573
0574
0575 inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
0576 {
0577
0578
0579
0580
0581
0582 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e, idxDEDX);
0583 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
0584 return x;
0585 }
0586
0587
0588
0589 inline
0590 G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e, G4double loge)
0591 {
0592
0593
0594
0595
0596
0597 G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->LogVectorValue(e,loge);
0598 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
0599 return x;
0600 }
0601
0602
0603
0604 inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
0605 {
0606 G4double x =
0607 fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e, idxIonisation);
0608 if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
0609 return x;
0610 }
0611
0612
0613
0614 inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
0615 {
0616
0617
0618
0619 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
0620 coupleIdxRange = currentCoupleIndex;
0621 fRangeEnergy = e;
0622 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->Value(e, idxRange);
0623 if (fRange < 0.0) { fRange = 0.0; }
0624 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
0625 }
0626
0627
0628
0629 return fRange;
0630 }
0631
0632 inline G4double
0633 G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e, G4double loge)
0634 {
0635
0636
0637
0638 if(currentCoupleIndex != coupleIdxRange || fRangeEnergy != e) {
0639 coupleIdxRange = currentCoupleIndex;
0640 fRangeEnergy = e;
0641 fRange = reduceFactor*((*theRangeTableForLoss)[basedCoupleIndex])->LogVectorValue(e, loge);
0642 if (fRange < 0.0) { fRange = 0.0; }
0643 else if (e < minKinEnergy) { fRange *= std::sqrt(e/minKinEnergy); }
0644 }
0645
0646
0647
0648 return fRange;
0649 }
0650
0651
0652
0653 inline G4double
0654 G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
0655 {
0656 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e, idxCSDA);
0657 if (x < 0.0) { x = 0.0; }
0658 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
0659 return x;
0660 }
0661
0662
0663
0664 inline G4double
0665 G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e,
0666 G4double loge)
0667 {
0668 G4double x = ((*theCSDARangeTable)[basedCoupleIndex])->LogVectorValue(e, loge);
0669 if (x < 0.0) { x = 0.0; }
0670 else if (e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
0671 return x;
0672 }
0673
0674
0675
0676 inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
0677 {
0678
0679
0680
0681 G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
0682 G4double rmin = v->Energy(0);
0683 G4double e = 0.0;
0684 if(r >= rmin) { e = v->Value(r, idxInverseRange); }
0685 else if(r > 0.0) {
0686 G4double x = r/rmin;
0687 e = minKinEnergy*x*x;
0688 }
0689 return e;
0690 }
0691
0692
0693
0694 inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
0695 {
0696 return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e, idxLambda);
0697 }
0698
0699
0700
0701 inline G4double
0702 G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e, G4double loge)
0703 {
0704 return fFactor*((*theLambdaTable)[basedCoupleIndex])->LogVectorValue(e, loge);
0705 }
0706
0707
0708
0709 inline G4double G4VEnergyLossProcess::LogScaledEkin(const G4Track& track)
0710 {
0711 return track.GetDynamicParticle()->GetLogKineticEnergy() + logMassRatio;
0712 }
0713
0714
0715
0716 inline G4double
0717 G4VEnergyLossProcess::GetDEDX(G4double kinEnergy,
0718 const G4MaterialCutsCouple* couple)
0719 {
0720 DefineMaterial(couple);
0721 return GetDEDXForScaledEnergy(kinEnergy*massRatio);
0722 }
0723
0724
0725
0726 inline G4double
0727 G4VEnergyLossProcess::GetDEDX(G4double kinEnergy,
0728 const G4MaterialCutsCouple* couple,
0729 G4double logKinEnergy)
0730 {
0731 DefineMaterial(couple);
0732 return GetDEDXForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
0733 }
0734
0735
0736
0737 inline G4double
0738 G4VEnergyLossProcess::GetRange(G4double kinEnergy,
0739 const G4MaterialCutsCouple* couple)
0740 {
0741 DefineMaterial(couple);
0742 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio);
0743 }
0744
0745
0746
0747 inline G4double
0748 G4VEnergyLossProcess::GetRange(G4double kinEnergy,
0749 const G4MaterialCutsCouple* couple,
0750 G4double logKinEnergy)
0751 {
0752 DefineMaterial(couple);
0753 return GetScaledRangeForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio);
0754 }
0755
0756
0757
0758 inline G4double
0759 G4VEnergyLossProcess::GetCSDARange(G4double kineticEnergy,
0760 const G4MaterialCutsCouple* couple)
0761 {
0762 DefineMaterial(couple);
0763 return (nullptr == theCSDARangeTable) ? DBL_MAX :
0764 GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
0765 }
0766
0767
0768
0769 inline G4double
0770 G4VEnergyLossProcess::GetKineticEnergy(G4double range,
0771 const G4MaterialCutsCouple* couple)
0772 {
0773 DefineMaterial(couple);
0774 return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
0775 }
0776
0777
0778
0779 inline G4double
0780 G4VEnergyLossProcess::GetLambda(G4double kinEnergy,
0781 const G4MaterialCutsCouple* couple)
0782 {
0783 DefineMaterial(couple);
0784 return (nullptr != theLambdaTable) ?
0785 GetLambdaForScaledEnergy(kinEnergy*massRatio) : 0.0;
0786 }
0787
0788
0789
0790 inline G4double
0791 G4VEnergyLossProcess::GetLambda(G4double kinEnergy,
0792 const G4MaterialCutsCouple* couple,
0793 G4double logKinEnergy)
0794 {
0795 DefineMaterial(couple);
0796 return (nullptr != theLambdaTable) ?
0797 GetLambdaForScaledEnergy(kinEnergy*massRatio, logKinEnergy+logMassRatio)
0798 : 0.0;
0799 }
0800
0801
0802
0803 inline void G4VEnergyLossProcess::SetFluctModel(G4VEmFluctuationModel* p)
0804 {
0805 fluctModel = p;
0806 }
0807
0808
0809
0810 inline G4VEmFluctuationModel* G4VEnergyLossProcess::FluctModel() const
0811 {
0812 return fluctModel;
0813 }
0814
0815
0816
0817 inline void G4VEnergyLossProcess::SetParticle(const G4ParticleDefinition* p)
0818 {
0819 particle = p;
0820 }
0821
0822
0823
0824 inline void
0825 G4VEnergyLossProcess::SetSecondaryParticle(const G4ParticleDefinition* p)
0826 {
0827 secondaryParticle = p;
0828 }
0829
0830
0831
0832 inline void
0833 G4VEnergyLossProcess::SetBaseParticle(const G4ParticleDefinition* p)
0834 {
0835 baseParticle = p;
0836 }
0837
0838
0839
0840 inline const G4ParticleDefinition* G4VEnergyLossProcess::Particle() const
0841 {
0842 return particle;
0843 }
0844
0845
0846
0847 inline const G4ParticleDefinition* G4VEnergyLossProcess::BaseParticle() const
0848 {
0849 return baseParticle;
0850 }
0851
0852
0853
0854 inline const G4ParticleDefinition*
0855 G4VEnergyLossProcess::SecondaryParticle() const
0856 {
0857 return secondaryParticle;
0858 }
0859
0860
0861
0862 inline void G4VEnergyLossProcess::SetLossFluctuations(G4bool val)
0863 {
0864 lossFluctuationFlag = val;
0865 actLossFluc = true;
0866 }
0867
0868
0869
0870 inline void G4VEnergyLossProcess::SetSpline(G4bool val)
0871 {
0872 spline = val;
0873 }
0874
0875
0876
0877 inline void G4VEnergyLossProcess::SetCrossSectionType(G4CrossSectionType val)
0878 {
0879 fXSType = val;
0880 }
0881
0882
0883
0884 inline G4CrossSectionType G4VEnergyLossProcess::CrossSectionType() const
0885 {
0886 return fXSType;
0887 }
0888
0889
0890
0891 inline G4bool G4VEnergyLossProcess::IsIonisationProcess() const
0892 {
0893 return isIonisation;
0894 }
0895
0896
0897
0898 inline G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions() const
0899 {
0900 return nSCoffRegions;
0901 }
0902
0903
0904
0905 inline G4double G4VEnergyLossProcess::MinKinEnergy() const
0906 {
0907 return minKinEnergy;
0908 }
0909
0910
0911
0912 inline G4double G4VEnergyLossProcess::MaxKinEnergy() const
0913 {
0914 return maxKinEnergy;
0915 }
0916
0917
0918
0919 inline G4double G4VEnergyLossProcess::CrossSectionBiasingFactor() const
0920 {
0921 return biasFactor;
0922 }
0923
0924
0925
0926 inline G4bool G4VEnergyLossProcess::TablesAreBuilt() const
0927 {
0928 return tablesAreBuilt;
0929 }
0930
0931
0932
0933 inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTable() const
0934 {
0935 return theDEDXTable;
0936 }
0937
0938
0939
0940 inline G4PhysicsTable* G4VEnergyLossProcess::DEDXunRestrictedTable() const
0941 {
0942 return theDEDXunRestrictedTable;
0943 }
0944
0945
0946
0947 inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTable() const
0948 {
0949 return theIonisationTable;
0950 }
0951
0952
0953
0954 inline G4PhysicsTable* G4VEnergyLossProcess::CSDARangeTable() const
0955 {
0956 return theCSDARangeTable;
0957 }
0958
0959
0960
0961 inline G4PhysicsTable* G4VEnergyLossProcess::RangeTableForLoss() const
0962 {
0963 return theRangeTableForLoss;
0964 }
0965
0966
0967
0968 inline G4PhysicsTable* G4VEnergyLossProcess::InverseRangeTable() const
0969 {
0970 return theInverseRangeTable;
0971 }
0972
0973
0974
0975 inline G4PhysicsTable* G4VEnergyLossProcess::LambdaTable() const
0976 {
0977 return theLambdaTable;
0978 }
0979
0980
0981
0982 inline G4bool G4VEnergyLossProcess::UseBaseMaterial() const
0983 {
0984 return baseMat;
0985 }
0986
0987
0988
0989 inline std::vector<G4double>*
0990 G4VEnergyLossProcess::EnergyOfCrossSectionMax() const
0991 {
0992 return theEnergyOfCrossSectionMax;
0993 }
0994
0995
0996
0997 inline std::vector<G4TwoPeaksXS*>* G4VEnergyLossProcess::TwoPeaksXS() const
0998 {
0999 return fXSpeaks;
1000 }
1001
1002
1003
1004 inline std::size_t G4VEnergyLossProcess::NumberOfModels() const
1005 {
1006 return numberOfModels;
1007 }
1008
1009
1010
1011 inline G4VEmModel* G4VEnergyLossProcess::EmModel(std::size_t index) const
1012 {
1013 return (index < emModels->size()) ? (*emModels)[index] : nullptr;
1014 }
1015
1016
1017
1018 inline G4VEmModel*
1019 G4VEnergyLossProcess::GetModelByIndex(std::size_t idx, G4bool ver) const
1020 {
1021 return modelManager->GetModel((G4int)idx, ver);
1022 }
1023
1024
1025
1026 #endif