source: trunk/source/processes/electromagnetic/standard/src/G4BetheHeitlerModel.cc @ 846

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 11.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: G4BetheHeitlerModel.cc,v 1.11 2007/05/22 17:34:36 vnivanch Exp $
27// GEANT4 tag $Name:  $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4BetheHeitlerModel
35//
36// Author:        Vladimir Ivanchenko on base of Michel Maire code
37//
38// Creation date: 15.03.2005
39//
40// Modifications:
41// 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
42// 24-06-05 Increase number of bins to 200 (V.Ivantchenko)
43// 16-11-05 replace shootBit() by G4UniformRand()  mma
44// 04-12-05 SetProposedKineticEnergy(0.) for the killed photon (mma)
45// 20-02-20 SelectRandomElement is called for any initial gamma energy
46//          in order to have selected element for polarized model (VI)
47//
48// Class Description:
49//
50// -------------------------------------------------------------------
51//
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
55#include "G4BetheHeitlerModel.hh"
56#include "G4Electron.hh"
57#include "G4Positron.hh"
58#include "G4Gamma.hh"
59#include "Randomize.hh"
60#include "G4DataVector.hh"
61#include "G4PhysicsLogVector.hh"
62#include "G4ParticleChangeForGamma.hh"
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66using namespace std;
67
68G4BetheHeitlerModel::G4BetheHeitlerModel(const G4ParticleDefinition*,
69                                         const G4String& nam)
70  : G4VEmModel(nam),
71    theCrossSectionTable(0),
72    nbins(200)
73{
74  theGamma    = G4Gamma::Gamma();
75  thePositron = G4Positron::Positron();
76  theElectron = G4Electron::Electron();
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
81G4BetheHeitlerModel::~G4BetheHeitlerModel()
82{
83  if(theCrossSectionTable) {
84    theCrossSectionTable->clearAndDestroy();
85    delete theCrossSectionTable;
86  }
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90
91void G4BetheHeitlerModel::Initialise(const G4ParticleDefinition*,
92                                     const G4DataVector&)
93{
94  if(pParticleChange)
95    fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
96  else
97    fParticleChange = new G4ParticleChangeForGamma();
98
99  if(theCrossSectionTable) {
100    theCrossSectionTable->clearAndDestroy();
101    delete theCrossSectionTable;
102  }
103
104  const G4ElementTable* theElementTable = G4Element::GetElementTable();
105  size_t nvect = G4Element::GetNumberOfElements();
106  theCrossSectionTable = new G4PhysicsTable(nvect);
107  G4PhysicsLogVector* ptrVector;
108  G4double emin = LowEnergyLimit();
109  G4double emax = HighEnergyLimit();
110  G4double e, value;
111
112  for(size_t j=0; j<nvect ; j++) { 
113
114    ptrVector  = new G4PhysicsLogVector(emin, emax, nbins);
115    G4double Z = (*theElementTable)[j]->GetZ();
116    G4int   iz = G4int(Z);
117    indexZ[iz] = j;
118 
119    for(G4int i=0; i<nbins; i++) {
120      e = ptrVector->GetLowEdgeEnergy( i ) ;
121      value = ComputeCrossSectionPerAtom(theGamma, e, Z); 
122      ptrVector->PutValue( i, value );
123    }
124
125    theCrossSectionTable->insert(ptrVector);
126  }
127}
128
129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130
131G4double G4BetheHeitlerModel::ComputeCrossSectionPerAtom(
132                                                   const G4ParticleDefinition*,
133                                              G4double GammaEnergy, G4double Z,
134                                              G4double, G4double, G4double)
135// Calculates the microscopic cross section in GEANT4 internal units.
136// A parametrized formula from L. Urban is used to estimate
137// the total cross section.
138// It gives a good description of the data from 1.5 MeV to 100 GeV.
139// below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass)
140//                                   *(GammaEnergy-2electronmass)
141{
142  static const G4double GammaEnergyLimit = 1.5*MeV;
143  G4double CrossSection = 0.0 ;
144  if ( Z < 1. ) return CrossSection;
145  if ( GammaEnergy <= 2.0*electron_mass_c2 ) return CrossSection;
146
147  static const G4double
148    a0= 8.7842e+2*microbarn, a1=-1.9625e+3*microbarn, a2= 1.2949e+3*microbarn,
149    a3=-2.0028e+2*microbarn, a4= 1.2575e+1*microbarn, a5=-2.8333e-1*microbarn;
150
151  static const G4double
152    b0=-1.0342e+1*microbarn, b1= 1.7692e+1*microbarn, b2=-8.2381   *microbarn,
153    b3= 1.3063   *microbarn, b4=-9.0815e-2*microbarn, b5= 2.3586e-3*microbarn;
154
155  static const G4double
156    c0=-4.5263e+2*microbarn, c1= 1.1161e+3*microbarn, c2=-8.6749e+2*microbarn,
157    c3= 2.1773e+2*microbarn, c4=-2.0467e+1*microbarn, c5= 6.5372e-1*microbarn;
158
159  G4double GammaEnergySave = GammaEnergy;
160  if (GammaEnergy < GammaEnergyLimit) GammaEnergy = GammaEnergyLimit ;
161
162  G4double X=log(GammaEnergy/electron_mass_c2), X2=X*X, X3=X2*X, X4=X3*X, X5=X4*X;
163
164  G4double F1 = a0 + a1*X + a2*X2 + a3*X3 + a4*X4 + a5*X5,
165           F2 = b0 + b1*X + b2*X2 + b3*X3 + b4*X4 + b5*X5,
166           F3 = c0 + c1*X + c2*X2 + c3*X3 + c4*X4 + c5*X5;     
167
168  CrossSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
169
170  if (GammaEnergySave < GammaEnergyLimit) {
171
172    X = (GammaEnergySave  - 2.*electron_mass_c2)
173      / (GammaEnergyLimit - 2.*electron_mass_c2);
174    CrossSection *= X*X;
175  }
176
177  if (CrossSection < 0.) CrossSection = 0.;
178  return CrossSection;
179}
180
181//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
182
183void G4BetheHeitlerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
184                                            const G4MaterialCutsCouple* couple,
185                                            const G4DynamicParticle* aDynamicGamma,
186                                            G4double,
187                                            G4double)
188// The secondaries e+e- energies are sampled using the Bethe - Heitler
189// cross sections with Coulomb correction.
190// A modified version of the random number techniques of Butcher & Messel
191// is used (Nuc Phys 20(1960),15).
192//
193// GEANT4 internal units.
194//
195// Note 1 : Effects due to the breakdown of the Born approximation at
196//          low energy are ignored.
197// Note 2 : The differential cross section implicitly takes account of
198//          pair creation in both nuclear and atomic electron fields.
199//          However triplet prodution is not generated.
200{
201  const G4Material* aMaterial = couple->GetMaterial();
202
203  G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
204  G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
205
206  G4double epsil ;
207  G4double epsil0 = electron_mass_c2/GammaEnergy ;
208  if(epsil0 > 1.0) return;
209
210  // do it fast if GammaEnergy < 2. MeV
211  static const G4double Egsmall=2.*MeV;
212
213  // select randomly one element constituing the material
214  const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
215
216  if (GammaEnergy < Egsmall) {
217
218    epsil = epsil0 + (0.5-epsil0)*G4UniformRand();
219
220  } else {
221    // now comes the case with GammaEnergy >= 2. MeV
222
223    // Extract Coulomb factor for this Element
224    G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
225    if (GammaEnergy > 50.*MeV) FZ += 8.*(anElement->GetfCoulomb());
226
227    // limits of the screening variable
228    G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
229    G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
230    G4double screenmin = min(4.*screenfac,screenmax);
231
232    // limits of the energy sampling
233    G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
234    G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
235
236    //
237    // sample the energy rate of the created electron (or positron)
238    //
239    //G4double epsil, screenvar, greject ;
240    G4double  screenvar, greject ;
241
242    G4double F10 = ScreenFunction1(screenmin) - FZ;
243    G4double F20 = ScreenFunction2(screenmin) - FZ;
244    G4double NormF1 = max(F10*epsilrange*epsilrange,0.); 
245    G4double NormF2 = max(1.5*F20,0.);
246
247    do {
248      if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) {
249        epsil = 0.5 - epsilrange*pow(G4UniformRand(), 0.333333);
250        screenvar = screenfac/(epsil*(1-epsil));
251        greject = (ScreenFunction1(screenvar) - FZ)/F10;
252             
253      } else { 
254        epsil = epsilmin + epsilrange*G4UniformRand();
255        screenvar = screenfac/(epsil*(1-epsil));
256        greject = (ScreenFunction2(screenvar) - FZ)/F20;
257      }
258
259    } while( greject < G4UniformRand() );
260
261  }   //  end of epsil sampling
262   
263  //
264  // fixe charges randomly
265  //
266
267  G4double ElectTotEnergy, PositTotEnergy;
268  if (G4UniformRand() > 0.5) {
269
270    ElectTotEnergy = (1.-epsil)*GammaEnergy;
271    PositTotEnergy = epsil*GammaEnergy;
272     
273  } else {
274   
275    PositTotEnergy = (1.-epsil)*GammaEnergy;
276    ElectTotEnergy = epsil*GammaEnergy;
277  }
278
279  //
280  // scattered electron (positron) angles. ( Z - axis along the parent photon)
281  //
282  //  universal distribution suggested by L. Urban
283  // (Geant3 manual (1993) Phys211),
284  //  derived from Tsai distribution (Rev Mod Phys 49,421(1977))
285
286  G4double u;
287  const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
288
289  if (9./(9.+d) >G4UniformRand()) u= - log(G4UniformRand()*G4UniformRand())/a1;
290  else                            u= - log(G4UniformRand()*G4UniformRand())/a2;
291
292  G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
293  G4double TetPo = u*electron_mass_c2/PositTotEnergy;
294  G4double Phi  = twopi * G4UniformRand();
295  G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
296  G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
297   
298  //
299  // kinematic of the created pair
300  //
301  // the electron and positron are assumed to have a symetric
302  // angular distribution with respect to the Z axis along the parent photon.
303
304  G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
305
306  G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
307  ElectDirection.rotateUz(GammaDirection);   
308
309  // create G4DynamicParticle object for the particle1 
310  G4DynamicParticle* aParticle1= new G4DynamicParticle(
311                     theElectron,ElectDirection,ElectKineEnergy);
312 
313  // the e+ is always created (even with Ekine=0) for further annihilation.
314
315  G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
316
317  G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
318  PositDirection.rotateUz(GammaDirection);   
319
320  // create G4DynamicParticle object for the particle2
321  G4DynamicParticle* aParticle2= new G4DynamicParticle(
322                      thePositron,PositDirection,PositKineEnergy);
323
324  // Fill output vector
325  fvect->push_back(aParticle1);
326  fvect->push_back(aParticle2);
327
328  // kill incident photon
329  fParticleChange->SetProposedKineticEnergy(0.);
330  fParticleChange->ProposeTrackStatus(fStopAndKill);   
331}
332
333//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
334
335
Note: See TracBrowser for help on using the repository browser.