Geant4  10.01.p01
G4KleinNishinaModel.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: G4KleinNishinaModel.cc 88979 2015-03-17 10:10:21Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4KleinNishinaModel
34 //
35 // Author: Vladimir Ivanchenko on base of G4KleinNishinaCompton
36 //
37 // Creation date: 13.06.2010
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // -------------------------------------------------------------------
44 //
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47 
48 #include "G4KleinNishinaModel.hh"
49 #include "G4PhysicalConstants.hh"
50 #include "G4SystemOfUnits.hh"
51 #include "G4Electron.hh"
52 #include "G4Gamma.hh"
53 #include "Randomize.hh"
54 #include "G4RandomDirection.hh"
55 #include "G4DataVector.hh"
57 #include "G4VAtomDeexcitation.hh"
58 #include "G4AtomicShells.hh"
59 #include "G4LossTableManager.hh"
60 #include "G4Log.hh"
61 #include "G4Exp.hh"
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64 
65 using namespace std;
66 
67 static const G4double
68  d1= 2.7965e-1*CLHEP::barn, d2=-1.8300e-1*CLHEP::barn,
69  d3= 6.7527 *CLHEP::barn, d4=-1.9798e+1*CLHEP::barn,
70  e1= 1.9756e-5*CLHEP::barn, e2=-1.0205e-2*CLHEP::barn,
71  e3=-7.3913e-2*CLHEP::barn, e4= 2.7079e-2*CLHEP::barn,
72  f1=-3.9178e-7*CLHEP::barn, f2= 6.8241e-5*CLHEP::barn,
73  f3= 6.0480e-5*CLHEP::barn, f4= 3.0274e-4*CLHEP::barn;
74 static const G4double dT0 = keV;
75 static const G4int nlooplim = 1000;
76 
78  : G4VEmModel(nam)
79 {
83  limitFactor = 4;
84  fProbabilities.resize(9,0.0);
85  SetDeexcitationFlag(true);
86  fParticleChange = 0;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93 {}
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96 
98  const G4DataVector& cuts)
99 {
101  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106 
108  G4VEmModel* masterModel)
109 {
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114 
115 G4double
117  G4double gammaEnergy,
118  G4double Z, G4double,
120 {
121  G4double xSection = 0.0 ;
122  if (gammaEnergy <= LowEnergyLimit()) { return xSection; }
123 
124  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
125 
126  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
127  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
128 
129  G4double T0 = 15.0*keV;
130  if (Z < 1.5) { T0 = 40.0*keV; }
131 
132  G4double X = max(gammaEnergy, T0) / electron_mass_c2;
133  xSection = p1Z*G4Log(1.+2.*X)/X
134  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
135 
136  // modification for low energy. (special case for Hydrogen)
137  if (gammaEnergy < T0) {
138  X = (T0+dT0) / electron_mass_c2 ;
139  G4double sigma = p1Z*G4Log(1.+2*X)/X
140  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
141  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
142  G4double c2 = 0.150;
143  if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
144  G4double y = G4Log(gammaEnergy/T0);
145  xSection *= G4Exp(-y*(c1+c2*y));
146  }
147 
148  if(xSection < 0.0) { xSection = 0.0; }
149  // G4cout << "e= " << GammaEnergy << " Z= " << Z
150  // << " cross= " << xSection << G4endl;
151  return xSection;
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
157  std::vector<G4DynamicParticle*>* fvect,
158  const G4MaterialCutsCouple* couple,
159  const G4DynamicParticle* aDynamicGamma,
160  G4double,
161  G4double)
162 {
163  // primary gamma
164  G4double energy = aDynamicGamma->GetKineticEnergy();
165 
166  // do nothing below the threshold
167  if(energy <= LowEnergyLimit()) { return; }
168 
169  G4ThreeVector direction = aDynamicGamma->GetMomentumDirection();
170 
171  // select atom
172  const G4Element* elm = SelectRandomAtom(couple, theGamma, energy);
173 
174  // select shell first
175  G4int nShells = elm->GetNbOfAtomicShells();
176  if(nShells > (G4int)fProbabilities.size()) { fProbabilities.resize(nShells); }
177  G4double totprob = 0.0;
178  G4int i;
179  for(i=0; i<nShells; ++i) {
180  //G4double bindingEnergy = elm->GetAtomicShell(i);
181  totprob += elm->GetNbOfShellElectrons(i);
182  //totprob += elm->GetNbOfShellElectrons(i)/(bindingEnergy*bindingEnergy);
183  fProbabilities[i] = totprob;
184  }
185 
186  // Loop on sampling
187  G4int nloop = 0;
188 
189  G4double bindingEnergy, ePotEnergy, eKinEnergy;
190  G4double gamEnergy0, gamEnergy1;
191 
192  //static const G4double eminus2 = 1.0 - G4Exp(-2.0);
193 
194  do {
195  ++nloop;
196  G4double xprob = totprob*rndmEngineMod->flat();
197 
198  // select shell
199  for(i=0; i<nShells; ++i) { if(xprob <= fProbabilities[i]) { break; } }
200 
201  bindingEnergy = elm->GetAtomicShell(i);
202  lv1.set(0.0,0.0,energy,energy);
203  /*
204  G4cout << "nShells= " << nShells << " i= " << i
205  << " Egamma= " << energy << " Ebind= " << bindingEnergy
206  << G4endl;
207  */
208  // for rest frame of the electron
209  G4double x = -G4Log(rndmEngineMod->flat());
210  eKinEnergy = bindingEnergy*x;
211  ePotEnergy = bindingEnergy*(1.0 + x);
212 
213  // for rest frame of the electron
214  G4double eTotMomentum = sqrt(eKinEnergy*(eKinEnergy + 2*electron_mass_c2));
215  G4double phi = rndmEngineMod->flat()*twopi;
216  G4double costet = 2*rndmEngineMod->flat() - 1;
217  G4double sintet = sqrt((1 - costet)*(1 + costet));
218  lv2.set(eTotMomentum*sintet*cos(phi),eTotMomentum*sintet*sin(phi),
219  eTotMomentum*costet,eKinEnergy + electron_mass_c2);
220  bst = lv2.boostVector();
221  lv1.boost(-bst);
222 
223  gamEnergy0 = lv1.e();
224 
225  // In the rest frame of the electron
226  // The scattered gamma energy is sampled according to Klein-Nishina formula
227  // The random number techniques of Butcher & Messel are used
228  // (Nuc Phys 20(1960),15).
229  G4double E0_m = gamEnergy0/electron_mass_c2;
230 
231  //G4cout << "Nloop= "<< nloop << " Ecm(keV)= " << gamEnergy0/keV << G4endl;
232  //
233  // sample the energy rate of the scattered gamma
234  //
235 
236  G4double epsilon, epsilonsq, onecost, sint2, greject ;
237 
238  G4double eps0 = 1./(1 + 2*E0_m);
239  G4double epsilon0sq = eps0*eps0;
240  G4double alpha1 = - G4Log(eps0);
241  G4double alpha2 = 0.5*(1 - epsilon0sq);
242 
243  do {
244  ++nloop;
245  // false interaction if too many iterations
246  if(nloop > nlooplim) { return; }
247 
248  if ( alpha1/(alpha1+alpha2) > rndmEngineMod->flat() ) {
249  epsilon = G4Exp(-alpha1*rndmEngineMod->flat()); // epsilon0**r
250  epsilonsq = epsilon*epsilon;
251 
252  } else {
253  epsilonsq = epsilon0sq + (1.- epsilon0sq)*rndmEngineMod->flat();
254  epsilon = sqrt(epsilonsq);
255  }
256 
257  onecost = (1.- epsilon)/(epsilon*E0_m);
258  sint2 = onecost*(2.-onecost);
259  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
260 
261  } while (greject < rndmEngineMod->flat());
262  gamEnergy1 = epsilon*gamEnergy0;
263 
264  // before scattering total 4-momentum in e- system
265  lv2.set(0.0,0.0,0.0,electron_mass_c2);
266  lv2 += lv1;
267 
268  //
269  // scattered gamma angles. ( Z - axis along the parent gamma)
270  //
271  if(sint2 < 0.0) { sint2 = 0.0; }
272  costet = 1. - onecost;
273  sintet = sqrt(sint2);
274  phi = twopi * rndmEngineMod->flat();
275 
276  // e- recoil
277  //
278  // in rest frame of the electron
279  G4ThreeVector gamDir = lv1.vect().unit();
280  G4ThreeVector v = G4ThreeVector(sintet*cos(phi),sintet*sin(phi),costet);
281  v.rotateUz(gamDir);
282  lv1.set(gamEnergy1*v.x(),gamEnergy1*v.y(),gamEnergy1*v.z(),gamEnergy1);
283  lv2 -= lv1;
284  //G4cout<<"Egam(keV)= " << lv1.e()/keV
285  // <<" Ee(keV)= " << (lv2.e()-electron_mass_c2)/keV << G4endl;
286  lv2.boost(bst);
287  eKinEnergy = lv2.e() - electron_mass_c2 - ePotEnergy;
288  //G4cout << "Nloop= " << nloop << " eKinEnergy= " << eKinEnergy << G4endl;
289 
290  } while ( eKinEnergy < 0.0 );
291 
292  //
293  // update G4VParticleChange for the scattered gamma
294  //
295 
296  lv1.boost(bst);
297  gamEnergy1 = lv1.e();
298  if(gamEnergy1 > lowestSecondaryEnergy) {
299  G4ThreeVector gamDirection1 = lv1.vect().unit();
300  gamDirection1.rotateUz(direction);
302  } else {
304  gamEnergy1 = 0.0;
305  }
307 
308  //
309  // kinematic of the scattered electron
310  //
311 
312  if(eKinEnergy > lowestSecondaryEnergy) {
313  G4ThreeVector eDirection = lv2.vect().unit();
314  eDirection.rotateUz(direction);
315  G4DynamicParticle* dp =
316  new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
317  fvect->push_back(dp);
318  } else { eKinEnergy = 0.0; }
319 
320  G4double edep = energy - gamEnergy1 - eKinEnergy;
321  G4double esec = 0.0;
322 
323  // sample deexcitation
324  //
325  if(fAtomDeexcitation) {
326  G4int index = couple->GetIndex();
328  G4int Z = G4lrint(elm->GetZ());
330  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
331  G4int nbefore = fvect->size();
332  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
333  G4int nafter = fvect->size();
334  //G4cout << "N1= " << nbefore << " N2= " << nafter << G4endl;
335  for (G4int j=nbefore; j<nafter; ++j) {
336  G4double e = ((*fvect)[j])->GetKineticEnergy();
337  if(esec + e > edep) {
338  // correct energy in order to have energy balance
339  e = edep - esec;
340  ((*fvect)[j])->SetKineticEnergy(e);
341  esec += e;
342  /*
343  G4cout << "### G4KleinNishinaModel Edep(eV)= " << edep/eV
344  << " Esec(eV)= " << esec/eV
345  << " E["<< j << "](eV)= " << e/eV
346  << " N= " << nafter
347  << " Z= " << Z << " shell= " << i
348  << " Ebind(keV)= " << bindingEnergy/keV
349  << " Eshell(keV)= " << shell->BindingEnergy()/keV
350  << G4endl;
351  */
352  // delete the rest of secondaries
353  for (G4int jj=nafter-1; jj>j; --jj) {
354  delete (*fvect)[jj];
355  fvect->pop_back();
356  }
357  break;
358  }
359  esec += e;
360  }
361  edep -= esec;
362  }
363  }
364  if(fabs(energy - gamEnergy1 - eKinEnergy - esec - edep) > eV) {
365  G4cout << "### G4KleinNishinaModel dE(eV)= "
366  << (energy - gamEnergy1 - eKinEnergy - esec - edep)/eV
367  << " shell= " << i
368  << " E(keV)= " << energy/keV
369  << " Ebind(keV)= " << bindingEnergy/keV
370  << " Eg(keV)= " << gamEnergy1/keV
371  << " Ee(keV)= " << eKinEnergy/keV
372  << " Esec(keV)= " << esec/keV
373  << " Edep(keV)= " << edep/keV
374  << G4endl;
375  }
376  // energy balance
377  if(edep > 0.0) {
379  }
380 }
381 
382 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
383 
static c2_factory< G4double > c2
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:622
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
static G4LossTableManager * Instance()
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
static const G4double e3
G4ParticleChangeForGamma * fParticleChange
static const G4double d1
static const G4double e4
G4double GetZ() const
Definition: G4Element.hh:131
G4bool IsMaster() const
Definition: G4VEmModel.hh:699
G4double a
Definition: TRTMaterials.hh:39
G4int GetNbOfShellElectrons(G4int index) const
Definition: G4Element.cc:388
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
static const G4double d4
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static const G4double f4
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4GLOB_DLL std::ostream G4cout
CLHEP::HepRandomEngine * rndmEngineMod
Definition: G4VEmModel.hh:400
std::vector< G4double > fProbabilities
const G4ThreeVector & GetMomentumDirection() const
static const G4double d2
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double flat()
Definition: G4AblaRandom.cc:47
static const G4double c1
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:783
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
static const G4double e1
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4KleinNishinaModel(const G4String &nam="KleinNishina")
static const double eV
Definition: G4SIunits.hh:194
static const G4double e2
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:791
G4ParticleDefinition * theElectron
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double energy(const ThreeVector &p, const G4double m)
static const G4int nlooplim
static const G4double f1
G4VAtomDeexcitation * fAtomDeexcitation
G4ParticleDefinition * theGamma
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
static const double keV
Definition: G4SIunits.hh:195
static const double barn
Definition: G4SIunits.hh:95
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:373
double G4double
Definition: G4Types.hh:76
static const G4double d3
void ProposeTrackStatus(G4TrackStatus status)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:762
static const G4double f2
G4double bindingEnergy(G4int A, G4int Z)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4AtomicShellEnumerator
static const G4double dT0
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:525
static const G4double f3
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134