Geant4 9.6.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// $Id$
27//
28
29// class description
30//
31// The class contains few (physical) quantities related to the Ionisation
32// process, for a material defined by its pointer G4Material*
33//
34
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
36
37// 09-07-98: data moved from G4Material (mma)
38// 09-03-01: copy constructor and assignement operator in public (mma)
39// 28-10-02: add setMeanExcitationEnergy (V.Ivanchenko)
40// 27-09-07: add computation of parameters for ions (V.Ivanchenko)
41// 04-03-08: add fBirks constant (mma)
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
51class G4Material; // forward declaration
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
55
56class G4IonisParamMat // with description
57{
58public:
59
61 virtual ~G4IonisParamMat();
62
63 //
64 // retrieval methods
65 //
66
67 // parameters for mean energy loss calculation:
68 G4double GetMeanExcitationEnergy() const {return fMeanExcitationEnergy;};
71 G4double GetLogMeanExcEnergy() const {return fLogMeanExcEnergy;};
72 G4double* GetShellCorrectionVector() const {return fShellCorrectionVector;};
73 G4double GetTaul() const {return fTaul;};
74
75 // parameters of the density correction:
76 G4double GetPlasmaEnergy() const {return fPlasmaEnergy;};
77 G4double GetAdjustmentFactor() const {return fAdjustmentFactor;};
78 G4double GetCdensity() const {return fCdensity;};
79 G4double GetMdensity() const {return fMdensity;};
80 G4double GetAdensity() const {return fAdensity;};
81 G4double GetX0density() const {return fX0density;};
82 G4double GetX1density() const {return fX1density;};
83 G4double GetD0density() const {return fD0density;};
84
85 // compute density correction as a function of the kinematic variable
86 // x = log10(beta*gamma)
89
90 // parameters of the energy loss fluctuation model:
91 G4double GetF1fluct() const {return fF1fluct;};
92 G4double GetF2fluct() const {return fF2fluct;};
93 G4double GetEnergy1fluct() const {return fEnergy1fluct;};
94 G4double GetLogEnergy1fluct() const {return fLogEnergy1fluct;};
95 G4double GetEnergy2fluct() const {return fEnergy2fluct;};
96 G4double GetLogEnergy2fluct() const {return fLogEnergy2fluct;};
97 G4double GetEnergy0fluct() const {return fEnergy0fluct;};
98 G4double GetRateionexcfluct() const {return fRateionexcfluct;};
99
100 // parameters for ion corrections computations
101 G4double GetZeffective() const {return fZeff;};
102 G4double GetFermiEnergy() const {return fFermiEnergy;};
103 G4double GetLFactor() const {return fLfactor;};
104 G4double GetInvA23() const {return fInvA23;};
105
106 // parameters for Birks attenuation:
107 void SetBirksConstant(G4double value) {fBirks = value;};
108 G4double GetBirksConstant() const {return fBirks;};
109
110 // parameters for average energy per ion
111 void SetMeanEnergyPerIonPair(G4double value) {fMeanEnergyPerIon = value;};
112 G4double GetMeanEnergyPerIonPair() const {return fMeanEnergyPerIon;};
113
114public: // without description
115
118 G4int operator==(const G4IonisParamMat&) const;
119 G4int operator!=(const G4IonisParamMat&) const;
120
121 G4IonisParamMat(__void__&);
122 // Fake default constructor for usage restricted to direct object
123 // persistency for clients requiring preallocation of memory for
124 // persistifiable objects.
125
126private:
127
128 // Compute mean parameters : ExcitationEnergy,Shell corretion vector ...
129 void ComputeMeanParameters();
130
131 // Compute parameters for the density effect
132 void ComputeDensityEffect();
133
134 // Compute parameters for the energy fluctuation model
135 void ComputeFluctModel();
136
137 // Compute parameters for ion parameterizations
138 void ComputeIonParameters();
139
140private:
141
142//
143// data members
144//
145 G4Material* fMaterial; // this material
146
147 // parameters for mean energy loss calculation
148 G4double fMeanExcitationEnergy; //
149 G4double fLogMeanExcEnergy; //
150 G4double* fShellCorrectionVector; // shell correction coefficients
151 G4double fTaul; // lower limit of Bethe-Bloch formula
152
153 // parameters of the density correction
154 G4double fCdensity; // mat.constant
155 G4double fMdensity; // exponent
156 G4double fAdensity; //
157 G4double fX0density; //
158 G4double fX1density; //
159 G4double fD0density;
160
161 G4double fPlasmaEnergy;
162 G4double fAdjustmentFactor;
163
164 // parameters of the energy loss fluctuation model
165 G4double fF1fluct;
166 G4double fF2fluct;
167 G4double fEnergy1fluct;
168 G4double fLogEnergy1fluct;
169 G4double fEnergy2fluct;
170 G4double fLogEnergy2fluct;
171 G4double fEnergy0fluct;
172 G4double fRateionexcfluct;
173
174 // parameters for ion corrections computations
175 G4double fZeff;
176 G4double fFermiEnergy;
177 G4double fLfactor;
178 G4double fInvA23;
179
180 // parameter for Birks attenuation
181 G4double fBirks;
182 // average energy per ion pair
183 G4double fMeanEnergyPerIon;
184
185 // static data created only once
186 static G4DensityEffectData* fDensityData;
187 G4double twoln10;
188};
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
191
193{
194 // x = log10(beta*gamma)
195
196 G4double y = 0.0;
197 if(x < fX0density) {
198 if(fD0density > 0.0) { y = fD0density*std::pow(10.,2*(x - fX0density)); }
199 } else if(x >= fX1density) { y = twoln10*x - fCdensity; }
200 else {y = twoln10*x - fCdensity + fAdensity*std::pow(fX1density - x, fMdensity);}
201 return y;
202}
203
204#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetF2fluct() const
G4double GetLogEnergy1fluct() const
G4double GetMdensity() const
G4IonisParamMat & operator=(const G4IonisParamMat &)
G4double GetX1density() const
G4double * GetShellCorrectionVector() const
G4double GetTaul() const
G4double GetAdjustmentFactor() const
G4double FindMeanExcitationEnergy(const G4String &chFormula)
G4double GetX0density() const
G4double GetCdensity() const
G4double GetLogEnergy2fluct() const
G4double GetMeanExcitationEnergy() const
G4double GetF1fluct() const
G4double GetMeanEnergyPerIonPair() const
static G4DensityEffectData * GetDensityEffectData()
G4double DensityCorrection(G4double x)
G4int operator==(const G4IonisParamMat &) const
G4double GetInvA23() const
G4double GetD0density() const
virtual ~G4IonisParamMat()
G4double GetEnergy2fluct() const
G4double GetZeffective() const
G4double GetAdensity() const
G4double GetEnergy1fluct() const
G4double GetFermiEnergy() const
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
G4int operator!=(const G4IonisParamMat &) const