Geant4  10.00.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 74309 2013-10-03 06:42:30Z 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*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
69  e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
70  f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
71 
73  : G4VEmModel(nam)
74 {
77  lowestGammaEnergy = 1.0*eV;
78  limitFactor = 4;
79  fProbabilities.resize(9,0.0);
80  SetDeexcitationFlag(true);
81  fParticleChange = 0;
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
88 {}
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93  const G4DataVector& cuts)
94 {
96  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
98 }
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101 
103  G4VEmModel* masterModel)
104 {
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109 
110 G4double
112  G4double GammaEnergy,
113  G4double Z, G4double,
115 {
116  G4double xSection = 0.0 ;
117  if ( Z < 0.9999 || GammaEnergy < 0.1*keV) { return xSection; }
118 
119  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
120 
121  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
122  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
123 
124  G4double T0 = 15.0*keV;
125  if (Z < 1.5) { T0 = 40.0*keV; }
126 
127  G4double X = max(GammaEnergy, T0) / electron_mass_c2;
128  xSection = p1Z*G4Log(1.+2.*X)/X
129  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
130 
131  // modification for low energy. (special case for Hydrogen)
132  if (GammaEnergy < T0) {
133  G4double dT0 = keV;
134  X = (T0+dT0) / electron_mass_c2 ;
135  G4double sigma = p1Z*G4Log(1.+2*X)/X
136  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
137  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
138  G4double c2 = 0.150;
139  if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
140  G4double y = G4Log(GammaEnergy/T0);
141  xSection *= G4Exp(-y*(c1+c2*y));
142  }
143 
144  if(xSection < 0.0) { xSection = 0.0; }
145  // G4cout << "e= " << GammaEnergy << " Z= " << Z
146  // << " cross= " << xSection << G4endl;
147  return xSection;
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151 
153  std::vector<G4DynamicParticle*>* fvect,
154  const G4MaterialCutsCouple* couple,
155  const G4DynamicParticle* aDynamicGamma,
156  G4double,
157  G4double)
158 {
159  // primary gamma
160  G4double energy = aDynamicGamma->GetKineticEnergy();
161  G4ThreeVector direction = aDynamicGamma->GetMomentumDirection();
162 
163  // select atom
164  const G4Element* elm = SelectRandomAtom(couple, theGamma, energy);
165 
166  // select shell first
167  G4int nShells = elm->GetNbOfAtomicShells();
168  if(nShells > (G4int)fProbabilities.size()) { fProbabilities.resize(nShells); }
169  G4double totprob = 0.0;
170  G4int i;
171  for(i=0; i<nShells; ++i) {
172  //G4double bindingEnergy = elm->GetAtomicShell(i);
173  totprob += elm->GetNbOfShellElectrons(i);
174  //totprob += elm->GetNbOfShellElectrons(i)/(bindingEnergy*bindingEnergy);
175  fProbabilities[i] = totprob;
176  }
177  //if(totprob == 0.0) { return; }
178 
179  // Loop on sampling
180  // const G4int nlooplim = 100;
181  //G4int nloop = 0;
182 
183  G4double bindingEnergy, ePotEnergy, eKinEnergy;
184  G4double gamEnergy0, gamEnergy1;
185 
186  //static const G4double eminus2 = 1.0 - G4Exp(-2.0);
187 
188  do {
189  //++nloop;
190  G4double xprob = totprob*G4UniformRand();
191 
192  // select shell
193  for(i=0; i<nShells; ++i) { if(xprob <= fProbabilities[i]) { break; } }
194 
195  bindingEnergy = elm->GetAtomicShell(i);
196  // ePotEnergy = bindingEnergy;
197  // gamEnergy0 = energy;
198  lv1.set(0.0,0.0,energy,energy);
199 
200  //G4cout << "nShells= " << nShells << " i= " << i
201  // << " Egamma= " << energy << " Ebind= " << bindingEnergy
202  // << " Elim= " << limitEnergy
203  // << G4endl;
204 
205  // for rest frame of the electron
206  G4double x = -G4Log(G4UniformRand());
207  eKinEnergy = bindingEnergy*x;
208  ePotEnergy = bindingEnergy*(1.0 + x);
209 
210  // for rest frame of the electron
211  G4double eTotMomentum = sqrt(eKinEnergy*(eKinEnergy + 2*electron_mass_c2));
212  G4double phi = G4UniformRand()*twopi;
213  G4double costet = 2*G4UniformRand() - 1;
214  G4double sintet = sqrt((1 - costet)*(1 + costet));
215  lv2.set(eTotMomentum*sintet*cos(phi),eTotMomentum*sintet*sin(phi),
216  eTotMomentum*costet,eKinEnergy + electron_mass_c2);
217  bst = lv2.boostVector();
218  lv1.boost(-bst);
219 
220  gamEnergy0 = lv1.e();
221 
222  // In the rest frame of the electron
223  // The scattered gamma energy is sampled according to Klein-Nishina formula
224  // The random number techniques of Butcher & Messel are used
225  // (Nuc Phys 20(1960),15).
226  G4double E0_m = gamEnergy0/electron_mass_c2;
227 
228  //
229  // sample the energy rate of the scattered gamma
230  //
231 
232  G4double epsilon, epsilonsq, onecost, sint2, greject ;
233 
234  G4double eps0 = 1./(1 + 2*E0_m);
235  G4double epsilon0sq = eps0*eps0;
236  G4double alpha1 = - G4Log(eps0);
237  G4double alpha2 = 0.5*(1 - epsilon0sq);
238 
239  do {
240  if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
241  epsilon = G4Exp(-alpha1*G4UniformRand()); // epsilon0**r
242  epsilonsq = epsilon*epsilon;
243 
244  } else {
245  epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
246  epsilon = sqrt(epsilonsq);
247  }
248 
249  onecost = (1.- epsilon)/(epsilon*E0_m);
250  sint2 = onecost*(2.-onecost);
251  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
252 
253  } while (greject < G4UniformRand());
254  gamEnergy1 = epsilon*gamEnergy0;
255 
256  // before scattering total 4-momentum in e- system
257  lv2.set(0.0,0.0,0.0,electron_mass_c2);
258  lv2 += lv1;
259 
260  //
261  // scattered gamma angles. ( Z - axis along the parent gamma)
262  //
263  if(sint2 < 0.0) { sint2 = 0.0; }
264  costet = 1. - onecost;
265  sintet = sqrt(sint2);
266  phi = twopi * G4UniformRand();
267 
268  // e- recoil
269  //
270  // in rest frame of the electron
271  G4ThreeVector gamDir = lv1.vect().unit();
272  G4ThreeVector v = G4ThreeVector(sintet*cos(phi),sintet*sin(phi),costet);
273  v.rotateUz(gamDir);
274  lv1.set(gamEnergy1*v.x(),gamEnergy1*v.y(),gamEnergy1*v.z(),gamEnergy1);
275  lv2 -= lv1;
276  //G4cout<<"Egam= "<<lv1.e()<<" Ee= "<< lv2.e()-electron_mass_c2 << G4endl;
277  lv2.boost(bst);
278  eKinEnergy = lv2.e() - electron_mass_c2 - ePotEnergy;
279  //G4cout << "eKinEnergy= " << eKinEnergy << G4endl;
280 
281  } while ( eKinEnergy < 0.0 );
282 
283  //
284  // update G4VParticleChange for the scattered gamma
285  //
286 
287  lv1.boost(bst);
288  gamEnergy1 = lv1.e();
289  if(gamEnergy1 > lowestGammaEnergy) {
290  G4ThreeVector gamDirection1 = lv1.vect().unit();
291  gamDirection1.rotateUz(direction);
293  } else {
295  gamEnergy1 = 0.0;
296  }
298 
299  //
300  // kinematic of the scattered electron
301  //
302 
303  if(eKinEnergy > lowestGammaEnergy) {
304  G4ThreeVector eDirection = lv2.vect().unit();
305  eDirection.rotateUz(direction);
306  G4DynamicParticle* dp =
307  new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
308  fvect->push_back(dp);
309  } else { eKinEnergy = 0.0; }
310 
311  G4double edep = energy - gamEnergy1 - eKinEnergy;
312  G4double esec = 0.0;
313 
314  // sample deexcitation
315  //
316  if(fAtomDeexcitation) {
317  G4int index = couple->GetIndex();
319  G4int Z = G4lrint(elm->GetZ());
321  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
322  size_t nbefore = fvect->size();
323  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
324  size_t nafter = fvect->size();
325  if(nafter > nbefore) {
326  for (size_t j=nbefore; j<nafter; ++j) {
327  G4double e = ((*fvect)[j])->GetKineticEnergy();
328  if(esec + e > edep) {
329  /*
330  G4cout << "### G4KleinNishinaModel Edep(eV)= " << edep/eV
331  << " Esec(eV)= " << esec/eV
332  << " E["<< j << "](eV)= " << e/eV
333  << " N= " << nafter
334  << " Z= " << Z << " shell= " << i
335  << " Ebind(keV)= " << bindingEnergy/keV
336  << " Eshell(keV)= " << shell->BindingEnergy()/keV
337  << G4endl;
338  */
339  for (size_t jj=j; jj<nafter; ++jj) { delete (*fvect)[jj]; }
340  for (size_t jj=j; jj<nafter; ++jj) { fvect->pop_back(); }
341  break;
342  }
343  esec += e;
344  }
345  }
346  edep -= esec;
347  }
348  }
349  if(fabs(energy - gamEnergy1 - eKinEnergy - esec - edep) > eV) {
350  G4cout << "### G4KleinNishinaModel dE(eV)= "
351  << (energy - gamEnergy1 - eKinEnergy - esec - edep)/eV
352  << " shell= " << i
353  << " E(keV)= " << energy/keV
354  << " Ebind(keV)= " << bindingEnergy/keV
355  << " Eg(keV)= " << gamEnergy1/keV
356  << " Ee(keV)= " << eKinEnergy/keV
357  << " Esec(keV)= " << esec/keV
358  << " Edep(keV)= " << edep/keV
359  << G4endl;
360  }
361  // energy balance
362  if(edep > 0.0) {
364  }
365 }
366 
367 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
368 
static c2_factory< G4double > c2
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:135
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:676
G4double a
Definition: TRTMaterials.hh:39
G4int GetNbOfShellElectrons(G4int index) const
Definition: G4Element.cc:380
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
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
std::vector< G4double > fProbabilities
const G4ThreeVector & GetMomentumDirection() const
static const G4double d2
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static const G4double c1
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
static const G4double e1
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
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:768
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 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:365
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:739
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
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
static const G4double f3
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121