Geant4  10.00.p01
G4BetheHeitlerModel.cc
Go to the documentation of this file.
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 74581 2013-10-15 12:03:25Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4BetheHeitlerModel
34 //
35 // Author: Vladimir Ivanchenko on base of Michel Maire code
36 //
37 // Creation date: 15.03.2005
38 //
39 // Modifications:
40 // 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
41 // 24-06-05 Increase number of bins to 200 (V.Ivantchenko)
42 // 16-11-05 replace shootBit() by G4UniformRand() mma
43 // 04-12-05 SetProposedKineticEnergy(0.) for the killed photon (mma)
44 // 20-02-07 SelectRandomElement is called for any initial gamma energy
45 // in order to have selected element for polarized model (VI)
46 // 25-10-10 Removed unused table, added element selector (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 "G4PhysicalConstants.hh"
57 #include "G4SystemOfUnits.hh"
58 #include "G4Electron.hh"
59 #include "G4Positron.hh"
60 #include "G4Gamma.hh"
61 #include "Randomize.hh"
63 #include "G4Pow.hh"
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
67 using namespace std;
68 
69 static const G4double GammaEnergyLimit = 1.5*MeV;
70 static const G4double Egsmall=2.*MeV;
71 static const G4double
72  a0= 8.7842e+2*microbarn, a1=-1.9625e+3*microbarn, a2= 1.2949e+3*microbarn,
73  a3=-2.0028e+2*microbarn, a4= 1.2575e+1*microbarn, a5=-2.8333e-1*microbarn;
74 
75 static const G4double
76  b0=-1.0342e+1*microbarn, b1= 1.7692e+1*microbarn, b2=-8.2381 *microbarn,
77  b3= 1.3063 *microbarn, b4=-9.0815e-2*microbarn, b5= 2.3586e-3*microbarn;
78 
79 static const G4double
80  c0=-4.5263e+2*microbarn, c1= 1.1161e+3*microbarn, c2=-8.6749e+2*microbarn,
81  c3= 2.1773e+2*microbarn, c4=-2.0467e+1*microbarn, c5= 6.5372e-1*microbarn;
82 
83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84 
86  const G4String& nam)
87  : G4VEmModel(nam)
88 {
89  fParticleChange = 0;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97 
99 {}
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102 
104  const G4DataVector& cuts)
105 {
107  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
111 
113  G4VEmModel* masterModel)
114 {
116 }
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119 
120 G4double
122  G4double GammaEnergy, G4double Z,
124 // Calculates the microscopic cross section in GEANT4 internal units.
125 // A parametrized formula from L. Urban is used to estimate
126 // the total cross section.
127 // It gives a good description of the data from 1.5 MeV to 100 GeV.
128 // below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass)
129 // *(GammaEnergy-2electronmass)
130 {
131  G4double xSection = 0.0 ;
132  if ( Z < 0.9 || GammaEnergy <= 2.0*electron_mass_c2 ) { return xSection; }
133 
134 
135  G4double GammaEnergySave = GammaEnergy;
136  if (GammaEnergy < GammaEnergyLimit) { GammaEnergy = GammaEnergyLimit; }
137 
138  G4double X=G4Log(GammaEnergy/electron_mass_c2), X2=X*X, X3=X2*X, X4=X3*X, X5=X4*X;
139 
140  G4double F1 = a0 + a1*X + a2*X2 + a3*X3 + a4*X4 + a5*X5,
141  F2 = b0 + b1*X + b2*X2 + b3*X3 + b4*X4 + b5*X5,
142  F3 = c0 + c1*X + c2*X2 + c3*X3 + c4*X4 + c5*X5;
143 
144  xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
145 
146  if (GammaEnergySave < GammaEnergyLimit) {
147 
148  X = (GammaEnergySave - 2.*electron_mass_c2)
149  / (GammaEnergyLimit - 2.*electron_mass_c2);
150  xSection *= X*X;
151  }
152 
153  if (xSection < 0.) { xSection = 0.; }
154  return xSection;
155 }
156 
157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158 
159 void G4BetheHeitlerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
160  const G4MaterialCutsCouple* couple,
161  const G4DynamicParticle* aDynamicGamma,
162  G4double,
163  G4double)
164 // The secondaries e+e- energies are sampled using the Bethe - Heitler
165 // cross sections with Coulomb correction.
166 // A modified version of the random number techniques of Butcher & Messel
167 // is used (Nuc Phys 20(1960),15).
168 //
169 // GEANT4 internal units.
170 //
171 // Note 1 : Effects due to the breakdown of the Born approximation at
172 // low energy are ignored.
173 // Note 2 : The differential cross section implicitly takes account of
174 // pair creation in both nuclear and atomic electron fields.
175 // However triplet prodution is not generated.
176 {
177  const G4Material* aMaterial = couple->GetMaterial();
178 
179  G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
180  G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
181 
182  G4double epsil ;
183  G4double epsil0 = electron_mass_c2/GammaEnergy ;
184  if(epsil0 > 1.0) { return; }
185 
186  // do it fast if GammaEnergy < Egsmall
187  // select randomly one element constituing the material
188  const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
189 
190  if (GammaEnergy < Egsmall) {
191 
192  epsil = epsil0 + (0.5-epsil0)*G4UniformRand();
193 
194  } else {
195  // now comes the case with GammaEnergy >= 2. MeV
196 
197  // Extract Coulomb factor for this Element
198  G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
199  if (GammaEnergy > 50.*MeV) { FZ += 8.*(anElement->GetfCoulomb()); }
200 
201  // limits of the screening variable
202  G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
203  G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
204  G4double screenmin = min(4.*screenfac,screenmax);
205 
206  // limits of the energy sampling
207  G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
208  G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
209 
210  //
211  // sample the energy rate of the created electron (or positron)
212  //
213  //G4double epsil, screenvar, greject ;
214  G4double screenvar, greject ;
215 
216  G4double F10 = ScreenFunction1(screenmin) - FZ;
217  G4double F20 = ScreenFunction2(screenmin) - FZ;
218  G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
219  G4double NormF2 = max(1.5*F20,0.);
220 
221  do {
222  if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) {
223  epsil = 0.5 - epsilrange*g4pow->A13(G4UniformRand());
224  screenvar = screenfac/(epsil*(1-epsil));
225  greject = (ScreenFunction1(screenvar) - FZ)/F10;
226 
227  } else {
228  epsil = epsilmin + epsilrange*G4UniformRand();
229  screenvar = screenfac/(epsil*(1-epsil));
230  greject = (ScreenFunction2(screenvar) - FZ)/F20;
231  }
232 
233  } while( greject < G4UniformRand() );
234 
235  } // end of epsil sampling
236 
237  //
238  // fixe charges randomly
239  //
240 
241  G4double ElectTotEnergy, PositTotEnergy;
242  if (G4UniformRand() > 0.5) {
243 
244  ElectTotEnergy = (1.-epsil)*GammaEnergy;
245  PositTotEnergy = epsil*GammaEnergy;
246 
247  } else {
248 
249  PositTotEnergy = (1.-epsil)*GammaEnergy;
250  ElectTotEnergy = epsil*GammaEnergy;
251  }
252 
253  //
254  // scattered electron (positron) angles. ( Z - axis along the parent photon)
255  //
256  // universal distribution suggested by L. Urban
257  // (Geant3 manual (1993) Phys211),
258  // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
259 
260  G4double u;
261  static const G4double aa1 = 0.625;
262  static const G4double aa2 = 1.875;
263  static const G4double d = 27. ;
264 
265  if (9./(9.+d) >G4UniformRand()) u= - G4Log(G4UniformRand()*G4UniformRand())/aa1;
266  else u= - G4Log(G4UniformRand()*G4UniformRand())/aa2;
267 
268  G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
269  G4double TetPo = u*electron_mass_c2/PositTotEnergy;
270  G4double Phi = twopi * G4UniformRand();
271  G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
272  G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
273 
274  //
275  // kinematic of the created pair
276  //
277  // the electron and positron are assumed to have a symetric
278  // angular distribution with respect to the Z axis along the parent photon.
279 
280  G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
281 
282  G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
283  ElectDirection.rotateUz(GammaDirection);
284 
285  // create G4DynamicParticle object for the particle1
286  G4DynamicParticle* aParticle1= new G4DynamicParticle(
287  theElectron,ElectDirection,ElectKineEnergy);
288 
289  // the e+ is always created (even with Ekine=0) for further annihilation.
290 
291  G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
292 
293  G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
294  PositDirection.rotateUz(GammaDirection);
295 
296  // create G4DynamicParticle object for the particle2
297  G4DynamicParticle* aParticle2= new G4DynamicParticle(
298  thePositron,PositDirection,PositKineEnergy);
299 
300  // Fill output vector
301  fvect->push_back(aParticle1);
302  fvect->push_back(aParticle2);
303 
304  // kill incident photon
307 }
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4Pow * GetInstance()
Definition: G4Pow.cc:53
static const double MeV
Definition: G4SIunits.hh:193
static const G4double a0
G4ParticleChangeForGamma * fParticleChange
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
static const G4double c2
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static const G4double a1
static const double microbarn
Definition: G4SIunits.hh:97
G4ParticleDefinition * thePositron
static const G4double a4
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4double ScreenFunction2(G4double ScreenVariable)
#define F20
static const G4double Egsmall
#define G4UniformRand()
Definition: Randomize.hh:87
static const G4double b3
G4double ScreenFunction1(G4double ScreenVariable)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
const G4ThreeVector & GetMomentumDirection() const
static const G4double c3
static const G4double c4
static const G4double b2
G4ParticleDefinition * theGamma
#define F10
static const G4double c1
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
G4BetheHeitlerModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitler")
G4double GetlogZ3() const
static const G4double a3
G4double G4Log(G4double x)
Definition: G4Log.hh:227
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4ParticleDefinition * theElectron
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
G4double A13(G4double A) const
Definition: G4Pow.hh:134
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static const G4double c0
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:198
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const G4double b1
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static const G4double c5
static const G4double a5
static const G4double GammaEnergyLimit
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ThreeVector G4ParticleMomentum
static const G4double b0
static const G4double b5
static const G4double b4
G4double GetZ3() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
static const G4double a2
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
const G4Material * GetMaterial() const