Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAUeharaScreenedRutherfordElasticModel Class Reference

#include <G4DNAUeharaScreenedRutherfordElasticModel.hh>

+ Inheritance diagram for G4DNAUeharaScreenedRutherfordElasticModel:

Public Member Functions

 G4DNAUeharaScreenedRutherfordElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAUeharaScreenedRutherfordElasticModel")
 
virtual ~G4DNAUeharaScreenedRutherfordElasticModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SelectFasterComputation (G4bool input)
 
void SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
void SelectHighEnergyLimit (G4double threshold)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 36 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAUeharaScreenedRutherfordElasticModel()

G4DNAUeharaScreenedRutherfordElasticModel::G4DNAUeharaScreenedRutherfordElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAUeharaScreenedRutherfordElasticModel" 
)

Definition at line 40 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

42 :
43G4VEmModel(nam), isInitialised(false)
44{
45 fpWaterDensity = 0;
46
47 // Switch between two final state models
48 intermediateEnergyLimit = 200. * eV;
49
50 SetLowEnergyLimit(9.*eV);
51
52 SetHighEnergyLimit(1.*MeV);
53
54 verboseLevel = 0;
55 // Verbosity scale:
56 // 0 = nothing
57 // 1 = warning for energy non-conservation
58 // 2 = details of energy budget
59 // 3 = calculation of cross sections, file openings, sampling of atoms
60 // 4 = entering in methods
61
62#ifdef UEHARA_VERBOSE
63 if (verboseLevel)
64 {
65 G4cout << "Screened Rutherford Elastic model is constructed " << G4endl
66 << "Energy range: "
67 << LowEnergyLimit() / eV << " eV - "
68 << HighEnergyLimit() / MeV << " MeV"
69 << G4endl;
70 }
71#endif
72
74
75 // Selection of computation method
76 // We do not recommend "true" usage with the current cumul. proba. settings
77 fasterCode = false;
78}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ ~G4DNAUeharaScreenedRutherfordElasticModel()

G4DNAUeharaScreenedRutherfordElasticModel::~G4DNAUeharaScreenedRutherfordElasticModel ( )
virtual

Definition at line 82 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

84{
85}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAUeharaScreenedRutherfordElasticModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 192 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

198{
199#ifdef UEHARA_VERBOSE
200 if (verboseLevel > 3)
201 {
202 G4cout
203 << "Calling CrossSectionPerVolume() of G4DNAUeharaScreenedRutherfordElasticModel"
204 << G4endl;
205 }
206#endif
207
208 // Calculate total cross section for model
209
210 G4double sigma = 0.;
211 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
212
213 G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
214 G4double n = ScreeningFactor(ekin,z);
215 G4double crossSection = RutherfordCrossSection(ekin, z);
216 sigma = pi * crossSection / (n * (n + 1.));
217
218#ifdef UEHARA_VERBOSE
219 if (verboseLevel > 2)
220 {
221 G4cout << "__________________________________" << G4endl;
222 G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO START"
223 << G4endl;
224 G4cout << "=== Kinetic energy(eV)=" << ekin/eV
225 << " particle : " << particleDefinition->GetParticleName() << G4endl;
226 G4cout << "=== Cross section per water molecule (cm^2)=" << sigma/cm/cm
227 << G4endl;
228 G4cout << "=== Cross section per water molecule (cm^-1)="
229 << sigma*waterDensity/(1./cm) << G4endl;
230 G4cout << "=== G4DNAUeharaScreenedRutherfordElasticModel - XS INFO END"
231 << G4endl;
232 }
233#endif
234
235 return sigma*waterDensity;
236}
double G4double
Definition: G4Types.hh:83
size_t GetIndex() const
Definition: G4Material.hh:258
const G4double pi

◆ GetKillBelowThreshold()

G4double G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 147 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

148{
150 errMsg << "*** WARNING : "
151 << "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold"
152 << "is deprecated, the returned value is nonsense";
153
155 "G4DNAUeharaScreenedRutherfordElasticModel::GetKillBelowThreshold",
156 "DEPRECATED", JustWarning, errMsg);
157
158 return -1;
159}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40

◆ Initialise()

void G4DNAUeharaScreenedRutherfordElasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 90 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

