File indexing completed on 2025-01-18 09:58:10
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
0055
0056
0057
0058
0059
0060
0061
0062
0063
0064 #ifndef G4eCoulombScatteringModel_h
0065 #define G4eCoulombScatteringModel_h 1
0066
0067 #include "G4VEmModel.hh"
0068 #include "globals.hh"
0069 #include "G4MaterialCutsCouple.hh"
0070 #include "G4WentzelOKandVIxSection.hh"
0071
0072 class G4ParticleChangeForGamma;
0073 class G4ParticleDefinition;
0074 class G4IonTable;
0075 class G4NistManager;
0076
0077 class G4eCoulombScatteringModel : public G4VEmModel
0078 {
0079
0080 public:
0081
0082 explicit G4eCoulombScatteringModel(G4bool combined = true);
0083
0084 ~G4eCoulombScatteringModel() override;
0085
0086 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0087
0088 void InitialiseLocal(const G4ParticleDefinition*,
0089 G4VEmModel* masterModel) override;
0090
0091 G4double ComputeCrossSectionPerAtom(
0092 const G4ParticleDefinition*,
0093 G4double kinEnergy,
0094 G4double Z,
0095 G4double A,
0096 G4double cut,
0097 G4double emax) override;
0098
0099 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0100 const G4MaterialCutsCouple*,
0101 const G4DynamicParticle*,
0102 G4double tmin,
0103 G4double maxEnergy) override;
0104
0105 G4double MinPrimaryEnergy(const G4Material*, const G4ParticleDefinition*,
0106 G4double) final;
0107
0108
0109 inline void SetLowEnergyThreshold(G4double val);
0110
0111
0112 inline void SetRecoilThreshold(G4double eth);
0113
0114
0115 inline void SetFixedCut(G4double);
0116
0117
0118 inline G4double GetFixedCut() const;
0119
0120
0121 G4eCoulombScatteringModel & operator=
0122 (const G4eCoulombScatteringModel &right) = delete;
0123 G4eCoulombScatteringModel(const G4eCoulombScatteringModel&) = delete;
0124
0125 protected:
0126
0127 inline void DefineMaterial(const G4MaterialCutsCouple*);
0128
0129 inline void SetupParticle(const G4ParticleDefinition*);
0130
0131 private:
0132
0133 G4IonTable* theIonTable;
0134 G4ParticleChangeForGamma* fParticleChange = nullptr;
0135 G4WentzelOKandVIxSection* wokvi;
0136 G4NistManager* fNistManager;
0137
0138 const std::vector<G4double>* pCuts = nullptr;
0139
0140 const G4ParticleDefinition* particle = nullptr;
0141 const G4ParticleDefinition* theProton;
0142
0143 const G4MaterialCutsCouple* currentCouple = nullptr;
0144 const G4Material* currentMaterial = nullptr;
0145 G4int currentMaterialIndex = 0;
0146
0147 G4double cosThetaMin = 1.0;
0148 G4double cosThetaMax = -1.0;
0149 G4double recoilThreshold = 0.0;
0150 G4double elecRatio = 0.0;
0151 G4double fixedCut = -1.0;
0152 G4double mass;
0153
0154 G4bool isCombined;
0155 };
0156
0157
0158
0159 inline
0160 void G4eCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup)
0161 {
0162 if(cup != currentCouple) {
0163 currentCouple = cup;
0164 currentMaterial = cup->GetMaterial();
0165 currentMaterialIndex = currentCouple->GetIndex();
0166 }
0167 }
0168
0169
0170
0171 inline
0172 void G4eCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
0173 {
0174
0175 if(p != particle) {
0176 particle = p;
0177 mass = particle->GetPDGMass();
0178 wokvi->SetupParticle(p);
0179 }
0180 }
0181
0182
0183
0184 inline void G4eCoulombScatteringModel::SetRecoilThreshold(G4double eth)
0185 {
0186 recoilThreshold = eth;
0187 }
0188
0189
0190
0191 inline void G4eCoulombScatteringModel::SetFixedCut(G4double val)
0192 {
0193 fixedCut = val;
0194 }
0195
0196
0197
0198 inline G4double G4eCoulombScatteringModel::GetFixedCut() const
0199 {
0200 return fixedCut;
0201 }
0202
0203
0204
0205 #endif