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
0027
0028 #ifndef G4DNAScreenedRutherfordElasticModel_h
0029 #define G4DNAScreenedRutherfordElasticModel_h 1
0030
0031 #include <CLHEP/Units/SystemOfUnits.h>
0032
0033 #include "G4VEmModel.hh"
0034 #include "G4ParticleChangeForGamma.hh"
0035 #include "G4ProductionCutsTable.hh"
0036 #include "G4NistManager.hh"
0037
0038 class G4DNAScreenedRutherfordElasticModel : public G4VEmModel
0039 {
0040 public:
0041 G4DNAScreenedRutherfordElasticModel(const G4ParticleDefinition* p = nullptr,
0042 const G4String& nam = "DNAScreenedRutherfordElasticModel");
0043
0044 ~G4DNAScreenedRutherfordElasticModel() override;
0045
0046 G4DNAScreenedRutherfordElasticModel& operator=
0047 (const G4DNAScreenedRutherfordElasticModel &right) = delete;
0048 G4DNAScreenedRutherfordElasticModel(const G4DNAScreenedRutherfordElasticModel&) = delete;
0049
0050 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0051
0052 G4double CrossSectionPerVolume(const G4Material* material,
0053 const G4ParticleDefinition* p,
0054 G4double ekin,
0055 G4double emin,
0056 G4double emax) override;
0057
0058 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0059 const G4MaterialCutsCouple*,
0060 const G4DynamicParticle*,
0061 G4double tmin,
0062 G4double maxEnergy) override;
0063
0064 inline void SetKillBelowThreshold (G4double threshold);
0065
0066 inline void SelectFasterComputation(G4bool input);
0067
0068 protected:
0069 G4ParticleChangeForGamma* fParticleChangeForGamma;
0070
0071 private:
0072 G4double lowEnergyLimit;
0073 G4double intermediateEnergyLimit;
0074 G4double highEnergyLimit;
0075
0076
0077 std::vector<G4double> betaCoeff;
0078 std::vector<G4double> deltaCoeff;
0079 std::vector<G4double> gamma035_10Coeff;
0080 std::vector<G4double> gamma10_100Coeff;
0081 std::vector<G4double> gamma100_200Coeff;
0082
0083
0084 const std::vector<G4double>* fpWaterDensity;
0085
0086 G4int verboseLevel;
0087
0088 G4bool isInitialised{false};
0089 G4bool fasterCode;
0090
0091
0092 G4double RutherfordCrossSection(G4double energy, G4double z);
0093 G4double ScreeningFactor(G4double energy, G4double z);
0094
0095
0096 G4double BrennerZaiderRandomizeCosTheta(G4double k);
0097 G4double CalculatePolynomial(G4double k, std::vector<G4double>& vec);
0098
0099
0100
0101 G4double ScreenedRutherfordRandomizeCosTheta(G4double k, G4double z);
0102
0103 };
0104
0105 inline void
0106 G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold (G4double)
0107 {
0108 G4ExceptionDescription errMsg;
0109 errMsg << "The method G4DNAScreenedRutherfordElasticModel::"
0110 "SetKillBelowThreshold is deprecated";
0111
0112 G4Exception("G4DNAScreenedRutherfordElasticModel::SetKillBelowThreshold",
0113 "deprecated",
0114 JustWarning,
0115 errMsg);
0116 }
0117
0118 inline void
0119 G4DNAScreenedRutherfordElasticModel::SelectFasterComputation(G4bool input)
0120 {
0121 fasterCode = input;
0122 }
0123
0124
0125
0126 #endif