93{
94#ifdef UEHARA_VERBOSE
95 if (verboseLevel > 3)
96 {
97 G4cout << "Calling G4DNAUeharaScreenedRutherfordElasticModel::Initialise()"
98 << G4endl;
99 }
100#endif
101
102 if(particle->GetParticleName() != "e-")
103 {
104 G4Exception("*** WARNING: the G4DNAUeharaScreenedRutherfordElasticModel is "
105 "not intented to be used with another particle than the electron",
106 "",FatalException,"") ;
107 }
108
109 // Energy limits
110
111 if(LowEnergyLimit() < 9.*CLHEP::eV)
112 {
113 G4Exception("*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel "
114 "class is not validated below 9 eV",
115 "",JustWarning,"") ;
116 }
117
118 if (HighEnergyLimit() > 10.*CLHEP::keV)
119 {
120 G4Exception("*** WARNING: the G4DNAUeharaScreenedRutherfordElasticModel "
121 "class is used above 10 keV",
122 "",JustWarning,"") ;
123 }
124
125#ifdef UEHARA_VERBOSE
126 if( verboseLevel>0 )
127 {
128 G4cout << "Screened Rutherford elastic model is initialized " << G4endl
129 << "Energy range: "
130 << LowEnergyLimit() / eV << " eV - "
131 << HighEnergyLimit() / MeV << " MeV"
132 << G4endl;
133 }
134#endif
135
136 if (isInitialised){ return; }
137
138 // Constants for final state by Brenner & Zaider
139 // Note: the instantiation must be placed after if (isInitialised)
140
141 betaCoeff=
142 {
143 7.51525,
144 -0.41912,
145 7.2017E-3,
146 -4.646E-5,
147 1.02897E-7};
148
149 deltaCoeff=
150 {
151 2.9612,
152 -0.26376,
153 4.307E-3,
154 -2.6895E-5,
155 5.83505E-8};
156
157 gamma035_10Coeff=
158 {
159 -1.7013,
160 -1.48284,
161 0.6331,
162 -0.10911,
163 8.358E-3,
164 -2.388E-4};
165
166 gamma10_100Coeff =
167 {
168 -3.32517,
169 0.10996,
170 -4.5255E-3,
171 5.8372E-5,
172 -2.4659E-7};
173
174 gamma100_200Coeff=
175 {
176 2.4775E-2,
177 -2.96264E-5,
178 -1.20655E-7};
179
180 // Initialize water density pointer
181 fpWaterDensity =
183 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
184
186 isInitialised = true;
187}
@ FatalException
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ SampleSecondaries()

void G4DNAUeharaScreenedRutherfordElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 295 of file G4DNAUeharaScreenedRutherfordElasticModel.cc.

301{
302#ifdef UEHARA_VERBOSE
303 if (verboseLevel > 3)
304 {
305 G4cout
306 << "Calling SampleSecondaries() of G4DNAUeharaScreenedRutherfordElasticModel"
307 << G4endl;
308 }
309#endif
310
311 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
312
313 G4double cosTheta = 0.;
314
315 if (electronEnergy0<intermediateEnergyLimit)
316 {
317#ifdef UEHARA_VERBOSE
318 if (verboseLevel > 3)
319 G4cout << "---> Using Brenner & Zaider model" << G4endl;
320#endif
321 cosTheta = BrennerZaiderRandomizeCosTheta(electronEnergy0);
322 }
323 else //if (electronEnergy0>=intermediateEnergyLimit)
324 {
325#ifdef UEHARA_VERBOSE
326 if (verboseLevel > 3)
327 G4cout << "---> Using Screened Rutherford model" << G4endl;
328#endif
329 G4double z = 7.42; // FROM PMB 37 (1992) 1841-1858 p1842
330 cosTheta = ScreenedRutherfordRandomizeCosTheta(electronEnergy0,z);
331 }
332
333 G4double phi = 2. * pi * G4UniformRand();
334
335 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
336 G4ThreeVector xVers = zVers.orthogonal();
337 G4ThreeVector yVers = zVers.cross(xVers);
338
339 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
340 G4double yDir = xDir;
341 xDir *= std::cos(phi);
342 yDir *= std::sin(phi);
343
344 G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
345
347
349}
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)

◆ SelectFasterComputation()

void G4DNAUeharaScreenedRutherfordElasticModel::SelectFasterComputation ( G4bool  input)
inline

Definition at line 109 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

111{
112 fasterCode = input;
113}

Referenced by G4EmDNAPhysics_option5::ConstructProcess().

◆ SelectHighEnergyLimit()

void G4DNAUeharaScreenedRutherfordElasticModel::SelectHighEnergyLimit ( G4double  threshold)
inline

Definition at line 119 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

121{
122 if(threshold > 10. * CLHEP::keV)
123 {
125 "*** WARNING : the G4DNAUeharaScreenedRutherfordElasticModel class is "
126 "used above 10 keV !",
127 "", JustWarning, "");
128 }
129
130 SetHighEnergyLimit(threshold);
131}

◆ SetKillBelowThreshold()

void G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 134 of file G4DNAUeharaScreenedRutherfordElasticModel.hh.

135{
137 errMsg << "*** WARNING : "
138 << "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold"
139 << "is deprecated, the kill threshold won't be taken into account";
140
142 "G4DNAUeharaScreenedRutherfordElasticModel::SetKillBelowThreshold",
143 "DEPRECATED", JustWarning, errMsg);
144}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAUeharaScreenedRutherfordElasticModel::fParticleChangeForGamma
protected

The documentation for this class was generated from the following files: