Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PenelopeComptonModel.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$
27 //
28 // Author: Luciano Pandola
29 //
30 // History:
31 // --------
32 // 15 Feb 2010 L Pandola Implementation
33 // 18 Mar 2010 L Pandola Removed GetAtomsPerMolecule(), now demanded
34 // to G4PenelopeOscillatorManager
35 // 01 Feb 2011 L Pandola Suppress fake energy-violation warning when Auger is active.
36 // Make sure that fluorescence/Auger is generated only if
37 // above threshold
38 // 24 May 2011 L Pandola Renamed (make v2008 as default Penelope)
39 // 10 Jun 2011 L Pandola Migrate atomic deexcitation interface
40 //
42 #include "G4PhysicalConstants.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4ParticleDefinition.hh"
45 #include "G4MaterialCutsCouple.hh"
46 #include "G4DynamicParticle.hh"
47 #include "G4VEMDataSet.hh"
48 #include "G4PhysicsTable.hh"
49 #include "G4PhysicsLogVector.hh"
51 #include "G4AtomicShell.hh"
52 #include "G4Gamma.hh"
53 #include "G4Electron.hh"
55 #include "G4PenelopeOscillator.hh"
56 #include "G4LossTableManager.hh"
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59 
60 
62  const G4String& nam)
63  :G4VEmModel(nam),fParticleChange(0),isInitialised(false),fAtomDeexcitation(0),
64  oscManager(0)
65 {
66  fIntrinsicLowEnergyLimit = 100.0*eV;
67  fIntrinsicHighEnergyLimit = 100.0*GeV;
68  // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
69  SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
70  //
72 
73  verboseLevel= 0;
74  // Verbosity scale:
75  // 0 = nothing
76  // 1 = warning for energy non-conservation
77  // 2 = details of energy budget
78  // 3 = calculation of cross sections, file openings, sampling of atoms
79  // 4 = entering in methods
80 
81  //Mark this model as "applicable" for atomic deexcitation
82  SetDeexcitationFlag(true);
83 
84  fTransitionManager = G4AtomicTransitionManager::Instance();
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
90 {;}
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93 
95  const G4DataVector&)
96 {
97  if (verboseLevel > 3)
98  G4cout << "Calling G4PenelopeComptonModel::Initialise()" << G4endl;
99 
100  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
101  //Issue warning if the AtomicDeexcitation has not been declared
102  if (!fAtomDeexcitation)
103  {
104  G4cout << G4endl;
105  G4cout << "WARNING from G4PenelopeComptonModel " << G4endl;
106  G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
107  G4cout << "any fluorescence/Auger emission." << G4endl;
108  G4cout << "Please make sure this is intended" << G4endl;
109  }
110 
111 
112  if (verboseLevel > 0)
113  {
114  G4cout << "Penelope Compton model v2008 is initialized " << G4endl
115  << "Energy range: "
116  << LowEnergyLimit() / keV << " keV - "
117  << HighEnergyLimit() / GeV << " GeV";
118  }
119 
120  if(isInitialised) return;
122  isInitialised = true;
123 
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127 
129  const G4ParticleDefinition* p,
131  G4double,
132  G4double)
133 {
134  // Penelope model v2008 to calculate the Compton scattering cross section:
135  // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167
136  //
137  // The cross section for Compton scattering is calculated according to the Klein-Nishina
138  // formula for energy > 5 MeV.
139  // For E < 5 MeV it is used a parametrization for the differential cross-section dSigma/dOmega,
140  // which is integrated numerically in cos(theta), G4PenelopeComptonModel::DifferentialCrossSection().
141  // The parametrization includes the J(p)
142  // distribution profiles for the atomic shells, that are tabulated from Hartree-Fock calculations
143  // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201
144  //
145  if (verboseLevel > 3)
146  G4cout << "Calling CrossSectionPerVolume() of G4PenelopeComptonModel" << G4endl;
147  SetupForMaterial(p, material, energy);
148 
149  //Retrieve the oscillator table for this material
150  G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableCompton(material);
151 
152  G4double cs = 0;
153 
154  if (energy < 5*MeV) //explicit calculation for E < 5 MeV
155  {
156  size_t numberOfOscillators = theTable->size();
157  for (size_t i=0;i<numberOfOscillators;i++)
158  {
159  G4PenelopeOscillator* theOsc = (*theTable)[i];
160  //sum contributions coming from each oscillator
161  cs += OscillatorTotalCrossSection(energy,theOsc);
162  }
163  }
164  else //use Klein-Nishina for E>5 MeV
165  cs = KleinNishinaCrossSection(energy,material);
166 
167  //cross sections are in units of pi*classic_electr_radius^2
169 
170  //Now, cs is the cross section *per molecule*, let's calculate the
171  //cross section per volume
172 
173  G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
174  G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
175 
176  if (verboseLevel > 3)
177  G4cout << "Material " << material->GetName() << " has " << atPerMol <<
178  "atoms per molecule" << G4endl;
179 
180  G4double moleculeDensity = 0.;
181 
182  if (atPerMol)
183  moleculeDensity = atomDensity/atPerMol;
184 
185  G4double csvolume = cs*moleculeDensity;
186 
187  if (verboseLevel > 2)
188  G4cout << "Compton mean free path at " << energy/keV << " keV for material " <<
189  material->GetName() << " = " << (1./csvolume)/mm << " mm" << G4endl;
190  return csvolume;
191 }
192 
193 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194 
195 //This is a dummy method. Never inkoved by the tracking, it just issues
196 //a warning if one tries to get Cross Sections per Atom via the
197 //G4EmCalculator.
199  G4double,
200  G4double,
201  G4double,
202  G4double,
203  G4double)
204 {
205  G4cout << "*** G4PenelopeComptonModel -- WARNING ***" << G4endl;
206  G4cout << "Penelope Compton model v2008 does not calculate cross section _per atom_ " << G4endl;
207  G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
208  G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
209  return 0;
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
213 
214 void G4PenelopeComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
215  const G4MaterialCutsCouple* couple,
216  const G4DynamicParticle* aDynamicGamma,
217  G4double,
218  G4double)
219 {
220 
221  // Penelope model v2008 to sample the Compton scattering final state.
222  // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167
223  // The model determines also the original shell from which the electron is expelled,
224  // in order to produce fluorescence de-excitation (from G4DeexcitationManager)
225  //
226  // The final state for Compton scattering is calculated according to the Klein-Nishina
227  // formula for energy > 5 MeV. In this case, the Doppler broadening is negligible and
228  // one can assume that the target electron is at rest.
229  // For E < 5 MeV it is used the parametrization for the differential cross-section dSigma/dOmega,
230  // to sample the scattering angle and the energy of the emerging electron, which is
231  // G4PenelopeComptonModel::DifferentialCrossSection(). The rejection method is
232  // used to sample cos(theta). The efficiency increases monotonically with photon energy and is
233  // nearly independent on the Z; typical values are 35%, 80% and 95% for 1 keV, 1 MeV and 10 MeV,
234  // respectively.
235  // The parametrization includes the J(p) distribution profiles for the atomic shells, that are
236  // tabulated
237  // from Hartree-Fock calculations from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201.
238  // Doppler broadening is included.
239  //
240 
241  if (verboseLevel > 3)
242  G4cout << "Calling SampleSecondaries() of G4PenelopeComptonModel" << G4endl;
243 
244  G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
245 
246  if (photonEnergy0 <= fIntrinsicLowEnergyLimit)
247  {
251  return ;
252  }
253 
254  G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
255  const G4Material* material = couple->GetMaterial();
256 
257  G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableCompton(material);
258 
259  const G4int nmax = 64;
260  G4double rn[nmax]={0.0};
261  G4double pac[nmax]={0.0};
262 
263  G4double S=0.0;
264  G4double epsilon = 0.0;
265  G4double cosTheta = 1.0;
266  G4double hartreeFunc = 0.0;
267  G4double oscStren = 0.0;
268  size_t numberOfOscillators = theTable->size();
269  size_t targetOscillator = 0;
270  G4double ionEnergy = 0.0*eV;
271 
272  G4double ek = photonEnergy0/electron_mass_c2;
273  G4double ek2 = 2.*ek+1.0;
274  G4double eks = ek*ek;
275  G4double ek1 = eks-ek2-1.0;
276 
277  G4double taumin = 1.0/ek2;
278  G4double a1 = std::log(ek2);
279  G4double a2 = a1+2.0*ek*(1.0+ek)/(ek2*ek2);
280 
281  G4double TST = 0;
282  G4double tau = 0.;
283 
284  //If the incoming photon is above 5 MeV, the quicker approach based on the
285  //pure Klein-Nishina formula is used
286  if (photonEnergy0 > 5*MeV)
287  {
288  do{
289  do{
290  if ((a2*G4UniformRand()) < a1)
291  tau = std::pow(taumin,G4UniformRand());
292  else
293  tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0));
294  //rejection function
295  TST = (1.0+tau*(ek1+tau*(ek2+tau*eks)))/(eks*tau*(1.0+tau*tau));
296  }while (G4UniformRand()> TST);
297  epsilon=tau;
298  cosTheta = 1.0 - (1.0-tau)/(ek*tau);
299 
300  //Target shell electrons
301  TST = oscManager->GetTotalZ(material)*G4UniformRand();
302  targetOscillator = numberOfOscillators-1; //last level
303  S=0.0;
304  G4bool levelFound = false;
305  for (size_t j=0;j<numberOfOscillators && !levelFound; j++)
306  {
307  S += (*theTable)[j]->GetOscillatorStrength();
308  if (S > TST)
309  {
310  targetOscillator = j;
311  levelFound = true;
312  }
313  }
314  //check whether the level is valid
315  ionEnergy = (*theTable)[targetOscillator]->GetIonisationEnergy();
316  }while((epsilon*photonEnergy0-photonEnergy0+ionEnergy) >0);
317  }
318  else //photonEnergy0 < 5 MeV
319  {
320  //Incoherent scattering function for theta=PI
321  G4double s0=0.0;
322  G4double pzomc=0.0;
323  G4double rni=0.0;
324  G4double aux=0.0;
325  for (size_t i=0;i<numberOfOscillators;i++)
326  {
327  ionEnergy = (*theTable)[i]->GetIonisationEnergy();
328  if (photonEnergy0 > ionEnergy)
329  {
330  G4double aux2 = photonEnergy0*(photonEnergy0-ionEnergy)*2.0;
331  hartreeFunc = (*theTable)[i]->GetHartreeFactor();
332  oscStren = (*theTable)[i]->GetOscillatorStrength();
333  pzomc = hartreeFunc*(aux2-electron_mass_c2*ionEnergy)/
334  (electron_mass_c2*std::sqrt(2.0*aux2+ionEnergy*ionEnergy));
335  if (pzomc > 0)
336  rni = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)*
337  (std::sqrt(0.5)+std::sqrt(2.0)*pzomc));
338  else
339  rni = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)*
340  (std::sqrt(0.5)-std::sqrt(2.0)*pzomc));
341  s0 += oscStren*rni;
342  }
343  }
344  //Sampling tau
345  G4double cdt1 = 0.;
346  do
347  {
348  if ((G4UniformRand()*a2) < a1)
349  tau = std::pow(taumin,G4UniformRand());
350  else
351  tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0));
352  cdt1 = (1.0-tau)/(ek*tau);
353  //Incoherent scattering function
354  S = 0.;
355  for (size_t i=0;i<numberOfOscillators;i++)
356  {
357  ionEnergy = (*theTable)[i]->GetIonisationEnergy();
358  if (photonEnergy0 > ionEnergy) //sum only on excitable levels
359  {
360  aux = photonEnergy0*(photonEnergy0-ionEnergy)*cdt1;
361  hartreeFunc = (*theTable)[i]->GetHartreeFactor();
362  oscStren = (*theTable)[i]->GetOscillatorStrength();
363  pzomc = hartreeFunc*(aux-electron_mass_c2*ionEnergy)/
364  (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy));
365  if (pzomc > 0)
366  rn[i] = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)*
367  (std::sqrt(0.5)+std::sqrt(2.0)*pzomc));
368  else
369  rn[i] = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)*
370  (std::sqrt(0.5)-std::sqrt(2.0)*pzomc));
371  S += oscStren*rn[i];
372  pac[i] = S;
373  }
374  else
375  pac[i] = S-1e-6;
376  }
377  //Rejection function
378  TST = S*(1.0+tau*(ek1+tau*(ek2+tau*eks)))/(eks*tau*(1.0+tau*tau));
379  }while ((G4UniformRand()*s0) > TST);
380 
381  cosTheta = 1.0 - cdt1;
382  G4double fpzmax=0.0,fpz=0.0;
383  G4double A=0.0;
384  //Target electron shell
385  do
386  {
387  do
388  {
389  TST = S*G4UniformRand();
390  targetOscillator = numberOfOscillators-1; //last level
391  G4bool levelFound = false;
392  for (size_t i=0;i<numberOfOscillators && !levelFound;i++)
393  {
394  if (pac[i]>TST)
395  {
396  targetOscillator = i;
397  levelFound = true;
398  }
399  }
400  A = G4UniformRand()*rn[targetOscillator];
401  hartreeFunc = (*theTable)[targetOscillator]->GetHartreeFactor();
402  oscStren = (*theTable)[targetOscillator]->GetOscillatorStrength();
403  if (A < 0.5)
404  pzomc = (std::sqrt(0.5)-std::sqrt(0.5-std::log(2.0*A)))/
405  (std::sqrt(2.0)*hartreeFunc);
406  else
407  pzomc = (std::sqrt(0.5-std::log(2.0-2.0*A))-std::sqrt(0.5))/
408  (std::sqrt(2.0)*hartreeFunc);
409  } while (pzomc < -1);
410 
411  // F(EP) rejection
412  G4double XQC = 1.0+tau*(tau-2.0*cosTheta);
413  G4double AF = std::sqrt(XQC)*(1.0+tau*(tau-cosTheta)/XQC);
414  if (AF > 0)
415  fpzmax = 1.0+AF*0.2;
416  else
417  fpzmax = 1.0-AF*0.2;
418  fpz = 1.0+AF*std::max(std::min(pzomc,0.2),-0.2);
419  }while ((fpzmax*G4UniformRand())>fpz);
420 
421  //Energy of the scattered photon
422  G4double T = pzomc*pzomc;
423  G4double b1 = 1.0-T*tau*tau;
424  G4double b2 = 1.0-T*tau*cosTheta;
425  if (pzomc > 0.0)
426  epsilon = (tau/b1)*(b2+std::sqrt(std::abs(b2*b2-b1*(1.0-T))));
427  else
428  epsilon = (tau/b1)*(b2-std::sqrt(std::abs(b2*b2-b1*(1.0-T))));
429  } //energy < 5 MeV
430 
431  //Ok, the kinematics has been calculated.
432  G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
433  G4double phi = twopi * G4UniformRand() ;
434  G4double dirx = sinTheta * std::cos(phi);
435  G4double diry = sinTheta * std::sin(phi);
436  G4double dirz = cosTheta ;
437 
438  // Update G4VParticleChange for the scattered photon
439  G4ThreeVector photonDirection1(dirx,diry,dirz);
440  photonDirection1.rotateUz(photonDirection0);
441  fParticleChange->ProposeMomentumDirection(photonDirection1) ;
442 
443  G4double photonEnergy1 = epsilon * photonEnergy0;
444 
445  if (photonEnergy1 > 0.)
446  fParticleChange->SetProposedKineticEnergy(photonEnergy1) ;
447  else
448  {
451  }
452 
453  //Create scattered electron
454  G4double diffEnergy = photonEnergy0*(1-epsilon);
455  ionEnergy = (*theTable)[targetOscillator]->GetIonisationEnergy();
456 
457  G4double Q2 =
458  photonEnergy0*photonEnergy0+photonEnergy1*(photonEnergy1-2.0*photonEnergy0*cosTheta);
459  G4double cosThetaE = 0.; //scattering angle for the electron
460 
461  if (Q2 > 1.0e-12)
462  cosThetaE = (photonEnergy0-photonEnergy1*cosTheta)/std::sqrt(Q2);
463  else
464  cosThetaE = 1.0;
465  G4double sinThetaE = std::sqrt(1-cosThetaE*cosThetaE);
466 
467  //Now, try to handle fluorescence
468  //Notice: merged levels are indicated with Z=0 and flag=30
469  G4int shFlag = (*theTable)[targetOscillator]->GetShellFlag();
470  G4int Z = (G4int) (*theTable)[targetOscillator]->GetParentZ();
471 
472  //initialize here, then check photons created by Atomic-Deexcitation, and the final state e-
474  const G4AtomicShell* shell = 0;
475 
476  //Real level
477  if (Z > 0 && shFlag<30)
478  {
479  shell = fTransitionManager->Shell(Z,shFlag-1);
480  bindingEnergy = shell->BindingEnergy();
481  }
482 
483  G4double ionEnergyInPenelopeDatabase = ionEnergy;
484  //protection against energy non-conservation
485  ionEnergy = std::max(bindingEnergy,ionEnergyInPenelopeDatabase);
486 
487  //subtract the excitation energy. If not emitted by fluorescence
488  //the ionization energy is deposited as local energy deposition
489  G4double eKineticEnergy = diffEnergy - ionEnergy;
490  G4double localEnergyDeposit = ionEnergy;
491  G4double energyInFluorescence = 0.; //testing purposes only
492  G4double energyInAuger = 0; //testing purposes
493 
494  if (eKineticEnergy < 0)
495  {
496  //It means that there was some problem/mismatch between the two databases.
497  //Try to make it work
498  //In this case available Energy (diffEnergy) < ionEnergy
499  //Full residual energy is deposited locally
500  localEnergyDeposit = diffEnergy;
501  eKineticEnergy = 0.0;
502  }
503 
504  //the local energy deposit is what remains: part of this may be spent for fluorescence.
505  //Notice: shell might be NULL (invalid!) if shFlag=30. Must be protected
506  //Now, take care of fluorescence, if required
507  if (fAtomDeexcitation && shell)
508  {
509  G4int index = couple->GetIndex();
510  if (fAtomDeexcitation->CheckDeexcitationActiveRegion(index))
511  {
512  size_t nBefore = fvect->size();
513  fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
514  size_t nAfter = fvect->size();
515 
516  if (nAfter > nBefore) //actual production of fluorescence
517  {
518  for (size_t j=nBefore;j<nAfter;j++) //loop on products
519  {
520  G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
521  localEnergyDeposit -= itsEnergy;
522  if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
523  energyInFluorescence += itsEnergy;
524  else if (((*fvect)[j])->GetParticleDefinition() == G4Electron::Definition())
525  energyInAuger += itsEnergy;
526 
527  }
528  }
529 
530  }
531  }
532 
533 
534  /*
535  if(DeexcitationFlag() && Z > 5 && shellId>0) {
536 
537  const G4ProductionCutsTable* theCoupleTable=
538  G4ProductionCutsTable::GetProductionCutsTable();
539 
540  size_t index = couple->GetIndex();
541  G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[index];
542  G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[index];
543 
544  // Generation of fluorescence
545  // Data in EADL are available only for Z > 5
546  // Protection to avoid generating photons in the unphysical case of
547  // shell binding energy > photon energy
548  if (localEnergyDeposit > cutg || localEnergyDeposit > cute)
549  {
550  G4DynamicParticle* aPhoton;
551  deexcitationManager.SetCutForSecondaryPhotons(cutg);
552  deexcitationManager.SetCutForAugerElectrons(cute);
553 
554  photonVector = deexcitationManager.GenerateParticles(Z,shellId);
555  if(photonVector)
556  {
557  size_t nPhotons = photonVector->size();
558  for (size_t k=0; k<nPhotons; k++)
559  {
560  aPhoton = (*photonVector)[k];
561  if (aPhoton)
562  {
563  G4double itsEnergy = aPhoton->GetKineticEnergy();
564  G4bool keepIt = false;
565  if (itsEnergy <= localEnergyDeposit)
566  {
567  //check if good!
568  if(aPhoton->GetDefinition() == G4Gamma::Gamma()
569  && itsEnergy >= cutg)
570  {
571  keepIt = true;
572  energyInFluorescence += itsEnergy;
573  }
574  if (aPhoton->GetDefinition() == G4Electron::Electron() &&
575  itsEnergy >= cute)
576  {
577  energyInAuger += itsEnergy;
578  keepIt = true;
579  }
580  }
581  //good secondary, register it
582  if (keepIt)
583  {
584  localEnergyDeposit -= itsEnergy;
585  fvect->push_back(aPhoton);
586  }
587  else
588  {
589  delete aPhoton;
590  (*photonVector)[k] = 0;
591  }
592  }
593  }
594  delete photonVector;
595  }
596  }
597  }
598  */
599 
600 
601  //Always produce explicitely the electron
602  G4DynamicParticle* electron = 0;
603 
604  G4double xEl = sinThetaE * std::cos(phi+pi);
605  G4double yEl = sinThetaE * std::sin(phi+pi);
606  G4double zEl = cosThetaE;
607  G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction
608  eDirection.rotateUz(photonDirection0);
609  electron = new G4DynamicParticle (G4Electron::Electron(),
610  eDirection,eKineticEnergy) ;
611  fvect->push_back(electron);
612 
613 
614  if (localEnergyDeposit < 0)
615  {
616  G4cout << "WARNING-"
617  << "G4PenelopeComptonModel::SampleSecondaries - Negative energy deposit"
618  << G4endl;
619  localEnergyDeposit=0.;
620  }
621  fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
622 
623  G4double electronEnergy = 0.;
624  if (electron)
625  electronEnergy = eKineticEnergy;
626  if (verboseLevel > 1)
627  {
628  G4cout << "-----------------------------------------------------------" << G4endl;
629  G4cout << "Energy balance from G4PenelopeCompton" << G4endl;
630  G4cout << "Incoming photon energy: " << photonEnergy0/keV << " keV" << G4endl;
631  G4cout << "-----------------------------------------------------------" << G4endl;
632  G4cout << "Scattered photon: " << photonEnergy1/keV << " keV" << G4endl;
633  G4cout << "Scattered electron " << electronEnergy/keV << " keV" << G4endl;
634  if (energyInFluorescence)
635  G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
636  if (energyInAuger)
637  G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
638  G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
639  G4cout << "Total final state: " << (photonEnergy1+electronEnergy+energyInFluorescence+
640  localEnergyDeposit+energyInAuger)/keV <<
641  " keV" << G4endl;
642  G4cout << "-----------------------------------------------------------" << G4endl;
643  }
644  if (verboseLevel > 0)
645  {
646  G4double energyDiff = std::fabs(photonEnergy1+
647  electronEnergy+energyInFluorescence+
648  localEnergyDeposit+energyInAuger-photonEnergy0);
649  if (energyDiff > 0.05*keV)
650  G4cout << "Warning from G4PenelopeCompton: problem with energy conservation: " <<
651  (photonEnergy1+electronEnergy+energyInFluorescence+energyInAuger+localEnergyDeposit)/keV <<
652  " keV (final) vs. " <<
653  photonEnergy0/keV << " keV (initial)" << G4endl;
654  }
655 }
656 
657 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
658 
659 G4double G4PenelopeComptonModel::DifferentialCrossSection(G4double cosTheta,G4double energy,
661 {
662  //
663  // Penelope model v2008. Single differential cross section *per electron*
664  // for photon Compton scattering by
665  // electrons in the given atomic oscillator, differential in the direction of the
666  // scattering photon. This is in units of pi*classic_electr_radius**2
667  //
668  // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167
669  // The parametrization includes the J(p) distribution profiles for the atomic shells,
670  // that are tabulated from Hartree-Fock calculations
671  // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201
672  //
673  G4double ionEnergy = osc->GetIonisationEnergy();
674  G4double harFunc = osc->GetHartreeFactor();
675 
676  const G4double k2 = std::sqrt(2.);
677  const G4double k1 = 1./k2;
678 
679  if (energy < ionEnergy)
680  return 0;
681 
682  //energy of the Compton line
683  G4double cdt1 = 1.0-cosTheta;
684  G4double EOEC = 1.0+(energy/electron_mass_c2)*cdt1;
685  G4double ECOE = 1.0/EOEC;
686 
687  //Incoherent scattering function (analytical profile)
688  G4double aux = energy*(energy-ionEnergy)*cdt1;
689  G4double Pzimax =
690  (aux - electron_mass_c2*ionEnergy)/(electron_mass_c2*std::sqrt(2*aux+ionEnergy*ionEnergy));
691  G4double sia = 0.0;
692  G4double x = harFunc*Pzimax;
693  if (x > 0)
694  sia = 1.0-0.5*std::exp(0.5-(k1+k2*x)*(k1+k2*x));
695  else
696  sia = 0.5*std::exp(0.5-(k1-k2*x)*(k1-k2*x));
697 
698  //1st order correction, integral of Pz times the Compton profile.
699  //Calculated approximately using a free-electron gas profile
700  G4double pf = 3.0/(4.0*harFunc);
701  if (std::fabs(Pzimax) < pf)
702  {
703  G4double QCOE2 = 1.0+ECOE*ECOE-2.0*ECOE*cosTheta;
704  G4double p2 = Pzimax*Pzimax;
705  G4double dspz = std::sqrt(QCOE2)*
706  (1.0+ECOE*(ECOE-cosTheta)/QCOE2)*harFunc
707  *0.25*(2*p2-(p2*p2)/(pf*pf)-(pf*pf));
708  sia += std::max(dspz,-1.0*sia);
709  }
710 
711  G4double XKN = EOEC+ECOE-1.0+cosTheta*cosTheta;
712 
713  //Differential cross section (per electron, in units of pi*classic_electr_radius**2)
714  G4double diffCS = ECOE*ECOE*XKN*sia;
715 
716  return diffCS;
717 }
718 
719 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
720 
721 G4double G4PenelopeComptonModel::OscillatorTotalCrossSection(G4double energy,G4PenelopeOscillator* osc)
722 {
723  //Total cross section (integrated) for the given oscillator in units of
724  //pi*classic_electr_radius^2
725 
726  //Integrate differential cross section for each oscillator
727  G4double stre = osc->GetOscillatorStrength();
728 
729  // here one uses the using the 20-point
730  // Gauss quadrature method with an adaptive bipartition scheme
731  const G4int npoints=10;
732  const G4int ncallsmax=20000;
733  const G4int nst=256;
734  static G4double Abscissas[10] = {7.652651133497334e-02,2.2778585114164508e-01,3.7370608871541956e-01,
735  5.1086700195082710e-01,6.3605368072651503e-01,7.4633190646015079e-01,
736  8.3911697182221882e-01,9.1223442825132591e-01,9.6397192727791379e-01,
737  9.9312859918509492e-01};
738  static G4double Weights[10] = {1.5275338713072585e-01,1.4917298647260375e-01,1.4209610931838205e-01,
739  1.3168863844917663e-01,1.1819453196151842e-01,1.0193011981724044e-01,
740  8.3276741576704749e-02,6.2672048334109064e-02,4.0601429800386941e-02,
741  1.7614007139152118e-02};
742 
743  G4double MaxError = 1e-5;
744  //Error control
745  G4double Ctol = std::min(std::max(MaxError,1e-13),1e-02);
746  G4double Ptol = 0.01*Ctol;
747  G4double Err=1e35;
748 
749  //Gauss integration from -1 to 1
750  G4double LowPoint = -1.0;
751  G4double HighPoint = 1.0;
752 
753  G4double h=HighPoint-LowPoint;
754  G4double sumga=0.0;
755  G4double a=0.5*(HighPoint-LowPoint);
756  G4double b=0.5*(HighPoint+LowPoint);
757  G4double c=a*Abscissas[0];
758  G4double d= Weights[0]*
759  (DifferentialCrossSection(b+c,energy,osc)+DifferentialCrossSection(b-c,energy,osc));
760  for (G4int i=2;i<=npoints;i++)
761  {
762  c=a*Abscissas[i-1];
763  d += Weights[i-1]*
764  (DifferentialCrossSection(b+c,energy,osc)+DifferentialCrossSection(b-c,energy,osc));
765  }
766  G4int icall = 2*npoints;
767  G4int LH=1;
768  G4double S[nst],x[nst],sn[nst],xrn[nst];
769  S[0]=d*a;
770  x[0]=LowPoint;
771 
772  G4bool loopAgain = true;
773 
774  //Adaptive bipartition scheme
775  do{
776  G4double h0=h;
777  h=0.5*h; //bipartition
778  G4double sumr=0;
779  G4int LHN=0;
780  G4double si,xa,xb,xc;
781  for (G4int i=1;i<=LH;i++){
782  si=S[i-1];
783  xa=x[i-1];
784  xb=xa+h;
785  xc=xa+h0;
786  a=0.5*(xb-xa);
787  b=0.5*(xb+xa);
788  c=a*Abscissas[0];
789  G4double dLocal = Weights[0]*
790  (DifferentialCrossSection(b+c,energy,osc)+DifferentialCrossSection(b-c,energy,osc));
791 
792  for (G4int j=1;j<npoints;j++)
793  {
794  c=a*Abscissas[j];
795  dLocal += Weights[j]*
796  (DifferentialCrossSection(b+c,energy,osc)+DifferentialCrossSection(b-c,energy,osc));
797  }
798  G4double s1=dLocal*a;
799  a=0.5*(xc-xb);
800  b=0.5*(xc+xb);
801  c=a*Abscissas[0];
802  dLocal=Weights[0]*
803  (DifferentialCrossSection(b+c,energy,osc)+DifferentialCrossSection(b-c,energy,osc));
804 
805  for (G4int j=1;j<npoints;j++)
806  {
807  c=a*Abscissas[j];
808  dLocal += Weights[j]*
809  (DifferentialCrossSection(b+c,energy,osc)+DifferentialCrossSection(b-c,energy,osc));
810  }
811  G4double s2=dLocal*a;
812  icall=icall+4*npoints;
813  G4double s12=s1+s2;
814  if (std::abs(s12-si)<std::max(Ptol*std::abs(s12),1e-35))
815  sumga += s12;
816  else
817  {
818  sumr += s12;
819  LHN += 2;
820  sn[LHN-1]=s2;
821  xrn[LHN-1]=xb;
822  sn[LHN-2]=s1;
823  xrn[LHN-2]=xa;
824  }
825 
826  if (icall>ncallsmax || LHN>nst)
827  {
828  G4cout << "G4PenelopeComptonModel: " << G4endl;
829  G4cout << "LowPoint: " << LowPoint << ", High Point: " << HighPoint << G4endl;
830  G4cout << "Tolerance: " << MaxError << G4endl;
831  G4cout << "Calls: " << icall << ", Integral: " << sumga << ", Error: " << Err << G4endl;
832  G4cout << "Number of open subintervals: " << LHN << G4endl;
833  G4cout << "WARNING: the required accuracy has not been attained" << G4endl;
834  loopAgain = false;
835  }
836  }
837  Err=std::abs(sumr)/std::max(std::abs(sumr+sumga),1e-35);
838  if (Err < Ctol || LHN == 0)
839  loopAgain = false; //end of cycle
840  LH=LHN;
841  for (G4int i=0;i<LH;i++)
842  {
843  S[i]=sn[i];
844  x[i]=xrn[i];
845  }
846  }while(Ctol < 1.0 && loopAgain);
847 
848 
849  G4double xs = stre*sumga;
850 
851  return xs;
852 }
853 
854 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
855 
856 G4double G4PenelopeComptonModel::KleinNishinaCrossSection(G4double energy,
857  const G4Material* material)
858 {
859  // use Klein-Nishina formula
860  // total cross section in units of pi*classic_electr_radius^2
861 
862  G4double cs = 0;
863 
864  G4double ek =energy/electron_mass_c2;
865  G4double eks = ek*ek;
866  G4double ek2 = 1.0+ek+ek;
867  G4double ek1 = eks-ek2-1.0;
868 
869  G4double t0 = 1.0/ek2;
870  G4double csl = 0.5*eks*t0*t0+ek2*t0+ek1*std::log(t0)-(1.0/t0);
871 
872  G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableCompton(material);
873 
874  for (size_t i=0;i<theTable->size();i++)
875  {
876  G4PenelopeOscillator* theOsc = (*theTable)[i];
877  G4double ionEnergy = theOsc->GetIonisationEnergy();
878  G4double tau=(energy-ionEnergy)/energy;
879  if (tau > t0)
880  {
881  G4double csu = 0.5*eks*tau*tau+ek2*tau+ek1*std::log(tau)-(1.0/tau);
882  G4double stre = theOsc->GetOscillatorStrength();
883 
884  cs += stre*(csu-csl);
885  }
886  }
887 
888  cs /= (ek*eks);
889 
890  return cs;
891 
892 }
893