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

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

update geant4.9.3 tag

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