Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4IonisParamMat.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// class description
27//
28// The class contains few (physical) quantities related to the Ionisation
29// process, for a material defined by its pointer G4Material*
30//
31// 09-07-98: data moved from G4Material (mma)
32// 09-03-01: copy constructor and assignement operator in public (mma)
33// 28-10-02: add setMeanExcitationEnergy (V.Ivanchenko)
34// 27-09-07: add computation of parameters for ions (V.Ivanchenko)
35// 04-03-08: add fBirks constant (mma)
36// 16-01-19, add exact computation of the density effect (M. Strait)
37
38#ifndef G4IonisParamMat_HH
39#define G4IonisParamMat_HH
40
42#include "G4ios.hh"
43#include "globals.hh"
44
45class G4Material;
47
49{
50 public:
55
56 // parameters for mean energy loss calculation:
57 inline G4double GetMeanExcitationEnergy() const { return fMeanExcitationEnergy; };
58
61
62 inline G4double GetLogMeanExcEnergy() const { return fLogMeanExcEnergy; };
63 inline G4double* GetShellCorrectionVector() const { return fShellCorrectionVector; };
64 inline G4double GetTaul() const { return fTaul; };
65
66 // parameters of the density correction:
67 inline G4double GetPlasmaEnergy() const { return fPlasmaEnergy; };
68 inline G4double GetAdjustmentFactor() const { return fAdjustmentFactor; };
69 inline G4double GetCdensity() const { return fCdensity; };
70 inline G4double GetMdensity() const { return fMdensity; };
71 inline G4double GetAdensity() const { return fAdensity; };
72 inline G4double GetX0density() const { return fX0density; };
73 inline G4double GetX1density() const { return fX1density; };
74 inline G4double GetD0density() const { return fD0density; };
75
76 // user defined density correction parameterisation
78 G4double cd, G4double md, G4double ad, G4double x0, G4double x1, G4double d0);
79
80 // defined density correction parameterisation via base material
81 void SetDensityEffectParameters(const G4Material* bmat);
82
84
86 {
87 return fDensityEffectCalc;
88 }
89
90 // compute density correction as a function of the kinematic variable
91 // x = log10(beta*gamma) using parameterisation of calculator
93 {
94 return (nullptr == fDensityEffectCalc) ? GetDensityCorrection(x)
95 : fDensityEffectCalc->ComputeDensityCorrection(x);
96 }
97
98 // use parameterisation
100
102
103 // parameters of the energy loss fluctuation model:
104 inline G4double GetF1fluct() const { return fF1fluct; };
105 inline G4double GetF2fluct() const { return fF2fluct; };
106 inline G4double GetEnergy1fluct() const { return fEnergy1fluct; };
107 inline G4double GetLogEnergy1fluct() const { return fLogEnergy1fluct; };
108 inline G4double GetEnergy2fluct() const { return fEnergy2fluct; };
109 inline G4double GetLogEnergy2fluct() const { return fLogEnergy2fluct; };
110 inline G4double GetEnergy0fluct() const { return fEnergy0fluct; };
111 inline G4double GetRateionexcfluct() const { return fRateionexcfluct; };
112
113 // parameters for ion corrections computations
114 inline G4double GetZeffective() const { return fZeff; };
115 inline G4double GetFermiEnergy() const { return fFermiEnergy; };
116 inline G4double GetLFactor() const { return fLfactor; };
117 inline G4double GetInvA23() const { return fInvA23; };
118
119 // parameters for Birks attenuation:
120 inline void SetBirksConstant(G4double value) { fBirks = value; };
121 inline G4double GetBirksConstant() const { return fBirks; };
122
123 // parameters for average energy per ion
124 inline void SetMeanEnergyPerIonPair(G4double value) { fMeanEnergyPerIon = value; };
125 inline G4double GetMeanEnergyPerIonPair() const { return fMeanEnergyPerIon; };
126
127 // parameter for sampling of positron annihilation at rest
128 inline void SetOrtoPositroniumFraction(G4double value) { fOrtoPositroniumFraction = value; };
129 inline G4double GetOrtoPositroniumFraction() const { return fOrtoPositroniumFraction; };
130
131
132 // operators
133 G4bool operator==(const G4IonisParamMat&) const = delete;
134 G4bool operator!=(const G4IonisParamMat&) const = delete;
135
136 private:
137 // Compute mean parameters : ExcitationEnergy,Shell corretion vector ...
138 void ComputeMeanParameters();
139
140 // Compute parameters for the density effect
141 void ComputeDensityEffectParameters();
142
143 // Compute parameters for the energy fluctuation model
144 void ComputeFluctModel();
145
146 // Compute parameters for ion parameterizations
147 void ComputeIonParameters();
148
149 //
150 // data members
151 //
152 const G4Material* fMaterial; // this material
153
154 G4DensityEffectCalculator* fDensityEffectCalc; // online calculator
155 G4double* fShellCorrectionVector; // shell correction coefficients
156
157 // parameters for mean energy loss calculation
158 G4double fMeanExcitationEnergy; //
159 G4double fLogMeanExcEnergy; //
160 G4double fTaul; // lower limit of Bethe-Bloch formula
161
162 // parameters of the density correction
163 G4double fCdensity; // mat.constant
164 G4double fMdensity; // exponent
165 G4double fAdensity; //
166 G4double fX0density; //
167 G4double fX1density; //
168 G4double fD0density;
169
170 G4double fPlasmaEnergy;
171 G4double fAdjustmentFactor;
172
173 // parameters of the energy loss fluctuation model
174 G4double fF1fluct;
175 G4double fF2fluct;
176 G4double fEnergy1fluct;
177 G4double fLogEnergy1fluct;
178 G4double fEnergy2fluct;
179 G4double fLogEnergy2fluct;
180 G4double fEnergy0fluct;
181 G4double fRateionexcfluct;
182
183 // parameters for ion corrections computations
184 G4double fZeff;
185 G4double fFermiEnergy;
186 G4double fLfactor;
187 G4double fInvA23;
188
189 // parameter for Birks attenuation
190 G4double fBirks;
191 // average energy per ion pair
192 G4double fMeanEnergyPerIon;
193 G4double twoln10;
194 // parameter for sampling of positron annihilation at rest
195 G4double fOrtoPositroniumFraction{0.035};
196
197 // static data created only once
198 static G4DensityEffectData* fDensityData;
199};
200
201#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4double GetF2fluct() const
G4double GetLogEnergy1fluct() const
G4double GetOrtoPositroniumFraction() const
G4bool operator==(const G4IonisParamMat &) const =delete
G4IonisParamMat & operator=(const G4IonisParamMat &)=delete
G4double GetMdensity() const
G4double GetX1density() const
G4double * GetShellCorrectionVector() const
G4double GetTaul() const
G4DensityEffectCalculator * GetDensityEffectCalculator() const
G4double GetAdjustmentFactor() const
G4double GetX0density() const
G4double DensityCorrection(G4double x) const
G4double GetCdensity() const
void SetOrtoPositroniumFraction(G4double value)
G4double GetLogEnergy2fluct() const
void SetDensityEffectParameters(G4double cd, G4double md, G4double ad, G4double x0, G4double x1, G4double d0)
G4double GetMeanExcitationEnergy() const
G4double GetF1fluct() const
G4double GetMeanEnergyPerIonPair() const
static G4DensityEffectData * GetDensityEffectData()
G4double FindMeanExcitationEnergy(const G4Material *) const
void ComputeDensityEffectOnFly(G4bool)
G4double GetInvA23() const
G4double GetD0density() const
G4double GetEnergy2fluct() const
G4double GetZeffective() const
G4IonisParamMat(const G4IonisParamMat &)=delete
G4IonisParamMat(const G4Material *)
G4double GetAdensity() const
G4double GetEnergy1fluct() const
G4double GetFermiEnergy() const
G4double GetDensityCorrection(G4double x) const
G4bool operator!=(const G4IonisParamMat &) const =delete
void SetMeanExcitationEnergy(G4double value)
G4double GetPlasmaEnergy() const
G4double GetLFactor() const
void SetBirksConstant(G4double value)
void SetMeanEnergyPerIonPair(G4double value)
G4double GetRateionexcfluct() const
G4double GetEnergy0fluct() const
G4double GetLogMeanExcEnergy() const
G4double GetBirksConstant() const