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

Last change on this file since 1316 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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