source: trunk/source/processes/electromagnetic/highenergy/src/G4mplIonisationModel.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

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