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