Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // $Id: G4IonisParamMat.hh 97248 2016-05-30 15:00:11Z 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 
62  G4IonisParamMat(const G4Material*);
64 
65  // parameters for mean energy loss calculation:
66  inline
67  G4double GetMeanExcitationEnergy() const {return fMeanExcitationEnergy;};
68 
71 
72  inline
73  G4double GetLogMeanExcEnergy() const {return fLogMeanExcEnergy;};
74  inline
75  G4double* GetShellCorrectionVector() const {return fShellCorrectionVector;};
76  inline
77  G4double GetTaul() const {return fTaul;};
78 
79  // parameters of the density correction:
80  inline
81  G4double GetPlasmaEnergy() const {return fPlasmaEnergy;};
82  inline
83  G4double GetAdjustmentFactor() const {return fAdjustmentFactor;};
84  inline
85  G4double GetCdensity() const {return fCdensity;};
86  inline
87  G4double GetMdensity() const {return fMdensity;};
88  inline
89  G4double GetAdensity() const {return fAdensity;};
90  inline
91  G4double GetX0density() const {return fX0density;};
92  inline
93  G4double GetX1density() const {return fX1density;};
94  inline
95  G4double GetD0density() const {return fD0density;};
96 
97  // compute density correction as a function of the kinematic variable
98  // x = log10(beta*gamma)
100 
102 
103  // parameters of the energy loss fluctuation model:
104  inline
105  G4double GetF1fluct() const {return fF1fluct;};
106  inline
107  G4double GetF2fluct() const {return fF2fluct;};
108  inline
109  G4double GetEnergy1fluct() const {return fEnergy1fluct;};
110  inline
111  G4double GetLogEnergy1fluct() const {return fLogEnergy1fluct;};
112  inline
113  G4double GetEnergy2fluct() const {return fEnergy2fluct;};
114  inline
115  G4double GetLogEnergy2fluct() const {return fLogEnergy2fluct;};
116  inline
117  G4double GetEnergy0fluct() const {return fEnergy0fluct;};
118  inline
119  G4double GetRateionexcfluct() const {return fRateionexcfluct;};
120 
121  // parameters for ion corrections computations
122  inline
123  G4double GetZeffective() const {return fZeff;};
124  inline
125  G4double GetFermiEnergy() const {return fFermiEnergy;};
126  inline
127  G4double GetLFactor() const {return fLfactor;};
128  inline
129  G4double GetInvA23() const {return fInvA23;};
130 
131  // parameters for Birks attenuation:
132  inline
133  void SetBirksConstant(G4double value) {fBirks = value;};
134  inline
135  G4double GetBirksConstant() const {return fBirks;};
136 
137  // parameters for average energy per ion
138  inline
139  void SetMeanEnergyPerIonPair(G4double value) {fMeanEnergyPerIon = value;};
140  inline
141  G4double GetMeanEnergyPerIonPair() const {return fMeanEnergyPerIon;};
142 
143  G4IonisParamMat(__void__&);
144  // Fake default constructor for usage restricted to direct object
145  // persistency for clients requiring preallocation of memory for
146  // persistifiable objects.
147 
148 private:
149 
150  // Compute mean parameters : ExcitationEnergy,Shell corretion vector ...
151  void ComputeMeanParameters();
152 
153  // Compute parameters for the density effect
154  void ComputeDensityEffect();
155 
156  // Compute parameters for the energy fluctuation model
157  void ComputeFluctModel();
158 
159  // Compute parameters for ion parameterizations
160  void ComputeIonParameters();
161 
162  // operators
163  G4IonisParamMat& operator=(const G4IonisParamMat&) = delete;
164  G4int operator==(const G4IonisParamMat&) const = delete;
165  G4int operator!=(const G4IonisParamMat&) const = delete;
166  G4IonisParamMat(const G4IonisParamMat&) = delete;
167 
168  //
169  // data members
170  //
171  const G4Material* fMaterial; // this material
172 
173  // parameters for mean energy loss calculation
174  G4double fMeanExcitationEnergy; //
175  G4double fLogMeanExcEnergy; //
176  G4double* fShellCorrectionVector; // shell correction coefficients
177  G4double fTaul; // lower limit of Bethe-Bloch formula
178 
179  // parameters of the density correction
180  G4double fCdensity; // mat.constant
181  G4double fMdensity; // exponent
182  G4double fAdensity; //
183  G4double fX0density; //
184  G4double fX1density; //
185  G4double fD0density;
186 
187  G4double fPlasmaEnergy;
188  G4double fAdjustmentFactor;
189 
190  // parameters of the energy loss fluctuation model
191  G4double fF1fluct;
192  G4double fF2fluct;
193  G4double fEnergy1fluct;
194  G4double fLogEnergy1fluct;
195  G4double fEnergy2fluct;
196  G4double fLogEnergy2fluct;
197  G4double fEnergy0fluct;
198  G4double fRateionexcfluct;
199 
200  // parameters for ion corrections computations
201  G4double fZeff;
202  G4double fFermiEnergy;
203  G4double fLfactor;
204  G4double fInvA23;
205 
206  // parameter for Birks attenuation
207  G4double fBirks;
208  // average energy per ion pair
209  G4double fMeanEnergyPerIon;
210 
211  // static data created only once
212  static G4DensityEffectData* fDensityData;
213  G4double twoln10;
214 };
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
217 
219 {
220  // x = log10(beta*gamma)
221  G4double y = 0.0;
222  if(x < fX0density) {
223  if(fD0density > 0.0) { y = fD0density*G4Exp(twoln10*(x - fX0density)); }
224  } else if(x >= fX1density) { y = twoln10*x - fCdensity; }
225  else {y = twoln10*x - fCdensity + fAdensity*G4Exp(G4Log(fX1density - x)*fMdensity);}
226  return y;
227 }
228 
229 #endif
G4double GetAdensity() const
G4double GetEnergy2fluct() const
G4double GetMeanEnergyPerIonPair() const
void SetMeanExcitationEnergy(G4double value)
void SetBirksConstant(G4double value)
G4double GetX1density() const
G4double GetLogEnergy2fluct() const
G4IonisParamMat(const G4Material *)
G4double GetTaul() const
tuple x
Definition: test.py:50
G4double GetLogMeanExcEnergy() const
int G4int
Definition: G4Types.hh:78
G4double GetZeffective() const
G4double GetBirksConstant() const
G4double GetFermiEnergy() const
G4double GetEnergy0fluct() const
G4double FindMeanExcitationEnergy(const G4Material *) const
const XML_Char int const XML_Char * value
Definition: expat.h:331
G4double GetInvA23() const
G4double GetAdjustmentFactor() const
G4double GetPlasmaEnergy() const
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 GetLogEnergy1fluct() const
G4double GetX0density() const
G4double GetRateionexcfluct() const
G4double DensityCorrection(G4double x)
G4double GetCdensity() const
G4double * GetShellCorrectionVector() const
G4double GetMeanExcitationEnergy() const
G4double GetF2fluct() const
G4double GetLFactor() const
double G4double
Definition: G4Types.hh:76
G4double GetD0density() const
G4double GetMdensity() const
G4double GetF1fluct() const
G4double GetEnergy1fluct() const
void SetMeanEnergyPerIonPair(G4double value)