Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 96934 2016-05-18 09:10: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 
68  : G4VEmModel(nam),
69  lv1(0.,0.,0.,0.),
70  lv2(0.,0.,0.,0.),
71  bst(0.,0.,0.)
72 {
76  limitFactor = 4;
77  fProbabilities.resize(9,0.0);
78  SetDeexcitationFlag(true);
79  fParticleChange = nullptr;
80  fAtomDeexcitation = nullptr;
81 }
82 
83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84 
86 {}
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89 
91  const G4DataVector& cuts)
92 {
93  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
94  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
95  if(nullptr == fParticleChange) {
97  }
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,
115 {
116  G4double xSection = 0.0 ;
117  if (gammaEnergy <= LowEnergyLimit()) { return xSection; }
118 
119  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
120 
121 static const G4double
122  d1= 2.7965e-1*CLHEP::barn, d2=-1.8300e-1*CLHEP::barn,
123  d3= 6.7527 *CLHEP::barn, d4=-1.9798e+1*CLHEP::barn,
124  e1= 1.9756e-5*CLHEP::barn, e2=-1.0205e-2*CLHEP::barn,
125  e3=-7.3913e-2*CLHEP::barn, e4= 2.7079e-2*CLHEP::barn,
126  f1=-3.9178e-7*CLHEP::barn, f2= 6.8241e-5*CLHEP::barn,
127  f3= 6.0480e-5*CLHEP::barn, f4= 3.0274e-4*CLHEP::barn;
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  static const G4double dT0 = keV;
141  if (gammaEnergy < T0) {
142  X = (T0+dT0) / electron_mass_c2 ;
143  G4double sigma = p1Z*G4Log(1.+2*X)/X
144  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
145  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
146  G4double c2 = 0.150;
147  if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
148  G4double y = G4Log(gammaEnergy/T0);
149  xSection *= G4Exp(-y*(c1+c2*y));
150  }
151 
152  if(xSection < 0.0) { xSection = 0.0; }
153  // G4cout << "e= " << GammaEnergy << " Z= " << Z
154  // << " cross= " << xSection << G4endl;
155  return xSection;
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
159 
161  std::vector<G4DynamicParticle*>* fvect,
162  const G4MaterialCutsCouple* couple,
163  const G4DynamicParticle* aDynamicGamma,
164  G4double,
165  G4double)
166 {
167  // primary gamma
168  G4double energy = aDynamicGamma->GetKineticEnergy();
169 
170  // do nothing below the threshold
171  if(energy <= LowEnergyLimit()) { return; }
172 
173  G4ThreeVector direction = aDynamicGamma->GetMomentumDirection();
174 
175  // select atom
176  const G4Element* elm = SelectRandomAtom(couple, theGamma, energy);
177 
178  // select shell first
179  G4int nShells = elm->GetNbOfAtomicShells();
180  if(nShells > (G4int)fProbabilities.size()) { fProbabilities.resize(nShells); }
181  G4double totprob = 0.0;
182  G4int i;
183  for(i=0; i<nShells; ++i) {
184  //G4double bindingEnergy = elm->GetAtomicShell(i);
185  totprob += elm->GetNbOfShellElectrons(i);
186  //totprob += elm->GetNbOfShellElectrons(i)/(bindingEnergy*bindingEnergy);
187  fProbabilities[i] = totprob;
188  }
189 
190  // Loop on sampling
191  static const G4int nlooplim = 1000;
192  G4int nloop = 0;
193 
194  G4double bindingEnergy, ePotEnergy, eKinEnergy;
195  G4double gamEnergy0, gamEnergy1;
196 
197  CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
198  G4double rndm[4];
199 
200  do {
201  ++nloop;
202 
203  // 4 random numbers to select e-
204  rndmEngineMod->flatArray(4, rndm);
205  G4double xprob = totprob*rndm[0];
206 
207  // select shell
208  for(i=0; i<nShells; ++i) { if(xprob <= fProbabilities[i]) { break; } }
209 
210  bindingEnergy = elm->GetAtomicShell(i);
211  lv1.set(0.0,0.0,energy,energy);
212  /*
213  G4cout << "nShells= " << nShells << " i= " << i
214  << " Egamma= " << energy << " Ebind= " << bindingEnergy
215  << G4endl;
216  */
217  // for rest frame of the electron
218  G4double x = -G4Log(rndm[1]);
219  eKinEnergy = bindingEnergy*x;
220  ePotEnergy = bindingEnergy*(1.0 + x);
221 
222  // for rest frame of the electron
223  G4double eTotMomentum = sqrt(eKinEnergy*(eKinEnergy + 2*electron_mass_c2));
224  G4double phi = rndm[2]*twopi;
225  G4double costet = 2*rndm[3] - 1;
226  G4double sintet = sqrt((1 - costet)*(1 + costet));
227  lv2.set(eTotMomentum*sintet*cos(phi),eTotMomentum*sintet*sin(phi),
228  eTotMomentum*costet,eKinEnergy + electron_mass_c2);
229  bst = lv2.boostVector();
230  lv1.boost(-bst);
231 
232  gamEnergy0 = lv1.e();
233 
234  // In the rest frame of the electron
235  // The scattered gamma energy is sampled according to Klein-Nishina formula
236  // The random number techniques of Butcher & Messel are used
237  // (Nuc Phys 20(1960),15).
238  G4double E0_m = gamEnergy0/electron_mass_c2;
239 
240  //G4cout << "Nloop= "<< nloop << " Ecm(keV)= " << gamEnergy0/keV << G4endl;
241  //
242  // sample the energy rate of the scattered gamma
243  //
244 
245  G4double epsilon, epsilonsq, onecost, sint2, greject ;
246 
247  G4double eps0 = 1./(1 + 2*E0_m);
248  G4double epsilon0sq = eps0*eps0;
249  G4double alpha1 = - G4Log(eps0);
250  G4double alpha2 = alpha1 + 0.5*(1 - epsilon0sq);
251 
252  do {
253  ++nloop;
254  // false interaction if too many iterations
255  if(nloop > nlooplim) { return; }
256 
257  // 3 random numbers to sample scattering
258  rndmEngineMod->flatArray(3, rndm);
259 
260  if ( alpha1 > alpha2*rndm[0] ) {
261  epsilon = G4Exp(-alpha1*rndm[1]); // epsilon0**r
262  epsilonsq = epsilon*epsilon;
263 
264  } else {
265  epsilonsq = epsilon0sq + (1.- epsilon0sq)*rndm[1];
266  epsilon = sqrt(epsilonsq);
267  }
268 
269  onecost = (1.- epsilon)/(epsilon*E0_m);
270  sint2 = onecost*(2.-onecost);
271  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
272 
273  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
274  } while (greject < rndm[2]);
275  gamEnergy1 = epsilon*gamEnergy0;
276 
277  // before scattering total 4-momentum in e- system
278  lv2.set(0.0,0.0,0.0,electron_mass_c2);
279  lv2 += lv1;
280 
281  //
282  // scattered gamma angles. ( Z - axis along the parent gamma)
283  //
284  if(sint2 < 0.0) { sint2 = 0.0; }
285  costet = 1. - onecost;
286  sintet = sqrt(sint2);
287  phi = twopi * rndmEngineMod->flat();
288 
289  // e- recoil
290  //
291  // in rest frame of the electron
292  G4ThreeVector gamDir = lv1.vect().unit();
293  G4ThreeVector v = G4ThreeVector(sintet*cos(phi),sintet*sin(phi),costet);
294  v.rotateUz(gamDir);
295  lv1.set(gamEnergy1*v.x(),gamEnergy1*v.y(),gamEnergy1*v.z(),gamEnergy1);
296  lv2 -= lv1;
297  //G4cout<<"Egam(keV)= " << lv1.e()/keV
298  // <<" Ee(keV)= " << (lv2.e()-electron_mass_c2)/keV << G4endl;
299  lv2.boost(bst);
300  eKinEnergy = lv2.e() - electron_mass_c2 - ePotEnergy;
301  //G4cout << "Nloop= " << nloop << " eKinEnergy= " << eKinEnergy << G4endl;
302 
303  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
304  } while ( eKinEnergy < 0.0 );
305 
306  //
307  // update G4VParticleChange for the scattered gamma
308  //
309 
310  lv1.boost(bst);
311  gamEnergy1 = lv1.e();
312  if(gamEnergy1 > lowestSecondaryEnergy) {
313  G4ThreeVector gamDirection1 = lv1.vect().unit();
314  gamDirection1.rotateUz(direction);
316  } else {
318  gamEnergy1 = 0.0;
319  }
321 
322  //
323  // kinematic of the scattered electron
324  //
325 
326  if(eKinEnergy > lowestSecondaryEnergy) {
327  G4ThreeVector eDirection = lv2.vect().unit();
328  eDirection.rotateUz(direction);
329  G4DynamicParticle* dp =
330  new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
331  fvect->push_back(dp);
332  } else { eKinEnergy = 0.0; }
333 
334  G4double edep = energy - gamEnergy1 - eKinEnergy;
335  G4double esec = 0.0;
336 
337  // sample deexcitation
338  //
339  if(fAtomDeexcitation) {
340  G4int index = couple->GetIndex();
341  if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
342  G4int Z = G4lrint(elm->GetZ());
344  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
345  G4int nbefore = fvect->size();
346  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
347  G4int nafter = fvect->size();
348  //G4cout << "N1= " << nbefore << " N2= " << nafter << G4endl;
349  for (G4int j=nbefore; j<nafter; ++j) {
350  G4double e = ((*fvect)[j])->GetKineticEnergy();
351  if(esec + e > edep) {
352  // correct energy in order to have energy balance
353  e = edep - esec;
354  ((*fvect)[j])->SetKineticEnergy(e);
355  esec += e;
356  /*
357  G4cout << "### G4KleinNishinaModel Edep(eV)= " << edep/eV
358  << " Esec(eV)= " << esec/eV
359  << " E["<< j << "](eV)= " << e/eV
360  << " N= " << nafter
361  << " Z= " << Z << " shell= " << i
362  << " Ebind(keV)= " << bindingEnergy/keV
363  << " Eshell(keV)= " << shell->BindingEnergy()/keV
364  << G4endl;
365  */
366  // delete the rest of secondaries (should not happens)
367  for (G4int jj=nafter-1; jj>j; --jj) {
368  delete (*fvect)[jj];
369  fvect->pop_back();
370  }
371  break;
372  }
373  esec += e;
374  }
375  edep -= esec;
376  }
377  }
378  if(std::abs(energy - gamEnergy1 - eKinEnergy - esec - edep) > eV) {
379  G4cout << "### G4KleinNishinaModel dE(eV)= "
380  << (energy - gamEnergy1 - eKinEnergy - esec - edep)/eV
381  << " shell= " << i
382  << " E(keV)= " << energy/keV
383  << " Ebind(keV)= " << bindingEnergy/keV
384  << " Eg(keV)= " << gamEnergy1/keV
385  << " Ee(keV)= " << eKinEnergy/keV
386  << " Esec(keV)= " << esec/keV
387  << " Edep(keV)= " << edep/keV
388  << G4endl;
389  }
390  // energy balance
391  if(edep > 0.0) {
393  }
394 }
395 
396 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
397 
Hep3Vector boostVector() const
static c2_factory< G4double > c2
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
static G4LossTableManager * Instance()
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:147
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
double x() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
G4ParticleChangeForGamma * fParticleChange
static const G4int nlooplim
const char * p
Definition: xmltok.h:285
static const G4double d2
G4double GetZ() const
Definition: G4Element.hh:131
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
tuple x
Definition: test.py:50
virtual double flat()=0
G4int GetNbOfShellElectrons(G4int index) const
Definition: G4Element.cc:382
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
Definition: G4SIunits.hh:76
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
tuple b
Definition: test.py:12
Hep3Vector vect() const
G4GLOB_DLL std::ostream G4cout
static constexpr double barn
Definition: SystemOfUnits.h:85
const G4ThreeVector & GetMomentumDirection() const
HepLorentzVector & boost(double, double, double)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
static constexpr double eV
Definition: G4SIunits.hh:215
float electron_mass_c2
Definition: hepunit.py:274
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:809
static const G4double d1
tuple v
Definition: test.py:18
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")
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:817
void set(double x, double y, double z, double t)
G4ParticleDefinition * theElectron
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel) override
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double energy(const ThreeVector &p, const G4double m)
Hep3Vector unit() const
double y() const
G4ParticleDefinition * theGamma
static const G4double T0[78]
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:367
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
tuple c
Definition: test.py:13
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:788
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
G4double bindingEnergy(G4int A, G4int Z)
virtual void flatArray(const int size, double *vect)=0
static constexpr double keV
Definition: G4SIunits.hh:216
double epsilon(double density, double temperature)
G4AtomicShellEnumerator
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:544
tuple c1
Definition: plottest35.py:14
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132