Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
27
28// class description
29//
30// The class contains few (physical) quantities related to the Ionisation
31// process, for a material defined by its pointer G4Material*
32//
33
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
35
36// 09-07-98: data moved from G4Material (mma)
37// 09-03-01: copy constructor and assignement operator in public (mma)
38// 28-10-02: add setMeanExcitationEnergy (V.Ivanchenko)
39// 27-09-07: add computation of parameters for ions (V.Ivanchenko)
40// 04-03-08: add fBirks constant (mma)
41// 16-01-19, add exact computation of the density effect (M. Strait)
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
44
45#ifndef G4IonisParamMat_HH
46#define G4IonisParamMat_HH
47
48#include "G4ios.hh"
49#include "globals.hh"
50#include "G4Log.hh"
51#include "G4Exp.hh"
52#include "G4Threading.hh"
53
54class G4Material;
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
59
61{
62public:
63
66
67 // parameters for mean energy loss calculation:
68 inline
70
73
74 inline
76 inline
78 inline
79 G4double GetTaul() const {return fTaul;};
80
81 // parameters of the density correction:
82 inline
84 inline
86 inline
87 G4double GetCdensity() const {return fCdensity;};
88 inline
89 G4double GetMdensity() const {return fMdensity;};
90 inline
91 G4double GetAdensity() const {return fAdensity;};
92 inline
93 G4double GetX0density() const {return fX0density;};
94 inline
95 G4double GetX1density() const {return fX1density;};
96 inline
97 G4double GetD0density() const {return fD0density;};
98
99 // user defined density correction parameterisation
101 G4double x0, G4double x1, G4double d0);
102
103 // defined density correction parameterisation via base material
104 void SetDensityEffectParameters(const G4Material* bmat);
105
107
108 // compute density correction as a function of the kinematic variable
109 // x = log10(beta*gamma)
111
113 { return fDensityEffectCalc; }
114
115 // use parameterisation
117
119
120 // parameters of the energy loss fluctuation model:
121 inline
122 G4double GetF1fluct() const {return fF1fluct;};
123 inline
124 G4double GetF2fluct() const {return fF2fluct;};
125 inline
127 inline
129 inline
131 inline
133 inline
135 inline
137
138 // parameters for ion corrections computations
139 inline
140 G4double GetZeffective() const {return fZeff;};
141 inline
143 inline
144 G4double GetLFactor() const {return fLfactor;};
145 inline
146 G4double GetInvA23() const {return fInvA23;};
147
148 // parameters for Birks attenuation:
149 inline
150 void SetBirksConstant(G4double value) {fBirks = value;};
151 inline
153
154 // parameters for average energy per ion
155 inline
157 inline
159
160 G4IonisParamMat(__void__&);
161 // Fake default constructor for usage restricted to direct object
162 // persistency for clients requiring preallocation of memory for
163 // persistifiable objects.
164
165 // operators
167 G4bool operator==(const G4IonisParamMat&) const = delete;
168 G4bool operator!=(const G4IonisParamMat&) const = delete;
170
171private:
172
173 // Compute mean parameters : ExcitationEnergy,Shell corretion vector ...
175
176 // Compute parameters for the density effect
178
179 // Compute parameters for the energy fluctuation model
180 void ComputeFluctModel();
181
182 // Compute parameters for ion parameterizations
184
185 //
186 // data members
187 //
188 const G4Material* fMaterial; // this material
189
190 G4DensityEffectCalculator* fDensityEffectCalc; // calculator of the density effect
191 G4double* fShellCorrectionVector; // shell correction coefficients
192
193 // parameters for mean energy loss calculation
196 G4double fTaul; // lower limit of Bethe-Bloch formula
197
198 // parameters of the density correction
199 G4double fCdensity; // mat.constant
200 G4double fMdensity; // exponent
205
208
209 // parameters of the energy loss fluctuation model
218
219 // parameters for ion corrections computations
224
225 // parameter for Birks attenuation
227 // average energy per ion pair
229
230 // static data created only once
233#ifdef G4MULTITHREADED
234 static G4Mutex ionisMutex;
235#endif
236};
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
239
241{
242 // x = log10(beta*gamma)
243 G4double y = 0.0;
244 if(x < fX0density) {
245 if(fD0density > 0.0) { y = fD0density*G4Exp(twoln10*(x - fX0density)); }
246 } else if(x >= fX1density) { y = twoln10*x - fCdensity; }
248 return y;
249}
250
251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
252
253#endif
static const G4double cd
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4double GetF2fluct() const
G4double GetLogEnergy1fluct() const
G4bool operator==(const G4IonisParamMat &) const =delete
G4double * fShellCorrectionVector
G4IonisParamMat & operator=(const G4IonisParamMat &)=delete
G4double GetMdensity() const
G4DensityEffectCalculator * fDensityEffectCalc
G4double GetX1density() const
G4double * GetShellCorrectionVector() const
G4double GetTaul() const
G4double GetAdjustmentFactor() const
G4double GetX0density() const
G4double GetCdensity() const
G4double GetLogEnergy2fluct() const
G4double GetDensityCorrection(G4double x)
void SetDensityEffectParameters(G4double cd, G4double md, G4double ad, G4double x0, G4double x1, G4double d0)
static G4DensityEffectData * fDensityData
G4double GetMeanExcitationEnergy() const
G4double GetF1fluct() const
G4double GetMeanEnergyPerIonPair() const
static G4DensityEffectData * GetDensityEffectData()
G4double fLogMeanExcEnergy
G4double DensityCorrection(G4double x)
G4double FindMeanExcitationEnergy(const G4Material *) const
void ComputeDensityEffectOnFly(G4bool)
void ComputeMeanParameters()
G4double GetInvA23() const
G4double fMeanExcitationEnergy
G4double GetD0density() const
G4double GetEnergy2fluct() const
G4double GetZeffective() const
G4IonisParamMat(const G4IonisParamMat &)=delete
G4IonisParamMat(const G4Material *)
G4double GetAdensity() const
void ComputeDensityEffectParameters()
G4double GetEnergy1fluct() const
G4double GetFermiEnergy() const
G4double fLogEnergy2fluct
G4bool operator!=(const G4IonisParamMat &) const =delete
void SetMeanExcitationEnergy(G4double value)
G4double fAdjustmentFactor
G4double GetPlasmaEnergy() const
G4double GetLFactor() const
G4double fRateionexcfluct
void SetBirksConstant(G4double value)
G4double fMeanEnergyPerIon
void SetMeanEnergyPerIonPair(G4double value)
G4DensityEffectCalculator * GetDensityEffectCalculator()
G4double GetRateionexcfluct() const
G4double GetEnergy0fluct() const
G4double GetLogMeanExcEnergy() const
G4double fLogEnergy1fluct
G4double GetBirksConstant() const
const G4Material * fMaterial