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

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

update

File size: 13.5 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: G4MollerBhabhaModel.cc,v 1.30 2007/05/22 17:34:36 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02 $
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),
76  particle(0),
77  isElectron(true),
78  twoln10(2.0*log(10.0)),
79  lowLimit(0.2*keV)
80{
81  theElectron = G4Electron::Electron();
82  if(p) SetParticle(p);
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86
87G4MollerBhabhaModel::~G4MollerBhabhaModel()
88{}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91
92void G4MollerBhabhaModel::SetParticle(const G4ParticleDefinition* p)
93{
94  particle = p;
95  if(p != theElectron) isElectron = false;
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99
100G4double G4MollerBhabhaModel::MinEnergyCut(const G4ParticleDefinition*,
101                                           const G4MaterialCutsCouple* couple)
102{
103  G4double electronDensity = couple->GetMaterial()->GetElectronDensity();
104  G4double Zeff  = electronDensity/couple->GetMaterial()->GetTotNbOfAtomsPerVolume();
105  return 0.25*sqrt(Zeff)*keV;
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109
110void G4MollerBhabhaModel::Initialise(const G4ParticleDefinition* p,
111                                     const G4DataVector&)
112{
113  if(!particle) SetParticle(p);
114  if(pParticleChange)
115    fParticleChange = reinterpret_cast<G4ParticleChangeForLoss*>
116                                                     (pParticleChange);
117  else
118    fParticleChange = new G4ParticleChangeForLoss();
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
123G4double G4MollerBhabhaModel::ComputeCrossSectionPerElectron(
124                                           const G4ParticleDefinition* p,
125                                                 G4double kineticEnergy,
126                                                 G4double cutEnergy,
127                                                 G4double maxEnergy)
128{
129  if(!particle) SetParticle(p);
130
131  G4double cross = 0.0;
132  G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
133  tmax = min(maxEnergy, tmax);
134
135  if(cutEnergy < tmax) {
136
137    G4double xmin  = cutEnergy/kineticEnergy;
138    G4double xmax  = tmax/kineticEnergy;
139    G4double gam   = kineticEnergy/electron_mass_c2 + 1.0;
140    G4double gamma2= gam*gam;
141    G4double beta2 = 1.0 - 1.0/gamma2;
142
143    //Moller (e-e-) scattering
144    if (isElectron) {
145
146      G4double g     = (2.0*gam - 1.0)/gamma2;
147      cross = ((xmax - xmin)*(1.0 - g + 1.0/(xmin*xmax)
148                              + 1.0/((1.0-xmin)*(1.0 - xmax)))
149            - g*log( xmax*(1.0 - xmin)/(xmin*(1.0 - xmax)) ) ) / beta2;
150
151    //Bhabha (e+e-) scattering
152    } else {
153
154      G4double y   = 1.0/(1.0 + gam);
155      G4double y2  = y*y;
156      G4double y12 = 1.0 - 2.0*y;
157      G4double b1  = 2.0 - y2;
158      G4double b2  = y12*(3.0 + y2);
159      G4double y122= y12*y12;
160      G4double b4  = y122*y12;
161      G4double b3  = b4 + y122;
162
163      cross = (xmax - xmin)*(1.0/(beta2*xmin*xmax) + b2
164            - 0.5*b3*(xmin + xmax)
165            + b4*(xmin*xmin + xmin*xmax + xmax*xmax)/3.0)
166            - b1*log(xmax/xmin);
167    }
168
169    cross *= twopi_mc2_rcl2/kineticEnergy;
170  }
171  return cross;
172}
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175
176G4double G4MollerBhabhaModel::ComputeCrossSectionPerAtom(
177                                           const G4ParticleDefinition* p,
178                                                 G4double kineticEnergy,
179                                                 G4double Z, G4double,
180                                                 G4double cutEnergy,
181                                                 G4double maxEnergy)
182{
183  G4double cross = Z*ComputeCrossSectionPerElectron
184                                         (p,kineticEnergy,cutEnergy,maxEnergy);
185  return cross;
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
189
190G4double G4MollerBhabhaModel::CrossSectionPerVolume(
191                                           const G4Material* material,
192                                           const G4ParticleDefinition* p,
193                                                 G4double kineticEnergy,
194                                                 G4double cutEnergy,
195                                                 G4double maxEnergy)
196{
197  G4double eDensity = material->GetElectronDensity();
198  G4double cross = eDensity*ComputeCrossSectionPerElectron
199                                         (p,kineticEnergy,cutEnergy,maxEnergy);
200  return cross;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
204
205G4double G4MollerBhabhaModel::ComputeDEDXPerVolume(
206                                          const G4Material* material,
207                                          const G4ParticleDefinition* p,
208                                                G4double kineticEnergy,
209                                                G4double cutEnergy)
210{
211  if(!particle) SetParticle(p);
212  // calculate the dE/dx due to the ionization by Seltzer-Berger formula
213 
214  G4double electronDensity = material->GetElectronDensity();
215  G4double Zeff  = electronDensity/material->GetTotNbOfAtomsPerVolume();
216  G4double th    = 0.25*sqrt(Zeff)*keV;
217  G4double tkin  = kineticEnergy;
218  G4bool   lowEnergy = false;
219  if (kineticEnergy < th) {
220    tkin = th;
221    lowEnergy = true;
222  }
223  G4double tau   = tkin/electron_mass_c2;
224  G4double gam   = tau + 1.0;
225  G4double gamma2= gam*gam;
226  G4double beta2 = 1. - 1./gamma2;
227  G4double bg2   = beta2*gamma2;
228
229  G4double eexc  = material->GetIonisation()->GetMeanExcitationEnergy();
230  eexc          /= electron_mass_c2;
231  G4double eexc2 = eexc*eexc; 
232
233  G4double d = min(cutEnergy, MaxSecondaryEnergy(p, tkin))/electron_mass_c2;
234  G4double dedx;
235
236  // electron
237  if (isElectron) {
238
239    dedx = log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
240         + log((tau-d)*d) + tau/(tau-d)
241         + (0.5*d*d + (2.0*tau + 1.)*log(1. - d/tau))/gamma2;
242   
243  //positron
244  } else {
245
246    G4double d2 = d*d*0.5;
247    G4double d3 = d2*d/1.5;
248    G4double d4 = d3*d*3.75;
249    G4double y  = 1.0/(1.0 + gam);
250    dedx = log(2.0*(tau + 2.0)/eexc2) + log(tau*d)
251         - beta2*(tau + 2.0*d - y*(3.0*d2
252         + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
253  } 
254
255  //density correction
256  G4double cden  = material->GetIonisation()->GetCdensity();
257  G4double mden  = material->GetIonisation()->GetMdensity();
258  G4double aden  = material->GetIonisation()->GetAdensity();
259  G4double x0den = material->GetIonisation()->GetX0density();
260  G4double x1den = material->GetIonisation()->GetX1density();
261  G4double x     = log(bg2)/twoln10;
262
263  if (x >= x0den) {
264    dedx -= twoln10*x - cden;
265    if (x < x1den) dedx -= aden*pow(x1den-x, mden);
266  } 
267
268  // now you can compute the total ionization loss
269  dedx *= twopi_mc2_rcl2*electronDensity/beta2;
270  if (dedx < 0.0) dedx = 0.0;
271
272  // lowenergy extrapolation
273
274  if (lowEnergy) {
275
276    if (kineticEnergy >= lowLimit) dedx *= sqrt(tkin/kineticEnergy);
277    else                           dedx *= sqrt(tkin*kineticEnergy)/lowLimit;
278
279  }
280  return dedx;
281}
282
283//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
284
285void G4MollerBhabhaModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
286                                            const G4MaterialCutsCouple*,
287                                            const G4DynamicParticle* dp,
288                                            G4double tmin,
289                                            G4double maxEnergy)
290{
291  G4double tmax = std::min(maxEnergy, MaxSecondaryKinEnergy(dp));
292  if(tmin >= tmax) return;
293
294  G4double kineticEnergy = dp->GetKineticEnergy();
295  G4double energy = kineticEnergy + electron_mass_c2;
296  G4double totalMomentum = sqrt(kineticEnergy*(energy + electron_mass_c2));
297  G4double xmin   = tmin/kineticEnergy;
298  G4double xmax   = tmax/kineticEnergy;
299  G4double gam    = energy/electron_mass_c2;
300  G4double gamma2 = gam*gam;
301  G4double beta2  = 1.0 - 1.0/gamma2;
302  G4double x, z, q, grej;
303
304  G4ThreeVector direction = dp->GetMomentumDirection();
305
306  //Moller (e-e-) scattering
307  if (isElectron) {
308
309    G4double g = (2.0*gam - 1.0)/gamma2;
310    G4double y = 1.0 - xmax;
311    grej = 1.0 - g*xmax + xmax*xmax*(1.0 - g + (1.0 - g*y)/(y*y));
312
313    do {
314      q = G4UniformRand();
315      x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
316      y = 1.0 - x;
317      z = 1.0 - g*x + x*x*(1.0 - g + (1.0 - g*y)/(y*y));
318      /*
319      if(z > grej) {
320        G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
321               << "Majorant " << grej << " < "
322               << z << " for x= " << x
323               << " e-e- scattering"
324               << G4endl;
325      }
326      */
327    } while(grej * G4UniformRand() > z);
328
329  //Bhabha (e+e-) scattering
330  } else {
331
332    G4double y   = 1.0/(1.0 + gam);
333    G4double y2  = y*y;
334    G4double y12 = 1.0 - 2.0*y;
335    G4double b1  = 2.0 - y2;
336    G4double b2  = y12*(3.0 + y2);
337    G4double y122= y12*y12;
338    G4double b4  = y122*y12;
339    G4double b3  = b4 + y122;
340
341    y     = xmax*xmax;
342    grej  = -xmin*b1;
343    grej += y*b2;
344    grej -= xmin*xmin*xmin*b3;
345    grej += y*y*b4;
346    grej *= beta2;
347    grej += 1.0;
348    do {
349      q  = G4UniformRand();
350      x  = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
351      z  = -x*b1;
352      y  = x*x;
353      z += y*b2;
354      y *= x;
355      z -= y*b3;
356      y *= x;
357      z += y*b4;
358      z *= beta2;
359      z += 1.0;
360      /*
361      if(z > grej) {
362        G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
363               << "Majorant " << grej << " < "
364               << z << " for x= " << x
365               << " e+e- scattering"
366               << G4endl;
367      }
368      */
369    } while(grej * G4UniformRand() > z);
370  }
371
372  G4double deltaKinEnergy = x * kineticEnergy;
373
374  G4double deltaMomentum =
375           sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
376  G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
377                                   (deltaMomentum * totalMomentum);
378  G4double sint = 1.0 - cost*cost;
379  if(sint > 0.0) sint = sqrt(sint);
380
381  G4double phi = twopi * G4UniformRand() ;
382
383  G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
384  deltaDirection.rotateUz(direction);
385
386  // primary change
387  kineticEnergy -= deltaKinEnergy;
388  fParticleChange->SetProposedKineticEnergy(kineticEnergy);
389
390  if(kineticEnergy > DBL_MIN) {
391    G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
392    direction = dir.unit();
393    fParticleChange->SetProposedMomentumDirection(direction);
394  }
395
396  // create G4DynamicParticle object for delta ray
397  G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
398                                                 deltaDirection,deltaKinEnergy);
399  vdp->push_back(delta);
400}
401
402//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.