Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAUeharaScreenedRutherfordElasticModel.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25
26#ifndef G4DNAUeharaScreenedRutherfordElasticModel_h
27#define G4DNAUeharaScreenedRutherfordElasticModel_h 1
28
30
31#include "G4VEmModel.hh"
34#include "G4NistManager.hh"
35
37{
38public:
40 const G4String& nam = "DNAUeharaScreenedRutherfordElasticModel");
41
43
44 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
45
47 const G4ParticleDefinition* p,
48 G4double ekin,
49 G4double emin,
50 G4double emax) override;
51
52 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
54 const G4DynamicParticle*,
55 G4double tmin,
56 G4double maxEnergy) override;
57
58 inline void SelectFasterComputation(G4bool input);
59
60 //---
61 // kept for backward compatibility
62 inline void SetKillBelowThreshold (G4double threshold);
64 inline void SelectHighEnergyLimit(G4double threshold);
65 //---
66
67 //
72
73private:
74
75 // -- Cross section
76 G4double RutherfordCrossSection(G4double energy, G4double z);
77 G4double ScreeningFactor(G4double energy, G4double z);
78
79 // -- Final state according to Brenner & Zaider
80 G4double BrennerZaiderRandomizeCosTheta(G4double k);
81 G4double CalculatePolynomial(G4double k,
82 std::vector<G4double>& vec);
83
84 // -- Final state according to Screened Rutherford
85 G4double ScreenedRutherfordRandomizeCosTheta(G4double k,
86 G4double z);
87
88protected:
90
91private:
92 G4double iLowEnergyLimit;
93 G4double intermediateEnergyLimit;
94 G4double iHighEnergyLimit;
95
96 // -- Brenner & Zaider
97 std::vector<G4double> betaCoeff;
98 std::vector<G4double> deltaCoeff;
99 std::vector<G4double> gamma035_10Coeff;
100 std::vector<G4double> gamma10_100Coeff;
101 std::vector<G4double> gamma100_200Coeff;
102
103 // -- Water density table
104 const std::vector<G4double>* fpWaterDensity = nullptr;
105
106 G4int verboseLevel;
107 G4bool isInitialised = false;
108 // Selection of computation method
109 // We do not recommend "true" usage with the current cumul. proba. settings
110 G4bool fasterCode = false;
111};
112
113
116{
117 fasterCode = input;
118}
119
120//---
121// kept for backward compatibility
122
123inline void
125 G4double threshold)
126{
127 if(threshold > 10. * CLHEP::keV)
128 {
130 "*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel class is "
131 "used above 10 keV !",
132 "", JustWarning, "");
133 }
134
135 SetHighEnergyLimit(threshold);
136}
137
138inline void
140{
142 errMsg << "*** WARNING : "
143 << "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold"
144 << "is deprecated, the kill threshold won't be taken into account";
145
147 "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold",
148 "DEPRECATED", JustWarning, errMsg);
149}
150
151inline G4double
153{
155 errMsg << "*** WARNING : "
156 << "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold"
157 << "is deprecated, the returned value is nonsense";
158
160 "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold",
161 "DEPRECATED", JustWarning, errMsg);
162
163 return -1;
164}
165//---
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168
169#endif
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
~G4DNAUeharaScreenedRutherfordElasticModel() override=default
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4DNAUeharaScreenedRutherfordElasticModel(const G4DNAUeharaScreenedRutherfordElasticModel &)=delete
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
G4DNAUeharaScreenedRutherfordElasticModel & operator=(const G4DNAUeharaScreenedRutherfordElasticModel &right)=delete
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746