File indexing completed on 2025-01-18 09:58:09
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 #ifndef G4DNAUeharaScreenedRutherfordElasticModel_h
0027 #define G4DNAUeharaScreenedRutherfordElasticModel_h 1
0028
0029 #include <CLHEP/Units/SystemOfUnits.h>
0030
0031 #include "G4VEmModel.hh"
0032 #include "G4ParticleChangeForGamma.hh"
0033 #include "G4ProductionCutsTable.hh"
0034 #include "G4NistManager.hh"
0035
0036 class G4DNAUeharaScreenedRutherfordElasticModel : public G4VEmModel
0037 {
0038 public:
0039 G4DNAUeharaScreenedRutherfordElasticModel(const G4ParticleDefinition* p = nullptr,
0040 const G4String& nam = "DNAUeharaScreenedRutherfordElasticModel");
0041
0042 ~G4DNAUeharaScreenedRutherfordElasticModel() override = default;
0043
0044 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0045
0046 G4double CrossSectionPerVolume(const G4Material* material,
0047 const G4ParticleDefinition* p,
0048 G4double ekin,
0049 G4double emin,
0050 G4double emax) override;
0051
0052 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0053 const G4MaterialCutsCouple*,
0054 const G4DynamicParticle*,
0055 G4double tmin,
0056 G4double maxEnergy) override;
0057
0058 inline void SelectFasterComputation(G4bool input);
0059
0060
0061
0062 inline void SetKillBelowThreshold (G4double threshold);
0063 inline G4double GetKillBelowThreshold();
0064 inline void SelectHighEnergyLimit(G4double threshold);
0065
0066
0067
0068 G4DNAUeharaScreenedRutherfordElasticModel&
0069 operator=(const G4DNAUeharaScreenedRutherfordElasticModel &right) = delete;
0070 G4DNAUeharaScreenedRutherfordElasticModel(
0071 const G4DNAUeharaScreenedRutherfordElasticModel&) = delete;
0072
0073 private:
0074
0075
0076 G4double RutherfordCrossSection(G4double energy, G4double z);
0077 G4double ScreeningFactor(G4double energy, G4double z);
0078
0079
0080 G4double BrennerZaiderRandomizeCosTheta(G4double k);
0081 G4double CalculatePolynomial(G4double k,
0082 std::vector<G4double>& vec);
0083
0084
0085 G4double ScreenedRutherfordRandomizeCosTheta(G4double k,
0086 G4double z);
0087
0088 protected:
0089 G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr;
0090
0091 private:
0092 G4double iLowEnergyLimit;
0093 G4double intermediateEnergyLimit;
0094 G4double iHighEnergyLimit;
0095
0096
0097 std::vector<G4double> betaCoeff;
0098 std::vector<G4double> deltaCoeff;
0099 std::vector<G4double> gamma035_10Coeff;
0100 std::vector<G4double> gamma10_100Coeff;
0101 std::vector<G4double> gamma100_200Coeff;
0102
0103
0104 const std::vector<G4double>* fpWaterDensity = nullptr;
0105
0106 G4int verboseLevel;
0107 G4bool isInitialised = false;
0108
0109
0110 G4bool fasterCode = false;
0111 };
0112
0113
0114 inline void G4DNAUeharaScreenedRutherfordElasticModel::
0115 SelectFasterComputation(G4bool input)
0116 {
0117 fasterCode = input;
0118 }
0119
0120
0121
0122
0123 inline void
0124 G4DNAUeharaScreenedRutherfordElasticModel::SelectHighEnergyLimit(
0125 G4double threshold)
0126 {
0127 if(threshold > 10. * CLHEP::keV)
0128 {
0129 G4Exception (
0130 "*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel class is "
0131 "used above 10 keV !",
0132 "", JustWarning, "");
0133 }
0134
0135 SetHighEnergyLimit(threshold);
0136 }
0137
0138 inline void
0139 G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold(G4double)
0140 {
0141 G4ExceptionDescription errMsg;
0142 errMsg << "*** WARNING : "
0143 << "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold"
0144 << "is deprecated, the kill threshold won't be taken into account";
0145
0146 G4Exception (
0147 "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold",
0148 "DEPRECATED", JustWarning, errMsg);
0149 }
0150
0151 inline G4double
0152 G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold()
0153 {
0154 G4ExceptionDescription errMsg;
0155 errMsg << "*** WARNING : "
0156 << "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold"
0157 << "is deprecated, the returned value is nonsense";
0158
0159 G4Exception (
0160 "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold",
0161 "DEPRECATED", JustWarning, errMsg);
0162
0163 return -1;
0164 }
0165
0166
0167
0168
0169 #endif