Geant4  10.01.p02
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 90583 2015-06-04 11:16:41Z 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  lv1(0.,0.,0.,0.),
80  lv2(0.,0.,0.,0.),
81  bst(0.,0.,0.)
82 {
86  limitFactor = 4;
87  fProbabilities.resize(9,0.0);
88  SetDeexcitationFlag(true);
89  fParticleChange = 0;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 
96 {}
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
101  const G4DataVector& cuts)
102 {
104  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109 
111  G4VEmModel* masterModel)
112 {
114 }
115 
116 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117 
118 G4double
120  G4double gammaEnergy,
121  G4double Z, G4double,
123 {
124  G4double xSection = 0.0 ;
125  if (gammaEnergy <= LowEnergyLimit()) { return xSection; }
126 
127  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
128 
129  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
130  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
131 
132  G4double T0 = 15.0*keV;
133  if (Z < 1.5) { T0 = 40.0*keV; }
134 
135  G4double X = max(gammaEnergy, T0) / electron_mass_c2;
136  xSection = p1Z*G4Log(1.+2.*X)/X
137  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
138 
139  // modification for low energy. (special case for Hydrogen)
140  if (gammaEnergy < T0) {
141  X = (T0+dT0) / electron_mass_c2 ;
142  G4double sigma = p1Z*G4Log(1.+2*X)/X
143  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
144  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
145  G4double c2 = 0.150;
146  if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
147  G4double y = G4Log(gammaEnergy/T0);
148  xSection *= G4Exp(-y*(c1+c2*y));
149  }
150 
151  if(xSection < 0.0) { xSection = 0.0; }
152  // G4cout << "e= " << GammaEnergy << " Z= " << Z
153  // << " cross= " << xSection << G4endl;
154  return xSection;
155 }
156 
157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
158 
160  std::vector<G4DynamicParticle*>* fvect,
161  const G4MaterialCutsCouple* couple,
162  const G4DynamicParticle* aDynamicGamma,
163  G4double,
164  G4double)
165 {
166  // primary gamma
167  G4double energy = aDynamicGamma->GetKineticEnergy();
168 
169  // do nothing below the threshold
170  if(energy <= LowEnergyLimit()) { return; }
171 
172  G4ThreeVector direction = aDynamicGamma->GetMomentumDirection();
173 
174  // select atom
175  const G4Element* elm = SelectRandomAtom(couple, theGamma, energy);
176 
177  // select shell first
178  G4int nShells = elm->GetNbOfAtomicShells();
179  if(nShells > (G4int)fProbabilities.size()) { fProbabilities.resize(nShells); }
180  G4double totprob = 0.0;
181  G4int i;
182  for(i=0; i<nShells; ++i) {
183  //G4double bindingEnergy = elm->GetAtomicShell(i);
184  totprob += elm->GetNbOfShellElectrons(i);
185  //totprob += elm->GetNbOfShellElectrons(i)/(bindingEnergy*bindingEnergy);
186  fProbabilities[i] = totprob;
187  }
188 
189  // Loop on sampling
190  G4int nloop = 0;
191 
192  G4double bindingEnergy, ePotEnergy, eKinEnergy;
193  G4double gamEnergy0, gamEnergy1;
194 
195  //static const G4double eminus2 = 1.0 - G4Exp(-2.0);
196 
197  rndmEngineMod = G4Random::getTheEngine();
198  do {
199  ++nloop;
200  G4double xprob = totprob*rndmEngineMod->flat();
201 
202  // select shell
203  for(i=0; i<nShells; ++i) { if(xprob <= fProbabilities[i]) { break; } }
204 
205  bindingEnergy = elm->GetAtomicShell(i);
206  lv1.set(0.0,0.0,energy,energy);
207  /*
208  G4cout << "nShells= " << nShells << " i= " << i
209  << " Egamma= " << energy << " Ebind= " << bindingEnergy
210  << G4endl;
211  */
212  // for rest frame of the electron
213  G4double x = -G4Log(rndmEngineMod->flat());
214  eKinEnergy = bindingEnergy*x;
215  ePotEnergy = bindingEnergy*(1.0 + x);
216 
217  // for rest frame of the electron
218  G4double eTotMomentum = sqrt(eKinEnergy*(eKinEnergy + 2*electron_mass_c2));
219  G4double phi = rndmEngineMod->flat()*twopi;
220  G4double costet = 2*rndmEngineMod->flat() - 1;
221  G4double sintet = sqrt((1 - costet)*(1 + costet));
222  lv2.set(eTotMomentum*sintet*cos(phi),eTotMomentum*sintet*sin(phi),
223  eTotMomentum*costet,eKinEnergy + electron_mass_c2);
224  bst = lv2.boostVector();
225  lv1.boost(-bst);
226 
227  gamEnergy0 = lv1.e();
228 
229  // In the rest frame of the electron
230  // The scattered gamma energy is sampled according to Klein-Nishina formula
231  // The random number techniques of Butcher & Messel are used
232  // (Nuc Phys 20(1960),15).
233  G4double E0_m = gamEnergy0/electron_mass_c2;
234 
235  //G4cout << "Nloop= "<< nloop << " Ecm(keV)= " << gamEnergy0/keV << G4endl;
236  //
237  // sample the energy rate of the scattered gamma
238  //
239 
240  G4double epsilon, epsilonsq, onecost, sint2, greject ;
241 
242  G4double eps0 = 1./(1 + 2*E0_m);
243  G4double epsilon0sq = eps0*eps0;
244  G4double alpha1 = - G4Log(eps0);
245  G4double alpha2 = 0.5*(1 - epsilon0sq);
246 
247  do {
248  ++nloop;
249  // false interaction if too many iterations
250  if(nloop > nlooplim) { return; }
251 
252  if ( alpha1/(alpha1+alpha2) > rndmEngineMod->flat() ) {
253  epsilon = G4Exp(-alpha1*rndmEngineMod->flat()); // epsilon0**r
254  epsilonsq = epsilon*epsilon;
255 
256  } else {
257  epsilonsq = epsilon0sq + (1.- epsilon0sq)*rndmEngineMod->flat();
258  epsilon = sqrt(epsilonsq);
259  }
260 
261  onecost = (1.- epsilon)/(epsilon*E0_m);
262  sint2 = onecost*(2.-onecost);
263  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
264 
265  } while (greject < rndmEngineMod->flat());
266  gamEnergy1 = epsilon*gamEnergy0;
267 
268  // before scattering total 4-momentum in e- system
269  lv2.set(0.0,0.0,0.0,electron_mass_c2);
270  lv2 += lv1;
271 
272  //
273  // scattered gamma angles. ( Z - axis along the parent gamma)
274  //
275  if(sint2 < 0.0) { sint2 = 0.0; }
276  costet = 1. - onecost;
277  sintet = sqrt(sint2);
278  phi = twopi * rndmEngineMod->flat();
279 
280  // e- recoil
281  //
282  // in rest frame of the electron
283  G4ThreeVector gamDir = lv1.vect().unit();
284  G4ThreeVector v = G4ThreeVector(sintet*cos(phi),sintet*sin(phi),costet);
285  v.rotateUz(gamDir);
286  lv1.set(gamEnergy1*v.x(),gamEnergy1*v.y(),gamEnergy1*v.z(),gamEnergy1);
287  lv2 -= lv1;
288  //G4cout<<"Egam(keV)= " << lv1.e()/keV
289  // <<" Ee(keV)= " << (lv2.e()-electron_mass_c2)/keV << G4endl;
290  lv2.boost(bst);
291  eKinEnergy = lv2.e() - electron_mass_c2 - ePotEnergy;
292  //G4cout << "Nloop= " << nloop << " eKinEnergy= " << eKinEnergy << G4endl;
293 
294  } while ( eKinEnergy < 0.0 );
295 
296  //
297  // update G4VParticleChange for the scattered gamma
298  //
299 
300  lv1.boost(bst);
301  gamEnergy1 = lv1.e();
302  if(gamEnergy1 > lowestSecondaryEnergy) {
303  G4ThreeVector gamDirection1 = lv1.vect().unit();
304  gamDirection1.rotateUz(direction);
306  } else {
308  gamEnergy1 = 0.0;
309  }
311 
312  //
313  // kinematic of the scattered electron
314  //
315 
316  if(eKinEnergy > lowestSecondaryEnergy) {
317  G4ThreeVector eDirection = lv2.vect().unit();
318  eDirection.rotateUz(direction);
319  G4DynamicParticle* dp =
320  new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
321  fvect->push_back(dp);
322  } else { eKinEnergy = 0.0; }
323 
324  G4double edep = energy - gamEnergy1 - eKinEnergy;
325  G4double esec = 0.0;
326 
327  // sample deexcitation
328  //
329  if(fAtomDeexcitation) {
330  G4int index = couple->GetIndex();
332  G4int Z = G4lrint(elm->GetZ());
334  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
335  G4int nbefore = fvect->size();
336  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
337  G4int nafter = fvect->size();
338  //G4cout << "N1= " << nbefore << " N2= " << nafter << G4endl;
339  for (G4int j=nbefore; j<nafter; ++j) {
340  G4double e = ((*fvect)[j])->GetKineticEnergy();
341  if(esec + e > edep) {
342  // correct energy in order to have energy balance
343  e = edep - esec;
344  ((*fvect)[j])->SetKineticEnergy(e);
345  esec += e;
346  /*
347  G4cout << "### G4KleinNishinaModel Edep(eV)= " << edep/eV
348  << " Esec(eV)= " << esec/eV
349  << " E["<< j << "](eV)= " << e/eV
350  << " N= " << nafter
351  << " Z= " << Z << " shell= " << i
352  << " Ebind(keV)= " << bindingEnergy/keV
353  << " Eshell(keV)= " << shell->BindingEnergy()/keV
354  << G4endl;
355  */
356  // delete the rest of secondaries (should not happens)
357  for (G4int jj=nafter-1; jj>j; --jj) {
358  delete (*fvect)[jj];
359  fvect->pop_back();
360  }
361  break;
362  }
363  esec += e;
364  }
365  edep -= esec;
366  }
367  }
368  if(fabs(energy - gamEnergy1 - eKinEnergy - esec - edep) > eV) {
369  G4cout << "### G4KleinNishinaModel dE(eV)= "
370  << (energy - gamEnergy1 - eKinEnergy - esec - edep)/eV
371  << " shell= " << i
372  << " E(keV)= " << energy/keV
373  << " Ebind(keV)= " << bindingEnergy/keV
374  << " Eg(keV)= " << gamEnergy1/keV
375  << " Ee(keV)= " << eKinEnergy/keV
376  << " Esec(keV)= " << esec/keV
377  << " Edep(keV)= " << edep/keV
378  << G4endl;
379  }
380  // energy balance
381  if(edep > 0.0) {
383  }
384 }
385 
386 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
387 
static c2_factory< G4double > c2
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:623
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:700
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:784
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:792
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:763
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:526
static const G4double f3
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134