Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 100399 2016-10-20 07:38:12Z 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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70 
72  const G4String& nam)
73  : G4VEmModel(nam)
74 {
75  fParticleChange = nullptr;
76  theGamma = G4Gamma::Gamma();
77  thePositron = G4Positron::Positron();
78  theElectron = G4Electron::Electron();
79  g4calc = G4Pow::GetInstance();
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
85 {}
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
90  const G4DataVector& cuts)
91 {
92  if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
93  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
97 
99  G4VEmModel* masterModel)
100 {
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
106 G4double
108  G4double GammaEnergy, G4double Z,
110 // Calculates the microscopic cross section in GEANT4 internal units.
111 // A parametrized formula from L. Urban is used to estimate
112 // the total cross section.
113 // It gives a good description of the data from 1.5 MeV to 100 GeV.
114 // below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass)
115 // *(GammaEnergy-2electronmass)
116 {
117  G4double xSection = 0.0 ;
118  if ( Z < 0.9 || GammaEnergy <= 2.0*electron_mass_c2 ) { return xSection; }
119 
120  static const G4double GammaEnergyLimit = 1.5*MeV;
121  static const G4double
122  a0= 8.7842e+2*microbarn, a1=-1.9625e+3*microbarn, a2= 1.2949e+3*microbarn,
123  a3=-2.0028e+2*microbarn, a4= 1.2575e+1*microbarn, a5=-2.8333e-1*microbarn;
124 
125  static const G4double
126  b0=-1.0342e+1*microbarn, b1= 1.7692e+1*microbarn, b2=-8.2381 *microbarn,
127  b3= 1.3063 *microbarn, b4=-9.0815e-2*microbarn, b5= 2.3586e-3*microbarn;
128 
129  static const G4double
130  c0=-4.5263e+2*microbarn, c1= 1.1161e+3*microbarn, c2=-8.6749e+2*microbarn,
131  c3= 2.1773e+2*microbarn, c4=-2.0467e+1*microbarn, c5= 6.5372e-1*microbarn;
132 
133  G4double GammaEnergySave = GammaEnergy;
134  if (GammaEnergy < GammaEnergyLimit) { GammaEnergy = GammaEnergyLimit; }
135 
136  G4double X=G4Log(GammaEnergy/electron_mass_c2), X2=X*X, X3=X2*X, X4=X3*X, X5=X4*X;
137 
138  G4double F1 = a0 + a1*X + a2*X2 + a3*X3 + a4*X4 + a5*X5,
139  F2 = b0 + b1*X + b2*X2 + b3*X3 + b4*X4 + b5*X5,
140  F3 = c0 + c1*X + c2*X2 + c3*X3 + c4*X4 + c5*X5;
141 
142  xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
143 
144  if (GammaEnergySave < GammaEnergyLimit) {
145 
146  X = (GammaEnergySave - 2.*electron_mass_c2)
147  / (GammaEnergyLimit - 2.*electron_mass_c2);
148  xSection *= X*X;
149  }
150 
151  xSection = std::max(xSection, 0.);
152  return xSection;
153 }
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156 
157 void G4BetheHeitlerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
158  const G4MaterialCutsCouple* couple,
159  const G4DynamicParticle* aDynamicGamma,
160  G4double,
161  G4double)
162 // The secondaries e+e- energies are sampled using the Bethe - Heitler
163 // cross sections with Coulomb correction.
164 // A modified version of the random number techniques of Butcher & Messel
165 // is used (Nuc Phys 20(1960),15).
166 //
167 // GEANT4 internal units.
168 //
169 // Note 1 : Effects due to the breakdown of the Born approximation at
170 // low energy are ignored.
171 // Note 2 : The differential cross section implicitly takes account of
172 // pair creation in both nuclear and atomic electron fields.
173 // However triplet prodution is not generated.
174 {
175  const G4Material* aMaterial = couple->GetMaterial();
176 
177  G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
178  G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
179 
180  G4double epsil ;
181  G4double epsil0 = electron_mass_c2/GammaEnergy ;
182  if(epsil0 > 1.0) { return; }
183 
184  // do it fast if GammaEnergy < Egsmall
185  // select randomly one element constituing the material
186  const G4Element* anElement =
187  SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
188 
189  CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
190 
191  static const G4double Egsmall=2.*MeV;
192  if (GammaEnergy < Egsmall) {
193 
194  epsil = epsil0 + (0.5-epsil0)*rndmEngine->flat();
195 
196  } else {
197  // now comes the case with GammaEnergy >= 2. MeV
198 
199  // Extract Coulomb factor for this Element
200  G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
201  if (GammaEnergy > 50.*MeV) { FZ += 8.*(anElement->GetfCoulomb()); }
202 
203  // limits of the screening variable
204  G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
205  G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
206  G4double screenmin = min(4.*screenfac,screenmax);
207 
208  // limits of the energy sampling
209  G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
210  G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
211 
212  //
213  // sample the energy rate of the created electron (or positron)
214  //
215  //G4double epsil, screenvar, greject ;
216  G4double screenvar, greject ;
217 
218  G4double F10 = ScreenFunction1(screenmin) - FZ;
219  G4double F20 = ScreenFunction2(screenmin) - FZ;
220  G4double NormF1 = max(F10*epsilrange*epsilrange,0.);
221  G4double NormF2 = max(1.5*F20,0.);
222 
223  do {
224  if ( NormF1/(NormF1+NormF2) > rndmEngine->flat()) {
225  epsil = 0.5 - epsilrange*g4calc->A13(rndmEngine->flat());
226  screenvar = screenfac/(epsil*(1-epsil));
227  greject = (ScreenFunction1(screenvar) - FZ)/F10;
228 
229  } else {
230  epsil = epsilmin + epsilrange*rndmEngine->flat();
231  screenvar = screenfac/(epsil*(1-epsil));
232  greject = (ScreenFunction2(screenvar) - FZ)/F20;
233  }
234 
235  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
236  } while( greject < rndmEngine->flat());
237 
238  } // end of epsil sampling
239 
240  //
241  // fixe charges randomly
242  //
243 
244  G4double ElectTotEnergy, PositTotEnergy;
245  if (rndmEngine->flat() > 0.5) {
246 
247  ElectTotEnergy = (1.-epsil)*GammaEnergy;
248  PositTotEnergy = epsil*GammaEnergy;
249 
250  } else {
251 
252  PositTotEnergy = (1.-epsil)*GammaEnergy;
253  ElectTotEnergy = epsil*GammaEnergy;
254  }
255 
256  //
257  // scattered electron (positron) angles. ( Z - axis along the parent photon)
258  //
259  // universal distribution suggested by L. Urban
260  // (Geant3 manual (1993) Phys211),
261  // derived from Tsai distribution (Rev Mod Phys 49,421(1977))
262 
263  G4double u= - G4Log(rndmEngine->flat()*rndmEngine->flat());
264 
265  if (9. > 36.*rndmEngine->flat()) { u *= 1.6; }
266  else { u *= 0.53333; }
267 
268  G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
269  G4double TetPo = u*electron_mass_c2/PositTotEnergy;
270  G4double Phi = twopi * rndmEngine->flat();
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
305  fParticleChange->SetProposedKineticEnergy(0.);
306  fParticleChange->ProposeTrackStatus(fStopAndKill);
307 }
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
const G4double a0
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double GetKineticEnergy() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
const char * p
Definition: xmltok.h:285
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
G4double GetfCoulomb() const
Definition: G4Element.hh:191
virtual double flat()=0
#define F20
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cut=0., G4double emax=DBL_MAX) override
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
#define F10
double flat()
Definition: G4AblaRandom.cc:47
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
G4BetheHeitlerModel(const G4ParticleDefinition *p=0, const G4String &nam="BetheHeitler")
G4double GetlogZ3() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809
G4double A13(G4double A) const
Definition: G4Pow.hh:132
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:199
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double GetZ3() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541
static constexpr double microbarn
Definition: G4SIunits.hh:107
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const