Geant4  10.01
G4hBetheBlochModel.cc
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 //
29 // GEANT4 Class file
30 //
31 //
32 // File name: G4hBetheBlochModel
33 //
34 // Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
35 //
36 // Creation date: 20 July 2000
37 //
38 // Modifications:
39 // 20/07/2000 V.Ivanchenko First implementation
40 // 03/10/2000 V.Ivanchenko clean up accoding to CodeWizard
41 //
42 // Class Description:
43 //
44 // Bethe-Bloch ionisation model
45 //
46 // Class Description: End
47 //
48 // -------------------------------------------------------------------
49 //
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51 
52 #include "G4hBetheBlochModel.hh"
53 
54 #include "globals.hh"
55 #include "G4PhysicalConstants.hh"
56 #include "G4SystemOfUnits.hh"
57 #include "G4DynamicParticle.hh"
58 #include "G4ParticleDefinition.hh"
59 #include "G4Material.hh"
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62 
64  : G4VLowEnergyModel(name),
65  lowEnergyLimit(1.*MeV),
66  highEnergyLimit(100.*GeV),
67  twoln10(2.*std::log(10.)),
68  bg2lim(0.0169),
69  taulim(8.4146e-3)
70 {;}
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73 
75 {;}
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
80  const G4Material* material)
81 {
82  G4double energy = particle->GetKineticEnergy() ;
83  G4double particleMass = particle->GetMass() ;
84 
85  G4double eloss = BetheBlochFormula(material,energy,particleMass) ;
86 
87  return eloss ;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93  const G4Material* material,
94  G4double kineticEnergy)
95 {
96  G4double particleMass = aParticle->GetPDGMass() ;
97  G4double eloss = BetheBlochFormula(material,kineticEnergy,particleMass) ;
98 
99  return eloss ;
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103 
105  const G4ParticleDefinition* ,
106  const G4Material* ) const
107 {
108  return highEnergyLimit ;
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112 
114  const G4ParticleDefinition* aParticle,
115  const G4Material* material) const
116 {
117  G4double taul = (material->GetIonisation()->GetTaul())*
118  (aParticle->GetPDGMass()) ;
119  return taul ;
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 
125  const G4ParticleDefinition* ) const
126 {
127  return highEnergyLimit ;
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 
133  const G4ParticleDefinition* ) const
134 {
135  return lowEnergyLimit ;
136 }
137 
138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139 
141  const G4Material* ) const
142 {
143  return true ;
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147 
149  const G4Material* ) const
150 {
151  return true ;
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
157  const G4Material* material,
158  G4double kineticEnergy,
159  G4double particleMass) const
160 {
161  // This member function is applied normally to proton/antiproton
162  G4double ionloss ;
163 
164  G4double rateMass = electron_mass_c2/particleMass ;
165 
166  G4double taul = material->GetIonisation()->GetTaul() ;
167  G4double tau = kineticEnergy/particleMass ; // tau is relative energy
168 
169  // It is not normal case for this function
170  // for low energy parametrisation have to be applied
171  if ( tau < taul ) tau = taul ;
172 
173  // some local variables
174 
175  G4double gamma,bg2,beta2,tmax,x,delta,sh ;
176  G4double electronDensity = material->GetElectronDensity();
177  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
178  G4double eexc2 = eexc*eexc ;
179  G4double cden = material->GetIonisation()->GetCdensity();
180  G4double mden = material->GetIonisation()->GetMdensity();
181  G4double aden = material->GetIonisation()->GetAdensity();
182  G4double x0den = material->GetIonisation()->GetX0density();
183  G4double x1den = material->GetIonisation()->GetX1density();
184  G4double* shellCorrectionVector =
186 
187  gamma = tau + 1.0 ;
188  bg2 = tau*(tau+2.0) ;
189  beta2 = bg2/(gamma*gamma) ;
190  tmax = 2.*electron_mass_c2*bg2/(1.+2.*gamma*rateMass+rateMass*rateMass) ;
191 
192  ionloss = std::log(2.0*electron_mass_c2*bg2*tmax/eexc2)-2.0*beta2 ;
193 
194  // density correction
195  x = std::log(bg2)/twoln10 ;
196  if ( x < x0den ) {
197  delta = 0.0 ;
198 
199  } else {
200  delta = twoln10*x - cden ;
201  if ( x < x1den ) delta += aden*std::pow((x1den-x),mden) ;
202  }
203 
204  // shell correction
205  sh = 0.0 ;
206  x = 1.0 ;
207 
208  if ( bg2 > bg2lim ) {
209  for (G4int k=0; k<=2; k++) {
210  x *= bg2 ;
211  sh += shellCorrectionVector[k]/x;
212  }
213 
214  } else {
215  for (G4int k=0; k<=2; k++) {
216  x *= bg2lim ;
217  sh += shellCorrectionVector[k]/x;
218  }
219  sh *= std::log(tau/taul)/std::log(taulim/taul) ;
220  }
221 
222  // now compute the total ionization loss
223 
224  ionloss -= delta + sh ;
225  ionloss *= twopi_mc2_rcl2*electronDensity/beta2 ;
226 
227  if ( ionloss < 0.0) ionloss = 0.0 ;
228 
229  return ionloss;
230 }
231 
232 
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
static const double MeV
Definition: G4SIunits.hh:193
G4double GetAdensity() const
G4double GetKineticEnergy() const
G4String name
Definition: TRTMaterials.hh:40
G4double GetX1density() const
G4double BetheBlochFormula(const G4Material *material, G4double kineticEnergy, G4double particleMass) const
G4double GetTaul() const
G4bool IsInCharge(const G4DynamicParticle *particle, const G4Material *material) const
int G4int
Definition: G4Types.hh:78
static const G4double bg2lim
G4double GetElectronDensity() const
Definition: G4Material.hh:215
G4double GetMass() const
bool G4bool
Definition: G4Types.hh:79
G4double HighEnergyLimit(const G4ParticleDefinition *aParticle, const G4Material *material) const
static const double GeV
Definition: G4SIunits.hh:196
G4double TheValue(const G4DynamicParticle *particle, const G4Material *material)
G4hBetheBlochModel(const G4String &name)
G4double GetX0density() const
G4double GetPDGMass() const
G4double GetCdensity() const
G4double energy(const ThreeVector &p, const G4double m)
G4double * GetShellCorrectionVector() const
G4double GetMeanExcitationEnergy() const
static const G4double twoln10
double G4double
Definition: G4Types.hh:76
G4double GetMdensity() const
G4double LowEnergyLimit(const G4ParticleDefinition *aParticle, const G4Material *material) const