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

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

update to last version 4.9.4

File size: 11.0 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: G4mplIonisationWithDeltaModel.cc,v 1.1 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: G4mplIonisationWithDeltaModel
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 "G4mplIonisationWithDeltaModel.hh"
58#include "Randomize.hh"
59#include "G4LossTableManager.hh"
60#include "G4ParticleChangeForLoss.hh"
61#include "G4Electron.hh"
62#include "G4DynamicParticle.hh"
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66using namespace std;
67
68G4mplIonisationWithDeltaModel::G4mplIonisationWithDeltaModel(G4double mCharge, const G4String& nam)
69 : G4VEmModel(nam),G4VEmFluctuationModel(nam),
70 magCharge(mCharge),
71 twoln10(log(100.0)),
72 betalow(0.01),
73 betalim(0.1),
74 beta2lim(betalim*betalim),
75 bg2lim(beta2lim*(1.0 + beta2lim))
76{
77 nmpl = G4int(abs(magCharge) * 2 * fine_structure_const + 0.5);
78 if(nmpl > 6) { nmpl = 6; }
79 else if(nmpl < 1) { nmpl = 1; }
80 pi_hbarc2_over_mc2 = pi * hbarc * hbarc / electron_mass_c2;
81 chargeSquare = magCharge * magCharge;
82 dedxlim = 45.*nmpl*nmpl*GeV*cm2/g;
83 fParticleChange = 0;
84 theElectron = G4Electron::Electron();
85 G4cout << "### Monopole ionisation model with d-electron production, Gmag= "
86 << magCharge/eplus << G4endl;
87 mass = 0.0;
88}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
92G4mplIonisationWithDeltaModel::~G4mplIonisationWithDeltaModel()
93{}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
97void
98G4mplIonisationWithDeltaModel::Initialise(const G4ParticleDefinition* p,
99 const G4DataVector&)
100{
101 monopole = p;
102 mass = monopole->GetPDGMass();
103 if(!fParticleChange) { fParticleChange = GetParticleChangeForLoss(); }
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
108G4double
109G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume(const G4Material* material,
110 const G4ParticleDefinition* p,
111 G4double kineticEnergy,
112 G4double maxEnergy)
113{
114 G4double tmax = MaxSecondaryEnergy(p,kineticEnergy);
115 G4double cutEnergy = std::min(tmax, maxEnergy);
116 G4double tau = kineticEnergy / mass;
117 G4double gam = tau + 1.0;
118 G4double bg2 = tau * (tau + 2.0);
119 G4double beta2 = bg2 / (gam * gam);
120 G4double beta = sqrt(beta2);
121
122 // low-energy asymptotic formula
123 G4double dedx = dedxlim*beta*material->GetDensity();
124
125 // above asymptotic
126 if(beta > betalow) {
127
128 // high energy
129 if(beta >= betalim) {
130 dedx = ComputeDEDXAhlen(material, bg2, cutEnergy);
131
132 } else {
133
134 G4double dedx1 = dedxlim*betalow*material->GetDensity();
135 G4double dedx2 = ComputeDEDXAhlen(material, bg2lim, cutEnergy);
136
137 // extrapolation between two formula
138 G4double kapa2 = beta - betalow;
139 G4double kapa1 = betalim - beta;
140 dedx = (kapa1*dedx1 + kapa2*dedx2)/(kapa1 + kapa2);
141 }
142 }
143 return dedx;
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147
148G4double
149G4mplIonisationWithDeltaModel::ComputeDEDXAhlen(const G4Material* material,
150 G4double bg2,
151 G4double cutEnergy)
152{
153 G4double eDensity = material->GetElectronDensity();
154 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
155
156 // Ahlen's formula for nonconductors, [1]p157, f(5.7)
157 G4double dedx =
158 0.5*(log(2.0 * electron_mass_c2 * bg2*cutEnergy / (eexc*eexc)) - 1.0);
159
160 // Kazama et al. cross-section correction
161 G4double k = 0.406;
162 if(nmpl > 1) { k = 0.346; }
163
164 // Bloch correction
165 const G4double B[7] = { 0.0, 0.248, 0.672, 1.022, 1.243, 1.464, 1.685};
166
167 dedx += 0.5 * k - B[nmpl];
168
169 // density effect correction
170 G4double x = log(bg2)/twoln10;
171 dedx -= material->GetIonisation()->DensityCorrection(x);
172
173 // now compute the total ionization loss
174 dedx *= pi_hbarc2_over_mc2 * eDensity * nmpl * nmpl;
175
176 if (dedx < 0.0) { dedx = 0; }
177 return dedx;
178}
179
180//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
181
182G4double
183G4mplIonisationWithDeltaModel::ComputeCrossSectionPerElectron(
184 const G4ParticleDefinition* p,
185 G4double kineticEnergy,
186 G4double cutEnergy,
187 G4double maxKinEnergy)
188{
189 G4double cross = 0.0;
190 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
191 G4double maxEnergy = min(tmax,maxKinEnergy);
192 if(cutEnergy < maxEnergy) {
193 cross = (1.0/cutEnergy - 1.0/maxEnergy)*twopi_mc2_rcl2*chargeSquare;
194 }
195 return cross;
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
199
200G4double
201G4mplIonisationWithDeltaModel::ComputeCrossSectionPerAtom(
202 const G4ParticleDefinition* p,
203 G4double kineticEnergy,
204 G4double Z, G4double,
205 G4double cutEnergy,
206 G4double maxEnergy)
207{
208 G4double cross =
209 Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
210 return cross;
211}
212
213//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
214
215void
216G4mplIonisationWithDeltaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
217 const G4MaterialCutsCouple*,
218 const G4DynamicParticle* dp,
219 G4double minKinEnergy,
220 G4double maxEnergy)
221{
222 G4double kineticEnergy = dp->GetKineticEnergy();
223 G4double tmax = MaxSecondaryEnergy(dp->GetDefinition(),kineticEnergy);
224
225 G4double maxKinEnergy = std::min(maxEnergy,tmax);
226 if(minKinEnergy >= maxKinEnergy) { return; }
227
228 //G4cout << "G4mplIonisationWithDeltaModel::SampleSecondaries: E(GeV)= "
229 // << kineticEnergy/GeV << " M(GeV)= " << mass/GeV
230 // << " tmin(MeV)= " << minKinEnergy/MeV << G4endl;
231
232 G4double totEnergy = kineticEnergy + mass;
233 G4double etot2 = totEnergy*totEnergy;
234 G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*mass)/etot2;
235
236 // sampling without nuclear size effect
237 G4double q = G4UniformRand();
238 G4double deltaKinEnergy = minKinEnergy*maxKinEnergy
239 /(minKinEnergy*(1.0 - q) + maxKinEnergy*q);
240
241 // delta-electron is produced
242 G4double totMomentum = totEnergy*sqrt(beta2);
243 G4double deltaMomentum =
244 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
245 G4double cost = deltaKinEnergy * (totEnergy + electron_mass_c2) /
246 (deltaMomentum * totMomentum);
247 if(cost > 1.0) { cost = 1.0; }
248
249 G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
250
251 G4double phi = twopi * G4UniformRand() ;
252
253 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost);
254 G4ThreeVector direction = dp->GetMomentumDirection();
255 deltaDirection.rotateUz(direction);
256
257 // create G4DynamicParticle object for delta ray
258 G4DynamicParticle* delta =
259 new G4DynamicParticle(theElectron,deltaDirection,deltaKinEnergy);
260
261 vdp->push_back(delta);
262
263 // Change kinematics of primary particle
264 kineticEnergy -= deltaKinEnergy;
265 G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
266 finalP = finalP.unit();
267
268 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
269 fParticleChange->SetProposedMomentumDirection(finalP);
270}
271
272//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
273
274G4double G4mplIonisationWithDeltaModel::SampleFluctuations(
275 const G4Material* material,
276 const G4DynamicParticle* dp,
277 G4double& tmax,
278 G4double& length,
279 G4double& meanLoss)
280{
281 G4double siga = Dispersion(material,dp,tmax,length);
282 G4double loss = meanLoss;
283 siga = sqrt(siga);
284 G4double twomeanLoss = meanLoss + meanLoss;
285
286 if(twomeanLoss < siga) {
287 G4double x;
288 do {
289 loss = twomeanLoss*G4UniformRand();
290 x = (loss - meanLoss)/siga;
291 } while (1.0 - 0.5*x*x < G4UniformRand());
292 } else {
293 do {
294 loss = G4RandGauss::shoot(meanLoss,siga);
295 } while (0.0 > loss || loss > twomeanLoss);
296 }
297 return loss;
298}
299
300//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
301
302G4double
303G4mplIonisationWithDeltaModel::Dispersion(const G4Material* material,
304 const G4DynamicParticle* dp,
305 G4double& tmax,
306 G4double& length)
307{
308 G4double siga = 0.0;
309 G4double tau = dp->GetKineticEnergy()/mass;
310 if(tau > 0.0) {
311 G4double electronDensity = material->GetElectronDensity();
312 G4double gam = tau + 1.0;
313 G4double invbeta2 = (gam*gam)/(tau * (tau+2.0));
314 siga = (invbeta2 - 0.5) * twopi_mc2_rcl2 * tmax * length
315 * electronDensity * chargeSquare;
316 }
317 return siga;
318}
319
320//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
321
322G4double
323G4mplIonisationWithDeltaModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
324 G4double kinEnergy)
325{
326 G4double tau = kinEnergy/mass;
327 return 2.0*electron_mass_c2*tau*(tau + 2.);
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.