File indexing completed on 2025-01-18 09:58:04
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 #ifndef G4ContinuousGainOfEnergy_h
0036 #define G4ContinuousGainOfEnergy_h 1
0037
0038 #include "globals.hh"
0039 #include "G4ProductionCutsTable.hh"
0040 #include "G4VContinuousProcess.hh"
0041
0042 class G4Material;
0043 class G4MaterialCutsCouple;
0044 class G4ParticleChange;
0045 class G4ParticleDefinition;
0046 class G4Step;
0047 class G4Track;
0048 class G4VEmModel;
0049 class G4VEnergyLossProcess;
0050
0051 class G4ContinuousGainOfEnergy : public G4VContinuousProcess
0052 {
0053 public:
0054 explicit G4ContinuousGainOfEnergy(const G4String& name = "EnergyGain",
0055 G4ProcessType type = fElectromagnetic);
0056
0057 ~G4ContinuousGainOfEnergy() override;
0058
0059 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
0060
0061 void SetLossFluctuations(G4bool val);
0062
0063 inline void SetDirectEnergyLossProcess(G4VEnergyLossProcess* aProcess)
0064 {
0065 fDirectEnergyLossProcess = aProcess;
0066 };
0067
0068 void SetDirectParticle(G4ParticleDefinition* p);
0069
0070 void ProcessDescription(std::ostream&) const override;
0071 void DumpInfo() const override { ProcessDescription(G4cout); };
0072
0073 G4ContinuousGainOfEnergy(G4ContinuousGainOfEnergy&) = delete;
0074 G4ContinuousGainOfEnergy& operator=(const G4ContinuousGainOfEnergy& right) =
0075 delete;
0076
0077 protected:
0078 G4double GetContinuousStepLimit(const G4Track& track,
0079 G4double previousStepSize,
0080 G4double currentMinimumStep,
0081 G4double& currentSafety) override;
0082
0083 private:
0084 void DefineMaterial(const G4MaterialCutsCouple* couple);
0085 void SetDynamicMassCharge(const G4Track& track, G4double energy);
0086
0087 const G4Material* fCurrentMaterial = nullptr;
0088 const G4MaterialCutsCouple* fCurrentCouple = nullptr;
0089
0090 G4VEmModel* fCurrentModel = nullptr;
0091 G4VEnergyLossProcess* fDirectEnergyLossProcess = nullptr;
0092 G4ParticleDefinition* fDirectPartDef = nullptr;
0093
0094 G4double fCurrentTcut = 0.;
0095 G4double fPreStepKinEnergy = 1.;
0096 G4double fLinLossLimit = 0.05;
0097 G4double fMassRatio = 1.;
0098
0099 size_t fCurrentCoupleIndex = 9999999;
0100
0101 G4bool fIsIon = false;
0102 G4bool fLossFluctuationFlag = true;
0103 G4bool fLossFluctuationArePossible = true;
0104 };
0105
0106
0107 inline void G4ContinuousGainOfEnergy::DefineMaterial(
0108 const G4MaterialCutsCouple* couple)
0109 {
0110 if(couple != fCurrentCouple)
0111 {
0112 fCurrentCouple = couple;
0113 fCurrentMaterial = couple->GetMaterial();
0114 fCurrentCoupleIndex = couple->GetIndex();
0115
0116 const std::vector<G4double>* aVec =
0117 G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(1);
0118 fCurrentTcut = (*aVec)[fCurrentCoupleIndex];
0119 }
0120 }
0121
0122 #endif