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

Last change on this file was 1340, checked in by garnier, 15 years ago

update ti head

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.8 2010/10/26 15:40:03 vnivanch Exp $
27// GEANT4 tag $Name: emhighenergy-V09-03-02 $
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 mass = 0.0;
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
87G4mplIonisationModel::~G4mplIonisationModel()
88{}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
92void G4mplIonisationModel::Initialise(const G4ParticleDefinition* p,
93 const G4DataVector&)
94{
95 monopole = p;
96 mass = monopole->GetPDGMass();
97 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101
102G4double G4mplIonisationModel::ComputeDEDXPerVolume(const G4Material* material,
103 const G4ParticleDefinition*,
104 G4double kineticEnergy,
105 G4double)
106{
107 G4double tau = kineticEnergy / mass;
108 G4double gam = tau + 1.0;
109 G4double bg2 = tau * (tau + 2.0);
110 G4double beta2 = bg2 / (gam * gam);
111 G4double beta = sqrt(beta2);
112
113 // low-energy asymptotic formula
114 G4double dedx = dedxlim*beta*material->GetDensity();
115
116 // above asymptotic
117 if(beta > betalow) {
118
119 // high energy
120 if(beta >= betalim) {
121 dedx = ComputeDEDXAhlen(material, bg2);
122
123 } else {
124
125 G4double dedx1 = dedxlim*betalow*material->GetDensity();
126 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim);
127
128 // extrapolation between two formula
129 G4double kapa2 = beta - betalow;
130 G4double kapa1 = betalim - beta;
131 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
132 }
133 }
134 return dedx;
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
139G4double G4mplIonisationModel::ComputeDEDXAhlen(const G4Material* material,
140 G4double bg2)
141{
142 G4double eDensity = material->GetElectronDensity();
143 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
144 G4double cden = material->GetIonisation()->GetCdensity();
145 G4double mden = material->GetIonisation()->GetMdensity();
146 G4double aden = material->GetIonisation()->GetAdensity();
147 G4double x0den = material->GetIonisation()->GetX0density();
148 G4double x1den = material->GetIonisation()->GetX1density();
149
150 // Ahlen's formula for nonconductors, [1]p157, f(5.7)
151 G4double dedx = log(2.0 * electron_mass_c2 * bg2 / eexc) - 0.5;
152
153 // Kazama et al. cross-section correction
154 G4double k = 0.406;
155 if(nmpl > 1) k = 0.346;
156
157 // Bloch correction
158 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
159
160 dedx += 0.5 * k - B[nmpl];
161
162 // density effect correction
163 G4double deltam;
164 G4double x = log(bg2) / twoln10;
165 if ( x >= x0den ) {
166 deltam = twoln10 * x - cden;
167 if ( x < x1den ) deltam += aden * pow((x1den-x), mden);
168 dedx -= 0.5 * deltam;
169 }
170
171 // now compute the total ionization loss
172 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
173
174 if (dedx < 0.0) dedx = 0;
175 return dedx;
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
179
180void G4mplIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
181 const G4MaterialCutsCouple*,
182 const G4DynamicParticle*,
183 G4double,
184 G4double)
185{}
186
187//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188
189G4double G4mplIonisationModel::SampleFluctuations(
190 const G4Material* material,
191 const G4DynamicParticle* dp,
192 G4double& tmax,
193 G4double& length,
194 G4double& meanLoss)
195{
196 G4double siga = Dispersion(material,dp,tmax,length);
197 G4double loss = meanLoss;
198 siga = sqrt(siga);
199 G4double twomeanLoss = meanLoss + meanLoss;
200
201 if(twomeanLoss < siga) {
202 G4double x;
203 do {
204 loss = twomeanLoss*G4UniformRand();
205 x = (loss - meanLoss)/siga;
206 } while (1.0 - 0.5*x*x < G4UniformRand());
207 } else {
208 do {
209 loss = G4RandGauss::shoot(meanLoss,siga);
210 } while (0.0 > loss || loss > twomeanLoss);
211 }
212 return loss;
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
216
217G4double G4mplIonisationModel::Dispersion(const G4Material* material,
218 const G4DynamicParticle* dp,
219 G4double& tmax,
220 G4double& length)
221{
222 G4double siga = 0.0;
223 G4double tau = dp->GetKineticEnergy()/mass;
224 if(tau > 0.0) {
225 G4double electronDensity = material->GetElectronDensity();
226 G4double gam = tau + 1.0;
227 G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
228 siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
229 * electronDensity * chargeSquare;
230 }
231 return siga;
232}
233
234//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.