source: trunk/source/processes/electromagnetic/standard/src/G4MollerBhabhaModel.cc@ 1199

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

update CVS release candidate geant4.9.3.01

File size: 13.5 KB
RevLine 
[819]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//
[1196]26// $Id: G4MollerBhabhaModel.cc,v 1.35 2009/11/09 19:16:13 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
[819]28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4MollerBhabhaModel
35//
36// Author: Vladimir Ivanchenko on base of Laszlo Urban code
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 13-11-02 Minor fix - use normalised direction (V.Ivanchenko)
43// 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (V.Ivanchenko)
44// 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
45// 27-01-03 Make models region aware (V.Ivanchenko)
46// 13-02-03 Add name (V.Ivanchenko)
47// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
48// 25-07-05 Add protection in calculation of recoil direction for the case
49// of complete energy transfer from e+ to e- (V.Ivanchenko)
50// 06-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
51// 15-05-06 Fix MinEnergyCut (V.Ivanchenko)
52//
53//
54// Class Description:
55//
56// Implementation of energy loss and delta-electron production by e+/e-
57//
58// -------------------------------------------------------------------
59//
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62
63#include "G4MollerBhabhaModel.hh"
64#include "G4Electron.hh"
65#include "G4Positron.hh"
66#include "Randomize.hh"
67#include "G4ParticleChangeForLoss.hh"
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71using namespace std;
72
73G4MollerBhabhaModel::G4MollerBhabhaModel(const G4ParticleDefinition* p,
74 const G4String& nam)
75 : G4VEmModel(nam),
[1055]76 particle(0),
77 isElectron(true),
78 twoln10(2.0*log(10.0)),
79 lowLimit(0.2*keV),
80 isInitialised(false)
[819]81{
82 theElectron = G4Electron::Electron();
83 if(p) SetParticle(p);
84}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
88G4MollerBhabhaModel::~G4MollerBhabhaModel()
89{}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
93G4double G4MollerBhabhaModel::MinEnergyCut(const G4ParticleDefinition*,
94 const G4MaterialCutsCouple* couple)
95{
96 G4double electronDensity = couple->GetMaterial()->GetElectronDensity();
97 G4double Zeff = electronDensity/couple->GetMaterial()->GetTotNbOfAtomsPerVolume();
98 return 0.25*sqrt(Zeff)*keV;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
[1055]103G4double G4MollerBhabhaModel::MaxSecondaryEnergy(const G4ParticleDefinition*,
104 G4double kinEnergy)
105{
106 G4double tmax = kinEnergy;
107 if(isElectron) tmax *= 0.5;
108 return tmax;
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
[819]113void G4MollerBhabhaModel::Initialise(const G4ParticleDefinition* p,
114 const G4DataVector&)
115{
116 if(!particle) SetParticle(p);
[1055]117 SetDeexcitationFlag(false);
118
119 if(isInitialised) return;
120
121 isInitialised = true;
122 fParticleChange = GetParticleChangeForLoss();
[819]123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
126
127G4double G4MollerBhabhaModel::ComputeCrossSectionPerElectron(
128 const G4ParticleDefinition* p,
129 G4double kineticEnergy,
130 G4double cutEnergy,
131 G4double maxEnergy)
132{
133 if(!particle) SetParticle(p);
134
135 G4double cross = 0.0;
136 G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
137 tmax = min(maxEnergy, tmax);
138
139 if(cutEnergy < tmax) {
140
141 G4double xmin = cutEnergy/kineticEnergy;
142 G4double xmax = tmax/kineticEnergy;
143 G4double gam = kineticEnergy/electron_mass_c2 + 1.0;
144 G4double gamma2= gam*gam;
145 G4double beta2 = 1.0 - 1.0/gamma2;
146
147 //Moller (e-e-) scattering
148 if (isElectron) {
149
150 G4double g = (2.0*gam - 1.0)/gamma2;
151 cross = ((xmax - xmin)*(1.0 - g + 1.0/(xmin*xmax)
152 + 1.0/((1.0-xmin)*(1.0 - xmax)))
153 - g*log( xmax*(1.0 - xmin)/(xmin*(1.0 - xmax)) ) ) / beta2;
154
155 //Bhabha (e+e-) scattering
156 } else {
157
158 G4double y = 1.0/(1.0 + gam);
159 G4double y2 = y*y;
160 G4double y12 = 1.0 - 2.0*y;
161 G4double b1 = 2.0 - y2;
162 G4double b2 = y12*(3.0 + y2);
163 G4double y122= y12*y12;
164 G4double b4 = y122*y12;
165 G4double b3 = b4 + y122;
166
167 cross = (xmax - xmin)*(1.0/(beta2*xmin*xmax) + b2
168 - 0.5*b3*(xmin + xmax)
169 + b4*(xmin*xmin + xmin*xmax + xmax*xmax)/3.0)
170 - b1*log(xmax/xmin);
171 }
172
173 cross *= twopi_mc2_rcl2/kineticEnergy;
174 }
175 return cross;
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179
180G4double G4MollerBhabhaModel::ComputeCrossSectionPerAtom(
181 const G4ParticleDefinition* p,
182 G4double kineticEnergy,
183 G4double Z, G4double,
184 G4double cutEnergy,
185 G4double maxEnergy)
186{
187 G4double cross = Z*ComputeCrossSectionPerElectron
188 (p,kineticEnergy,cutEnergy,maxEnergy);
189 return cross;
190}
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193
194G4double G4MollerBhabhaModel::CrossSectionPerVolume(
195 const G4Material* material,
196 const G4ParticleDefinition* p,
197 G4double kineticEnergy,
198 G4double cutEnergy,
199 G4double maxEnergy)
200{
201 G4double eDensity = material->GetElectronDensity();
202 G4double cross = eDensity*ComputeCrossSectionPerElectron
203 (p,kineticEnergy,cutEnergy,maxEnergy);
204 return cross;
205}
206
207//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
208
209G4double G4MollerBhabhaModel::ComputeDEDXPerVolume(
210 const G4Material* material,
211 const G4ParticleDefinition* p,
212 G4double kineticEnergy,
213 G4double cutEnergy)
214{
215 if(!particle) SetParticle(p);
216 // calculate the dE/dx due to the ionization by Seltzer-Berger formula
217
218 G4double electronDensity = material->GetElectronDensity();
219 G4double Zeff = electronDensity/material->GetTotNbOfAtomsPerVolume();
220 G4double th = 0.25*sqrt(Zeff)*keV;
221 G4double tkin = kineticEnergy;
222 G4bool lowEnergy = false;
223 if (kineticEnergy < th) {
224 tkin = th;
225 lowEnergy = true;
226 }
227 G4double tau = tkin/electron_mass_c2;
228 G4double gam = tau + 1.0;
229 G4double gamma2= gam*gam;
230 G4double beta2 = 1. - 1./gamma2;
231 G4double bg2 = beta2*gamma2;
232
233 G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
234 eexc /= electron_mass_c2;
235 G4double eexc2 = eexc*eexc;
236
237 G4double d = min(cutEnergy, MaxSecondaryEnergy(p, tkin))/electron_mass_c2;
238 G4double dedx;
239
240 // electron
241 if (isElectron) {
242
243 dedx = log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
244 + log((tau-d)*d) + tau/(tau-d)
245 + (0.5*d*d + (2.0*tau + 1.)*log(1. - d/tau))/gamma2;
246
247 //positron
248 } else {
249
250 G4double d2 = d*d*0.5;
251 G4double d3 = d2*d/1.5;
252 G4double d4 = d3*d*3.75;
253 G4double y = 1.0/(1.0 + gam);
254 dedx = log(2.0*(tau + 2.0)/eexc2) + log(tau*d)
255 - beta2*(tau + 2.0*d - y*(3.0*d2
256 + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
257 }
258
259 //density correction
[1196]260 //G4double cden = material->GetIonisation()->GetCdensity();
261 //G4double mden = material->GetIonisation()->GetMdensity();
262 //G4double aden = material->GetIonisation()->GetAdensity();
263 //G4double x0den = material->GetIonisation()->GetX0density();
264 //G4double x1den = material->GetIonisation()->GetX1density();
[819]265 G4double x = log(bg2)/twoln10;
266
[1196]267 //if (x >= x0den) {
268 // dedx -= twoln10*x - cden;
269 // if (x < x1den) dedx -= aden*pow(x1den-x, mden);
270 //}
271 dedx -= material->GetIonisation()->DensityCorrection(x);
[819]272
273 // now you can compute the total ionization loss
274 dedx *= twopi_mc2_rcl2*electronDensity/beta2;
275 if (dedx < 0.0) dedx = 0.0;
276
277 // lowenergy extrapolation
278
279 if (lowEnergy) {
280
281 if (kineticEnergy >= lowLimit) dedx *= sqrt(tkin/kineticEnergy);
282 else dedx *= sqrt(tkin*kineticEnergy)/lowLimit;
283
284 }
285 return dedx;
286}
287
288//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
289
290void G4MollerBhabhaModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
291 const G4MaterialCutsCouple*,
292 const G4DynamicParticle* dp,
293 G4double tmin,
294 G4double maxEnergy)
295{
[1055]296 G4double kineticEnergy = dp->GetKineticEnergy();
297 G4double tmax = kineticEnergy;
298 if(isElectron) tmax *= 0.5;
299 if(maxEnergy < tmax) tmax = maxEnergy;
[819]300 if(tmin >= tmax) return;
301
302 G4double energy = kineticEnergy + electron_mass_c2;
303 G4double totalMomentum = sqrt(kineticEnergy*(energy + electron_mass_c2));
304 G4double xmin = tmin/kineticEnergy;
305 G4double xmax = tmax/kineticEnergy;
306 G4double gam = energy/electron_mass_c2;
307 G4double gamma2 = gam*gam;
308 G4double beta2 = 1.0 - 1.0/gamma2;
309 G4double x, z, q, grej;
310
311 G4ThreeVector direction = dp->GetMomentumDirection();
312
313 //Moller (e-e-) scattering
314 if (isElectron) {
315
316 G4double g = (2.0*gam - 1.0)/gamma2;
317 G4double y = 1.0 - xmax;
318 grej = 1.0 - g*xmax + xmax*xmax*(1.0 - g + (1.0 - g*y)/(y*y));
319
320 do {
321 q = G4UniformRand();
322 x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
323 y = 1.0 - x;
324 z = 1.0 - g*x + x*x*(1.0 - g + (1.0 - g*y)/(y*y));
325 /*
326 if(z > grej) {
327 G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
328 << "Majorant " << grej << " < "
329 << z << " for x= " << x
330 << " e-e- scattering"
331 << G4endl;
332 }
333 */
334 } while(grej * G4UniformRand() > z);
335
336 //Bhabha (e+e-) scattering
337 } else {
338
339 G4double y = 1.0/(1.0 + gam);
340 G4double y2 = y*y;
341 G4double y12 = 1.0 - 2.0*y;
342 G4double b1 = 2.0 - y2;
343 G4double b2 = y12*(3.0 + y2);
344 G4double y122= y12*y12;
345 G4double b4 = y122*y12;
346 G4double b3 = b4 + y122;
347
[1055]348 y = xmax*xmax;
349 grej = 1.0 + (y*y*b4 - xmin*xmin*xmin*b3 + y*b2 - xmin*b1)*beta2;
[819]350 do {
[1055]351 q = G4UniformRand();
352 x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
353 y = x*x;
354 z = 1.0 + (y*y*b4 - x*y*b3 + y*b2 - x*b1)*beta2;
[819]355 /*
356 if(z > grej) {
357 G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
358 << "Majorant " << grej << " < "
359 << z << " for x= " << x
360 << " e+e- scattering"
361 << G4endl;
362 }
363 */
364 } while(grej * G4UniformRand() > z);
365 }
366
367 G4double deltaKinEnergy = x * kineticEnergy;
368
369 G4double deltaMomentum =
370 sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
371 G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
372 (deltaMomentum * totalMomentum);
373 G4double sint = 1.0 - cost*cost;
374 if(sint > 0.0) sint = sqrt(sint);
375
376 G4double phi = twopi * G4UniformRand() ;
377
378 G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
379 deltaDirection.rotateUz(direction);
380
381 // primary change
382 kineticEnergy -= deltaKinEnergy;
383 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
384
385 if(kineticEnergy > DBL_MIN) {
386 G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
387 direction = dir.unit();
388 fParticleChange->SetProposedMomentumDirection(direction);
389 }
390
391 // create G4DynamicParticle object for delta ray
392 G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
[1055]393 deltaDirection,deltaKinEnergy);
[819]394 vdp->push_back(delta);
395}
396
397//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.