Geant4  10.02
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: G4IonisParamMat.hh 93568 2015-10-26 14:52:36Z gcosmo $
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 #include "G4Log.hh"
51 #include "G4Exp.hh"
52 
53 class G4Material; // forward declaration
55 
56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
57 
58 class G4IonisParamMat // with description
59 {
60 public:
61 
63  virtual ~G4IonisParamMat();
64 
65  //
66  // retrieval methods
67  //
68 
69  // parameters for mean energy loss calculation:
70  inline
72 
74  G4double FindMeanExcitationEnergy(const G4String& chFormula);
75 
76  inline
78  inline
80  inline
81  G4double GetTaul() const {return fTaul;};
82 
83  // parameters of the density correction:
84  inline
86  inline
88  inline
89  G4double GetCdensity() const {return fCdensity;};
90  inline
91  G4double GetMdensity() const {return fMdensity;};
92  inline
93  G4double GetAdensity() const {return fAdensity;};
94  inline
95  G4double GetX0density() const {return fX0density;};
96  inline
97  G4double GetX1density() const {return fX1density;};
98  inline
99  G4double GetD0density() const {return fD0density;};
100 
101  // compute density correction as a function of the kinematic variable
102  // x = log10(beta*gamma)
104 
106 
107  // parameters of the energy loss fluctuation model:
108  inline
109  G4double GetF1fluct() const {return fF1fluct;};
110  inline
111  G4double GetF2fluct() const {return fF2fluct;};
112  inline
114  inline
116  inline
118  inline
120  inline
122  inline
124 
125  // parameters for ion corrections computations
126  inline
127  G4double GetZeffective() const {return fZeff;};
128  inline
130  inline
131  G4double GetLFactor() const {return fLfactor;};
132  inline
133  G4double GetInvA23() const {return fInvA23;};
134 
135  // parameters for Birks attenuation:
136  inline
137  void SetBirksConstant(G4double value) {fBirks = value;};
138  inline
139  G4double GetBirksConstant() const {return fBirks;};
140 
141  // parameters for average energy per ion
142  inline
144  inline
146 
147  G4IonisParamMat(__void__&);
148  // Fake default constructor for usage restricted to direct object
149  // persistency for clients requiring preallocation of memory for
150  // persistifiable objects.
151 
152 private:
153 
154  // Compute mean parameters : ExcitationEnergy,Shell corretion vector ...
155  void ComputeMeanParameters();
156 
157  // Compute parameters for the density effect
158  void ComputeDensityEffect();
159 
160  // Compute parameters for the energy fluctuation model
161  void ComputeFluctModel();
162 
163  // Compute parameters for ion parameterizations
164  void ComputeIonParameters();
165 
166  // operators
168  G4int operator==(const G4IonisParamMat&) const;
169  G4int operator!=(const G4IonisParamMat&) const;
171 
172  //
173  // data members
174  //
175  G4Material* fMaterial; // this material
176 
177  // parameters for mean energy loss calculation
180  G4double* fShellCorrectionVector; // shell correction coefficients
181  G4double fTaul; // lower limit of Bethe-Bloch formula
182 
183  // parameters of the density correction
184  G4double fCdensity; // mat.constant
185  G4double fMdensity; // exponent
190 
193 
194  // parameters of the energy loss fluctuation model
203 
204  // parameters for ion corrections computations
209 
210  // parameter for Birks attenuation
212  // average energy per ion pair
214 
215  // static data created only once
218 };
219 
220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
221 
223 {
224  // x = log10(beta*gamma)
225 
226  G4double y = 0.0;
227  if(x < fX0density) {
228  if(fD0density > 0.0) { y = fD0density*G4Exp(twoln10*(x - fX0density)); }
229  } else if(x >= fX1density) { y = twoln10*x - fCdensity; }
230  else {y = twoln10*x - fCdensity + fAdensity*G4Exp(G4Log(fX1density - x)*fMdensity);}
231  return y;
232 }
233 
234 #endif
G4double fLogEnergy2fluct
G4double GetAdensity() const
G4double GetEnergy2fluct() const
G4double GetMeanEnergyPerIonPair() const
G4double FindMeanExcitationEnergy(const G4String &chFormula)
void SetMeanExcitationEnergy(G4double value)
void SetBirksConstant(G4double value)
G4double GetX1density() const
G4double fAdjustmentFactor
G4double GetLogEnergy2fluct() const
G4double GetTaul() const
G4double GetLogMeanExcEnergy() const
int G4int
Definition: G4Types.hh:78
G4double GetZeffective() const
G4IonisParamMat & operator=(const G4IonisParamMat &)
G4double GetBirksConstant() const
G4double GetFermiEnergy() const
G4double GetEnergy0fluct() const
G4double fRateionexcfluct
G4double GetInvA23() const
virtual ~G4IonisParamMat()
G4double GetAdjustmentFactor() const
G4double GetPlasmaEnergy() const
G4double fLogMeanExcEnergy
G4double fMeanEnergyPerIon
G4double fLogEnergy1fluct
static G4DensityEffectData * GetDensityEffectData()
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double fMeanExcitationEnergy
G4double * fShellCorrectionVector
G4double GetLogEnergy1fluct() const
G4double GetX0density() const
G4double GetRateionexcfluct() const
G4double DensityCorrection(G4double x)
G4double GetCdensity() const
G4Material * fMaterial
G4IonisParamMat(G4Material *)
G4double * GetShellCorrectionVector() const
const G4double x[NPOINTSGL]
G4int operator!=(const G4IonisParamMat &) const
G4double GetMeanExcitationEnergy() const
G4double GetF2fluct() const
G4double GetLFactor() const
double G4double
Definition: G4Types.hh:76
G4double GetD0density() const
G4double GetMdensity() const
G4int operator==(const G4IonisParamMat &) const
G4double GetF1fluct() const
static G4DensityEffectData * fDensityData
void ComputeMeanParameters()
G4double GetEnergy1fluct() const
void SetMeanEnergyPerIonPair(G4double value)