Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4mplIonisationModel.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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4mplIonisationModel
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 06.09.2005
38//
39// Modifications:
40// 12.08.2007 Changing low energy approximation and extrapolation.
41// Small bug fixing and refactoring (M. Vladymyrov)
42// 13.11.2007 Use low-energy asymptotic from [3] (V.Ivanchenko)
43//
44//
45// -------------------------------------------------------------------
46// References
47// [1] Steven P. Ahlen: Energy loss of relativistic heavy ionizing particles,
48// S.P. Ahlen, Rev. Mod. Phys 52(1980), p121
49// [2] K.A. Milton arXiv:hep-ex/0602040
50// [3] S.P. Ahlen and K. Kinoshita, Phys. Rev. D26 (1982) 2347
51
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55
57#include "Randomize.hh"
59#include "G4SystemOfUnits.hh"
60#include "G4LossTableManager.hh"
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
65using namespace std;
66
69 magCharge(mCharge),
70 twoln10(log(100.0)),
71 betalow(0.01),
72 betalim(0.1),
73 beta2lim(betalim*betalim),
74 bg2lim(beta2lim*(1.0 + beta2lim))
75{
76 nmpl = G4int(abs(magCharge) * 2 * fine_structure_const + 0.5);
77 if(nmpl > 6) { nmpl = 6; }
78 else if(nmpl < 1) { nmpl = 1; }
79 pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
80 chargeSquare = magCharge * magCharge;
81 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
82 fParticleChange = 0;
83 monopole = 0;
84 mass = 0.0;
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88
90{}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
95{
96 monopole = p;
97 mass = monopole->GetPDGMass();
98 G4double emin =
99 std::min(LowEnergyLimit(),0.1*mass*(1/sqrt(1 - betalow*betalow) - 1));
100 G4double emax =
101 std::max(HighEnergyLimit(),10*mass*(1/sqrt(1 - beta2lim) - 1));
102 SetLowEnergyLimit(emin);
103 SetHighEnergyLimit(emax);
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
109 const G4DataVector&)
110{
111 if(!monopole) { SetParticle(p); }
112 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
113}
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
118 const G4ParticleDefinition* p,
119 G4double kineticEnergy,
120 G4double)
121{
122 if(!monopole) { SetParticle(p); }
123 G4double tau = kineticEnergy / mass;
124 G4double gam = tau + 1.0;
125 G4double bg2 = tau * (tau + 2.0);
126 G4double beta2 = bg2 / (gam * gam);
127 G4double beta = sqrt(beta2);
128
129 // low-energy asymptotic formula
130 G4double dedx = dedxlim*beta*material->GetDensity();
131
132 // above asymptotic
133 if(beta > betalow) {
134
135 // high energy
136 if(beta >= betalim) {
137 dedx = ComputeDEDXAhlen(material, bg2);
138
139 } else {
140
141 G4double dedx1 = dedxlim*betalow*material->GetDensity();
142 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim);
143
144 // extrapolation between two formula
145 G4double kapa2 = beta - betalow;
146 G4double kapa1 = betalim - beta;
147 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
148 }
149 }
150 return dedx;
151}
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154
155G4double G4mplIonisationModel::ComputeDEDXAhlen(const G4Material* material,
156 G4double bg2)
157{
158 G4double eDensity = material->GetElectronDensity();
159 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
160 G4double cden = material->GetIonisation()->GetCdensity();
161 G4double mden = material->GetIonisation()->GetMdensity();
162 G4double aden = material->GetIonisation()->GetAdensity();
163 G4double x0den = material->GetIonisation()->GetX0density();
164 G4double x1den = material->GetIonisation()->GetX1density();
165
166 // Ahlen's formula for nonconductors, [1]p157, f(5.7)
167 G4double dedx = log(2.0 * electron_mass_c2 * bg2 / eexc) - 0.5;
168
169 // Kazama et al. cross-section correction
170 G4double k = 0.406;
171 if(nmpl > 1) k = 0.346;
172
173 // Bloch correction
174 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
175
176 dedx += 0.5 * k - B[nmpl];
177
178 // density effect correction
179 G4double deltam;
180 G4double x = log(bg2) / twoln10;
181 if ( x >= x0den ) {
182 deltam = twoln10 * x - cden;
183 if ( x < x1den ) deltam += aden * pow((x1den-x), mden);
184 dedx -= 0.5 * deltam;
185 }
186
187 // now compute the total ionization loss
188 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
189
190 if (dedx < 0.0) dedx = 0;
191 return dedx;
192}
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
195
196void G4mplIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
198 const G4DynamicParticle*,
199 G4double,
200 G4double)
201{}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
204
206 const G4Material* material,
207 const G4DynamicParticle* dp,
208 G4double& tmax,
209 G4double& length,
210 G4double& meanLoss)
211{
212 G4double siga = Dispersion(material,dp,tmax,length);
213 G4double loss = meanLoss;
214 siga = sqrt(siga);
215 G4double twomeanLoss = meanLoss + meanLoss;
216
217 if(twomeanLoss < siga) {
218 G4double x;
219 do {
220 loss = twomeanLoss*G4UniformRand();
221 x = (loss - meanLoss)/siga;
222 } while (1.0 - 0.5*x*x < G4UniformRand());
223 } else {
224 do {
225 loss = G4RandGauss::shoot(meanLoss,siga);
226 } while (0.0 > loss || loss > twomeanLoss);
227 }
228 return loss;
229}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
232
234 const G4DynamicParticle* dp,
235 G4double& tmax,
236 G4double& length)
237{
238 G4double siga = 0.0;
239 G4double tau = dp->GetKineticEnergy()/mass;
240 if(tau > 0.0) {
241 G4double electronDensity = material->GetElectronDensity();
242 G4double gam = tau + 1.0;
243 G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
244 siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
245 * electronDensity * chargeSquare;
246 }
247 return siga;
248}
249
250//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetKineticEnergy() const
G4double GetMdensity() const
G4double GetX1density() const
G4double GetX0density() const
G4double GetCdensity() const
G4double GetMeanExcitationEnergy() const
G4double GetAdensity() const
G4double GetDensity() const
Definition: G4Material.hh:179
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
G4double GetElectronDensity() const
Definition: G4Material.hh:216
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95
G4mplIonisationModel(G4double mCharge, const G4String &nam="mplIonisation")
virtual G4double SampleFluctuations(const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length, G4double &meanLoss)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
void SetParticle(const G4ParticleDefinition *p)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)