Geant4_10
G4PenelopeIonisationModel.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: G4PenelopeIonisationModel.cc 75573 2013-11-04 11:48:15Z gcosmo $
27 //
28 // Author: Luciano Pandola
29 //
30 // History:
31 // --------
32 // 27 Jul 2010 L Pandola First complete implementation
33 // 18 Jan 2011 L.Pandola Stricter check on production of sub-treshold delta-rays.
34 // Should never happen now
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 // 25 May 2011 L Pandola Renamed (make v2008 as default Penelope)
39 // 26 Jan 2012 L Pandola Migration of AtomicDeexcitation to the new interface
40 // 09 Mar 2012 L Pandola Moved the management and calculation of
41 // cross sections to a separate class. Use a different method to
42 // get normalized shell cross sections
43 // 07 Oct 2013 L. Pandola Migration to MT
44 //
45 
47 #include "G4PhysicalConstants.hh"
48 #include "G4SystemOfUnits.hh"
49 #include "G4ParticleDefinition.hh"
50 #include "G4MaterialCutsCouple.hh"
51 #include "G4ProductionCutsTable.hh"
52 #include "G4DynamicParticle.hh"
54 #include "G4AtomicShell.hh"
55 #include "G4Gamma.hh"
56 #include "G4Electron.hh"
57 #include "G4Positron.hh"
59 #include "G4PenelopeOscillator.hh"
61 #include "G4PhysicsFreeVector.hh"
62 #include "G4PhysicsLogVector.hh"
63 #include "G4LossTableManager.hh"
65 #include "G4AutoLock.hh"
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68 namespace { G4Mutex PenelopeIonisationModelMutex = G4MUTEX_INITIALIZER; }
69 
71  const G4String& nam)
72  :G4VEmModel(nam),fParticleChange(0),fParticle(0),isInitialised(false),
73  fAtomDeexcitation(0),kineticEnergy1(0.*eV),
74  cosThetaPrimary(1.0),energySecondary(0.*eV),
75  cosThetaSecondary(0.0),targetOscillator(-1),
76  theCrossSectionHandler(0),fLocalTable(false)
77 {
78  fIntrinsicLowEnergyLimit = 100.0*eV;
79  fIntrinsicHighEnergyLimit = 100.0*GeV;
80  // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
81  SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
82  nBins = 200;
83 
84  if (part)
85  SetParticle(part);
86 
87  //
89  //
90  verboseLevel= 0;
91 
92  // Verbosity scale:
93  // 0 = nothing
94  // 1 = warning for energy non-conservation
95  // 2 = details of energy budget
96  // 3 = calculation of cross sections, file openings, sampling of atoms
97  // 4 = entering in methods
98 
99  // Atomic deexcitation model activated by default
100  SetDeexcitationFlag(true);
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
104 
106 {
107  if (IsMaster() || fLocalTable)
108  {
109  if (theCrossSectionHandler)
110  delete theCrossSectionHandler;
111  }
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
115 
117  const G4DataVector& theCuts)
118 {
119  if (verboseLevel > 3)
120  G4cout << "Calling G4PenelopeIonisationModel::Initialise()" << G4endl;
121 
122  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
123  //Issue warning if the AtomicDeexcitation has not been declared
124  if (!fAtomDeexcitation)
125  {
126  G4cout << G4endl;
127  G4cout << "WARNING from G4PenelopeIonisationModel " << G4endl;
128  G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
129  G4cout << "any fluorescence/Auger emission." << G4endl;
130  G4cout << "Please make sure this is intended" << G4endl;
131  }
132  SetParticle(particle);
133 
134  //Only the master model creates/manages the tables. All workers get the
135  //pointer to the table, and use it as readonly
136  if (IsMaster() && particle == fParticle)
137  {
138 
139  //Set the number of bins for the tables. 20 points per decade
140  nBins = (size_t) (20*std::log10(HighEnergyLimit()/LowEnergyLimit()));
141  nBins = std::max(nBins,(size_t)100);
142 
143  //Clear and re-build the tables
144  if (theCrossSectionHandler)
145  {
146  delete theCrossSectionHandler;
147  theCrossSectionHandler = 0;
148  }
149  theCrossSectionHandler = new G4PenelopeIonisationXSHandler(nBins);
150  theCrossSectionHandler->SetVerboseLevel(verboseLevel);
151 
152  //Build tables for all materials
153  G4ProductionCutsTable* theCoupleTable =
155  for (size_t i=0;i<theCoupleTable->GetTableSize();i++)
156  {
157  const G4Material* theMat =
158  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
159  //Forces the building of the cross section tables
160  theCrossSectionHandler->BuildXSTable(theMat,theCuts.at(i),particle,
161  IsMaster());
162 
163  }
164 
165  if (verboseLevel > 2) {
166  G4cout << "Penelope Ionisation model v2008 is initialized " << G4endl
167  << "Energy range: "
168  << LowEnergyLimit() / keV << " keV - "
169  << HighEnergyLimit() / GeV << " GeV. Using "
170  << nBins << " bins."
171  << G4endl;
172  }
173  }
174 
175  if(isInitialised)
176  return;
178  isInitialised = true;
179 
180  return;
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184 
186  G4VEmModel *masterModel)
187 {
188  if (verboseLevel > 3)
189  G4cout << "Calling G4PenelopeIonisationModel::InitialiseLocal()" << G4endl;
190 
191  //
192  //Check that particle matches: one might have multiple master models (e.g.
193  //for e+ and e-).
194  //
195  if (part == fParticle)
196  {
197  //Get the const table pointers from the master to the workers
198  const G4PenelopeIonisationModel* theModel =
199  static_cast<G4PenelopeIonisationModel*> (masterModel);
200 
201  //Copy pointers to the data tables
202  theCrossSectionHandler = theModel->theCrossSectionHandler;
203 
204  //copy data
205  nBins = theModel->nBins;
206 
207  //Same verbosity for all workers, as the master
208  verboseLevel = theModel->verboseLevel;
209  }
210 
211  return;
212 }
213 
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
216 
218  const G4ParticleDefinition*
219  theParticle,
221  G4double cutEnergy,
222  G4double)
223 {
224  // Penelope model v2008 to calculate the cross section for inelastic collisions above the
225  // threshold. It makes use of the Generalised Oscillator Strength (GOS) model from
226  // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567
227  //
228  // The total cross section is calculated analytically by taking
229  // into account the atomic oscillators coming into the play for a given threshold.
230  //
231  // For incident e- the maximum energy allowed for the delta-rays is energy/2.
232  // because particles are undistinghishable.
233  //
234  // The contribution is splitted in three parts: distant longitudinal collisions,
235  // distant transverse collisions and close collisions. Each term is described by
236  // its own analytical function.
237  // Fermi density correction is calculated analytically according to
238  // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1
239  //
240  if (verboseLevel > 3)
241  G4cout << "Calling CrossSectionPerVolume() of G4PenelopeIonisationModel" << G4endl;
242 
243  SetupForMaterial(theParticle, material, energy);
244 
245  G4double totalCross = 0.0;
246  G4double crossPerMolecule = 0.;
247 
248  //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
249  //not invoked
250  if (!theCrossSectionHandler)
251  {
252  //create a **thread-local** version of the table. Used only for G4EmCalculator and
253  //Unit Tests
254  fLocalTable = true;
255  theCrossSectionHandler = new G4PenelopeIonisationXSHandler(nBins);
256  }
257 
258  const G4PenelopeCrossSection* theXS =
259  theCrossSectionHandler->GetCrossSectionTableForCouple(theParticle,
260  material,
261  cutEnergy);
262  if (!theXS)
263  {
264  //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
265  //not filled up. This can happen in a UnitTest or via G4EmCalculator
267  ed << "Unable to retrieve the cross section table for " << theParticle->GetParticleName() <<
268  " in " << material->GetName() << ", cut = " << cutEnergy/keV << " keV " << G4endl;
269  ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
270  G4Exception("G4PenelopeIonisationModel::CrossSectionPerVolume()",
271  "em2038",JustWarning,ed);
272  //protect file reading via autolock
273  G4AutoLock lock(&PenelopeIonisationModelMutex);
274  theCrossSectionHandler->BuildXSTable(material,cutEnergy,theParticle);
275  lock.unlock();
276  //now it should be ok
277  theXS =
278  theCrossSectionHandler->GetCrossSectionTableForCouple(theParticle,
279  material,
280  cutEnergy);
281  }
282 
283 
284  if (theXS)
285  crossPerMolecule = theXS->GetHardCrossSection(energy);
286 
287  G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
288  G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
289 
290  if (verboseLevel > 3)
291  G4cout << "Material " << material->GetName() << " has " << atPerMol <<
292  "atoms per molecule" << G4endl;
293 
294  G4double moleculeDensity = 0.;
295  if (atPerMol)
296  moleculeDensity = atomDensity/atPerMol;
297 
298  G4double crossPerVolume = crossPerMolecule*moleculeDensity;
299 
300  if (verboseLevel > 2)
301  {
302  G4cout << "G4PenelopeIonisationModel " << G4endl;
303  G4cout << "Mean free path for delta emission > " << cutEnergy/keV << " keV at " <<
304  energy/keV << " keV = " << (1./crossPerVolume)/mm << " mm" << G4endl;
305  if (theXS)
306  totalCross = (theXS->GetTotalCrossSection(energy))*moleculeDensity;
307  G4cout << "Total free path for ionisation (no threshold) at " <<
308  energy/keV << " keV = " << (1./totalCross)/mm << " mm" << G4endl;
309  }
310  return crossPerVolume;
311 }
312 
313 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
314 
315 //This is a dummy method. Never inkoved by the tracking, it just issues
316 //a warning if one tries to get Cross Sections per Atom via the
317 //G4EmCalculator.
319  G4double,
320  G4double,
321  G4double,
322  G4double,
323  G4double)
324 {
325  G4cout << "*** G4PenelopeIonisationModel -- WARNING ***" << G4endl;
326  G4cout << "Penelope Ionisation model v2008 does not calculate cross section _per atom_ " << G4endl;
327  G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
328  G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
329  return 0;
330 }
331 
332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
333 
335  const G4ParticleDefinition* theParticle,
336  G4double kineticEnergy,
337  G4double cutEnergy)
338 {
339  // Penelope model v2008 to calculate the stopping power for soft inelastic collisions
340  // below the threshold. It makes use of the Generalised Oscillator Strength (GOS)
341  // model from
342  // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567
343  //
344  // The stopping power is calculated analytically using the dsigma/dW cross
345  // section from the GOS models, which includes separate contributions from
346  // distant longitudinal collisions, distant transverse collisions and
347  // close collisions. Only the atomic oscillators that come in the play
348  // (according to the threshold) are considered for the calculation.
349  // Differential cross sections have a different form for e+ and e-.
350  //
351  // Fermi density correction is calculated analytically according to
352  // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1
353 
354  if (verboseLevel > 3)
355  G4cout << "Calling ComputeDEDX() of G4PenelopeIonisationModel" << G4endl;
356 
357  //Either Initialize() was not called, or we are in a slave and InitializeLocal() was
358  //not invoked
359  if (!theCrossSectionHandler)
360  {
361  //create a **thread-local** version of the table. Used only for G4EmCalculator and
362  //Unit Tests
363  fLocalTable = true;
364  theCrossSectionHandler = new G4PenelopeIonisationXSHandler(nBins);
365  }
366 
367  const G4PenelopeCrossSection* theXS =
368  theCrossSectionHandler->GetCrossSectionTableForCouple(theParticle,material,
369  cutEnergy);
370 
371  if (!theXS)
372  {
373  //If we are here, it means that Initialize() was inkoved, but the MaterialTable was
374  //not filled up. This can happen in a UnitTest or via G4EmCalculator
376  ed << "Unable to retrieve the cross section table for " << theParticle->GetParticleName() <<
377  " in " << material->GetName() << ", cut = " << cutEnergy/keV << " keV " << G4endl;
378  ed << "This can happen only in Unit Tests or via G4EmCalculator" << G4endl;
379  G4Exception("G4PenelopeIonisationModel::ComputeDEDXPerVolume()",
380  "em2038",JustWarning,ed);
381  //protect file reading via autolock
382  G4AutoLock lock(&PenelopeIonisationModelMutex);
383  theCrossSectionHandler->BuildXSTable(material,cutEnergy,theParticle);
384  lock.unlock();
385  //now it should be ok
386  theXS =
387  theCrossSectionHandler->GetCrossSectionTableForCouple(theParticle,
388  material,
389  cutEnergy);
390  }
391 
392  G4double sPowerPerMolecule = 0.0;
393  if (theXS)
394  sPowerPerMolecule = theXS->GetSoftStoppingPower(kineticEnergy);
395 
396 
397  G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
398  G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
399 
400  G4double moleculeDensity = 0.;
401  if (atPerMol)
402  moleculeDensity = atomDensity/atPerMol;
403 
404  G4double sPowerPerVolume = sPowerPerMolecule*moleculeDensity;
405 
406  if (verboseLevel > 2)
407  {
408  G4cout << "G4PenelopeIonisationModel " << G4endl;
409  G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
410  kineticEnergy/keV << " keV = " <<
411  sPowerPerVolume/(keV/mm) << " keV/mm" << G4endl;
412  }
413  return sPowerPerVolume;
414 }
415 
416 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
417 
419  const G4MaterialCutsCouple*)
420 {
421  return fIntrinsicLowEnergyLimit;
422 }
423 
424 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
425 
426 void G4PenelopeIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
427  const G4MaterialCutsCouple* couple,
428  const G4DynamicParticle* aDynamicParticle,
429  G4double cutE, G4double)
430 {
431  // Penelope model v2008 to sample the final state following an hard inelastic interaction.
432  // It makes use of the Generalised Oscillator Strength (GOS) model from
433  // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567
434  //
435  // The GOS model is used to calculate the individual cross sections for all
436  // the atomic oscillators coming in the play, taking into account the three
437  // contributions (distant longitudinal collisions, distant transverse collisions and
438  // close collisions). Then the target shell and the interaction channel are
439  // sampled. Final state of the delta-ray (energy, angle) are generated according
440  // to the analytical distributions (dSigma/dW) for the selected interaction
441  // channels.
442  // For e-, the maximum energy for the delta-ray is initialEnergy/2. (because
443  // particles are indistinghusbable), while it is the full initialEnergy for
444  // e+.
445  // The efficiency of the random sampling algorithm (e.g. for close collisions)
446  // decreases when initial and cutoff energy increase (e.g. 87% for 10-keV primary
447  // and 1 keV threshold, 99% for 10-MeV primary and 10-keV threshold).
448  // Differential cross sections have a different form for e+ and e-.
449  //
450  // WARNING: The model provides an _average_ description of inelastic collisions.
451  // Anyway, the energy spectrum associated to distant excitations of a given
452  // atomic shell is approximated as a single resonance. The simulated energy spectra
453  // show _unphysical_ narrow peaks at energies that are multiple of the shell
454  // resonance energies. The spurious speaks are automatically smoothed out after
455  // multiple inelastic collisions.
456  //
457  // The model determines also the original shell from which the delta-ray is expelled,
458  // in order to produce fluorescence de-excitation (from G4DeexcitationManager)
459  //
460  // Fermi density correction is calculated analytically according to
461  // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1
462 
463  if (verboseLevel > 3)
464  G4cout << "Calling SamplingSecondaries() of G4PenelopeIonisationModel" << G4endl;
465 
466  G4double kineticEnergy0 = aDynamicParticle->GetKineticEnergy();
467  const G4ParticleDefinition* theParticle = aDynamicParticle->GetDefinition();
468 
469  if (kineticEnergy0 <= fIntrinsicLowEnergyLimit)
470  {
473  return ;
474  }
475 
476  const G4Material* material = couple->GetMaterial();
477  const G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableIonisation(material);
478 
479  G4ParticleMomentum particleDirection0 = aDynamicParticle->GetMomentumDirection();
480 
481  //Initialise final-state variables. The proper values will be set by the methods
482  // SampleFinalStateElectron() and SampleFinalStatePositron()
483  kineticEnergy1=kineticEnergy0;
484  cosThetaPrimary=1.0;
485  energySecondary=0.0;
486  cosThetaSecondary=1.0;
487  targetOscillator = -1;
488 
489  if (theParticle == G4Electron::Electron())
490  SampleFinalStateElectron(material,cutE,kineticEnergy0);
491  else if (theParticle == G4Positron::Positron())
492  SampleFinalStatePositron(material,cutE,kineticEnergy0);
493  else
494  {
496  ed << "Invalid particle " << theParticle->GetParticleName() << G4endl;
497  G4Exception("G4PenelopeIonisationModel::SamplingSecondaries()",
498  "em0001",FatalException,ed);
499 
500  }
501  if (energySecondary == 0) return;
502 
503  if (verboseLevel > 3)
504  {
505  G4cout << "G4PenelopeIonisationModel::SamplingSecondaries() for " <<
506  theParticle->GetParticleName() << G4endl;
507  G4cout << "Final eKin = " << kineticEnergy1 << " keV" << G4endl;
508  G4cout << "Final cosTheta = " << cosThetaPrimary << G4endl;
509  G4cout << "Delta-ray eKin = " << energySecondary << " keV" << G4endl;
510  G4cout << "Delta-ray cosTheta = " << cosThetaSecondary << G4endl;
511  G4cout << "Oscillator: " << targetOscillator << G4endl;
512  }
513 
514  //Update the primary particle
515  G4double sint = std::sqrt(1. - cosThetaPrimary*cosThetaPrimary);
516  G4double phiPrimary = twopi * G4UniformRand();
517  G4double dirx = sint * std::cos(phiPrimary);
518  G4double diry = sint * std::sin(phiPrimary);
519  G4double dirz = cosThetaPrimary;
520 
521  G4ThreeVector electronDirection1(dirx,diry,dirz);
522  electronDirection1.rotateUz(particleDirection0);
523 
524  if (kineticEnergy1 > 0)
525  {
526  fParticleChange->ProposeMomentumDirection(electronDirection1);
527  fParticleChange->SetProposedKineticEnergy(kineticEnergy1);
528  }
529  else
531 
532 
533  //Generate the delta ray
534  G4double ionEnergyInPenelopeDatabase =
535  (*theTable)[targetOscillator]->GetIonisationEnergy();
536 
537  //Now, try to handle fluorescence
538  //Notice: merged levels are indicated with Z=0 and flag=30
539  G4int shFlag = (*theTable)[targetOscillator]->GetShellFlag();
540  G4int Z = (G4int) (*theTable)[targetOscillator]->GetParentZ();
541 
542  //initialize here, then check photons created by Atomic-Deexcitation, and the final state e-
545  //G4int shellId = 0;
546 
547  const G4AtomicShell* shell = 0;
548  //Real level
549  if (Z > 0 && shFlag<30)
550  {
551  shell = transitionManager->Shell(Z,shFlag-1);
552  bindingEnergy = shell->BindingEnergy();
553  //shellId = shell->ShellId();
554  }
555 
556  //correct the energySecondary to account for the fact that the Penelope
557  //database of ionisation energies is in general (slightly) different
558  //from the fluorescence database used in Geant4.
559  energySecondary += ionEnergyInPenelopeDatabase-bindingEnergy;
560 
561  G4double localEnergyDeposit = bindingEnergy;
562  //testing purposes only
563  G4double energyInFluorescence = 0;
564  G4double energyInAuger = 0;
565 
566  if (energySecondary < 0)
567  {
568  //It means that there was some problem/mismatch between the two databases.
569  //In this case, the available energy is ok to excite the level according
570  //to the Penelope database, but not according to the Geant4 database
571  //Full residual energy is deposited locally
572  localEnergyDeposit += energySecondary;
573  energySecondary = 0.0;
574  }
575 
576  //Notice: shell might be NULL (invalid!) if shFlag=30. Must be protected
577  //Now, take care of fluorescence, if required
578  if (fAtomDeexcitation && shell)
579  {
580  G4int index = couple->GetIndex();
581  if (fAtomDeexcitation->CheckDeexcitationActiveRegion(index))
582  {
583  size_t nBefore = fvect->size();
584  fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
585  size_t nAfter = fvect->size();
586 
587  if (nAfter > nBefore) //actual production of fluorescence
588  {
589  for (size_t j=nBefore;j<nAfter;j++) //loop on products
590  {
591  G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
592  localEnergyDeposit -= itsEnergy;
593  if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
594  energyInFluorescence += itsEnergy;
595  else if (((*fvect)[j])->GetParticleDefinition() == G4Electron::Definition())
596  energyInAuger += itsEnergy;
597  }
598  }
599  }
600  }
601 
602  // Generate the delta ray --> to be done only if above cut
603  if (energySecondary > cutE)
604  {
606  G4double sinThetaE = std::sqrt(1.-cosThetaSecondary*cosThetaSecondary);
607  G4double phiEl = phiPrimary+pi; //pi with respect to the primary electron/positron
608  G4double xEl = sinThetaE * std::cos(phiEl);
609  G4double yEl = sinThetaE * std::sin(phiEl);
610  G4double zEl = cosThetaSecondary;
611  G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction
612  eDirection.rotateUz(particleDirection0);
613  electron = new G4DynamicParticle (G4Electron::Electron(),
614  eDirection,energySecondary) ;
615  fvect->push_back(electron);
616  }
617  else
618  {
619  localEnergyDeposit += energySecondary;
620  energySecondary = 0;
621  }
622 
623  if (localEnergyDeposit < 0)
624  {
625  G4cout << "WARNING-"
626  << "G4PenelopeIonisationModel::SampleSecondaries - Negative energy deposit"
627  << G4endl;
628  localEnergyDeposit=0.;
629  }
630  fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
631 
632  if (verboseLevel > 1)
633  {
634  G4cout << "-----------------------------------------------------------" << G4endl;
635  G4cout << "Energy balance from G4PenelopeIonisation" << G4endl;
636  G4cout << "Incoming primary energy: " << kineticEnergy0/keV << " keV" << G4endl;
637  G4cout << "-----------------------------------------------------------" << G4endl;
638  G4cout << "Outgoing primary energy: " << kineticEnergy1/keV << " keV" << G4endl;
639  G4cout << "Delta ray " << energySecondary/keV << " keV" << G4endl;
640  if (energyInFluorescence)
641  G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
642  if (energyInAuger)
643  G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
644  G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
645  G4cout << "Total final state: " << (energySecondary+energyInFluorescence+kineticEnergy1+
646  localEnergyDeposit+energyInAuger)/keV <<
647  " keV" << G4endl;
648  G4cout << "-----------------------------------------------------------" << G4endl;
649  }
650 
651  if (verboseLevel > 0)
652  {
653  G4double energyDiff = std::fabs(energySecondary+energyInFluorescence+kineticEnergy1+
654  localEnergyDeposit+energyInAuger-kineticEnergy0);
655  if (energyDiff > 0.05*keV)
656  G4cout << "Warning from G4PenelopeIonisation: problem with energy conservation: " <<
657  (energySecondary+energyInFluorescence+kineticEnergy1+localEnergyDeposit+energyInAuger)/keV <<
658  " keV (final) vs. " <<
659  kineticEnergy0/keV << " keV (initial)" << G4endl;
660  }
661 
662 }
663 
664 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
665 void G4PenelopeIonisationModel::SampleFinalStateElectron(const G4Material* mat,
666  G4double cutEnergy,
667  G4double kineticEnergy)
668 {
669  // This method sets the final ionisation parameters
670  // kineticEnergy1, cosThetaPrimary (= updates of the primary e-)
671  // energySecondary, cosThetaSecondary (= info of delta-ray)
672  // targetOscillator (= ionised oscillator)
673  //
674  // The method implements SUBROUTINE EINa of Penelope
675  //
676 
677  G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableIonisation(mat);
678  size_t numberOfOscillators = theTable->size();
679  const G4PenelopeCrossSection* theXS =
680  theCrossSectionHandler->GetCrossSectionTableForCouple(G4Electron::Electron(),mat,
681  cutEnergy);
682  G4double delta = theCrossSectionHandler->GetDensityCorrection(mat,kineticEnergy);
683 
684  // Selection of the active oscillator
685  G4double TST = G4UniformRand();
686  targetOscillator = numberOfOscillators-1; //initialization, last oscillator
687  G4double XSsum = 0.;
688 
689  for (size_t i=0;i<numberOfOscillators-1;i++)
690  {
691  /* testing purposes
692  G4cout << "sampling: " << i << " " << XSsum << " " << TST << " " <<
693  theXS->GetShellCrossSection(i,kineticEnergy) << " " <<
694  theXS->GetNormalizedShellCrossSection(i,kineticEnergy) << " " <<
695  mat->GetName() << G4endl;
696  */
697  XSsum += theXS->GetNormalizedShellCrossSection(i,kineticEnergy);
698 
699  if (XSsum > TST)
700  {
701  targetOscillator = (G4int) i;
702  break;
703  }
704  }
705 
706 
707  if (verboseLevel > 3)
708  {
709  G4cout << "SampleFinalStateElectron: sampled oscillator #" << targetOscillator << "." << G4endl;
710  G4cout << "Ionisation energy: " << (*theTable)[targetOscillator]->GetIonisationEnergy()/eV <<
711  " eV " << G4endl;
712  G4cout << "Resonance energy: : " << (*theTable)[targetOscillator]->GetResonanceEnergy()/eV << " eV "
713  << G4endl;
714  }
715  //Constants
716  G4double rb = kineticEnergy + 2.0*electron_mass_c2;
717  G4double gam = 1.0+kineticEnergy/electron_mass_c2;
718  G4double gam2 = gam*gam;
719  G4double beta2 = (gam2-1.0)/gam2;
720  G4double amol = ((gam-1.0)/gam)*((gam-1.0)/gam);
721 
722  //Partial cross section of the active oscillator
723  G4double resEne = (*theTable)[targetOscillator]->GetResonanceEnergy();
724  G4double invResEne = 1.0/resEne;
725  G4double ionEne = (*theTable)[targetOscillator]->GetIonisationEnergy();
726  G4double cutoffEne = (*theTable)[targetOscillator]->GetCutoffRecoilResonantEnergy();
727  G4double XHDL = 0.;
728  G4double XHDT = 0.;
729  G4double QM = 0.;
730  G4double cps = 0.;
731  G4double cp = 0.;
732 
733  //Distant excitations
734  if (resEne > cutEnergy && resEne < kineticEnergy)
735  {
736  cps = kineticEnergy*rb;
737  cp = std::sqrt(cps);
738  G4double XHDT0 = std::max(std::log(gam2)-beta2-delta,0.);
739  if (resEne > 1.0e-6*kineticEnergy)
740  {
741  G4double cpp = std::sqrt((kineticEnergy-resEne)*(kineticEnergy-resEne+2.0*electron_mass_c2));
742  QM = std::sqrt((cp-cpp)*(cp-cpp)+electron_mass_c2*electron_mass_c2)-electron_mass_c2;
743  }
744  else
745  {
746  QM = resEne*resEne/(beta2*2.0*electron_mass_c2);
747  QM *= (1.0-QM*0.5/electron_mass_c2);
748  }
749  if (QM < cutoffEne)
750  {
751  XHDL = std::log(cutoffEne*(QM+2.0*electron_mass_c2)/(QM*(cutoffEne+2.0*electron_mass_c2)))
752  *invResEne;
753  XHDT = XHDT0*invResEne;
754  }
755  else
756  {
757  QM = cutoffEne;
758  XHDL = 0.;
759  XHDT = 0.;
760  }
761  }
762  else
763  {
764  QM = cutoffEne;
765  cps = 0.;
766  cp = 0.;
767  XHDL = 0.;
768  XHDT = 0.;
769  }
770 
771  //Close collisions
772  G4double EE = kineticEnergy + ionEne;
773  G4double wmaxc = 0.5*EE;
774  G4double wcl = std::max(cutEnergy,cutoffEne);
775  G4double rcl = wcl/EE;
776  G4double XHC = 0.;
777  if (wcl < wmaxc)
778  {
779  G4double rl1 = 1.0-rcl;
780  G4double rrl1 = 1.0/rl1;
781  XHC = (amol*(0.5-rcl)+1.0/rcl-rrl1+
782  (1.0-amol)*std::log(rcl*rrl1))/EE;
783  }
784 
785  //Total cross section per molecule for the active shell, in cm2
786  G4double XHTOT = XHC + XHDL + XHDT;
787 
788  //very small cross section, do nothing
789  if (XHTOT < 1.e-14*barn)
790  {
791  kineticEnergy1=kineticEnergy;
792  cosThetaPrimary=1.0;
793  energySecondary=0.0;
794  cosThetaSecondary=1.0;
795  targetOscillator = numberOfOscillators-1;
796  return;
797  }
798 
799  //decide which kind of interaction we'll have
800  TST = XHTOT*G4UniformRand();
801 
802 
803  // Hard close collision
804  G4double TS1 = XHC;
805 
806  if (TST < TS1)
807  {
808  G4double A = 5.0*amol;
809  G4double ARCL = A*0.5*rcl;
810  G4double rk=0.;
811  G4bool loopAgain = false;
812  do
813  {
814  loopAgain = false;
815  G4double fb = (1.0+ARCL)*G4UniformRand();
816  if (fb < 1)
817  rk = rcl/(1.0-fb*(1.0-(rcl+rcl)));
818  else
819  rk = rcl + (fb-1.0)*(0.5-rcl)/ARCL;
820  G4double rk2 = rk*rk;
821  G4double rkf = rk/(1.0-rk);
822  G4double phi = 1.0+rkf*rkf-rkf+amol*(rk2+rkf);
823  if (G4UniformRand()*(1.0+A*rk2) > phi)
824  loopAgain = true;
825  }while(loopAgain);
826  //energy and scattering angle (primary electron)
827  G4double deltaE = rk*EE;
828 
829  kineticEnergy1 = kineticEnergy - deltaE;
830  cosThetaPrimary = std::sqrt(kineticEnergy1*rb/(kineticEnergy*(rb-deltaE)));
831  //energy and scattering angle of the delta ray
832  energySecondary = deltaE - ionEne; //subtract ionisation energy
833  cosThetaSecondary= std::sqrt(deltaE*rb/(kineticEnergy*(deltaE+2.0*electron_mass_c2)));
834  if (verboseLevel > 3)
835  G4cout << "SampleFinalStateElectron: sampled close collision " << G4endl;
836  return;
837  }
838 
839  //Hard distant longitudinal collisions
840  TS1 += XHDL;
841  G4double deltaE = resEne;
842  kineticEnergy1 = kineticEnergy - deltaE;
843 
844  if (TST < TS1)
845  {
846  G4double QS = QM/(1.0+QM*0.5/electron_mass_c2);
847  G4double Q = QS/(std::pow((QS/cutoffEne)*(1.0+cutoffEne*0.5/electron_mass_c2),G4UniformRand())
848  - (QS*0.5/electron_mass_c2));
849  G4double QTREV = Q*(Q+2.0*electron_mass_c2);
850  G4double cpps = kineticEnergy1*(kineticEnergy1+2.0*electron_mass_c2);
851  cosThetaPrimary = (cpps+cps-QTREV)/(2.0*cp*std::sqrt(cpps));
852  if (cosThetaPrimary > 1.)
853  cosThetaPrimary = 1.0;
854  //energy and emission angle of the delta ray
855  energySecondary = deltaE - ionEne;
856  cosThetaSecondary = 0.5*(deltaE*(kineticEnergy+rb-deltaE)+QTREV)/std::sqrt(cps*QTREV);
857  if (cosThetaSecondary > 1.0)
858  cosThetaSecondary = 1.0;
859  if (verboseLevel > 3)
860  G4cout << "SampleFinalStateElectron: sampled distant longitudinal collision " << G4endl;
861  return;
862  }
863 
864  //Hard distant transverse collisions
865  cosThetaPrimary = 1.0;
866  //energy and emission angle of the delta ray
867  energySecondary = deltaE - ionEne;
868  cosThetaSecondary = 0.5;
869  if (verboseLevel > 3)
870  G4cout << "SampleFinalStateElectron: sampled distant transverse collision " << G4endl;
871 
872  return;
873 }
874 
875 
876 
877 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
878 void G4PenelopeIonisationModel::SampleFinalStatePositron(const G4Material* mat,
879  G4double cutEnergy,
880  G4double kineticEnergy)
881 {
882  // This method sets the final ionisation parameters
883  // kineticEnergy1, cosThetaPrimary (= updates of the primary e-)
884  // energySecondary, cosThetaSecondary (= info of delta-ray)
885  // targetOscillator (= ionised oscillator)
886  //
887  // The method implements SUBROUTINE PINa of Penelope
888  //
889 
890  G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableIonisation(mat);
891  size_t numberOfOscillators = theTable->size();
892  const G4PenelopeCrossSection* theXS =
893  theCrossSectionHandler->GetCrossSectionTableForCouple(G4Positron::Positron(),mat,
894  cutEnergy);
895  G4double delta = theCrossSectionHandler->GetDensityCorrection(mat,kineticEnergy);
896 
897  // Selection of the active oscillator
898  G4double TST = G4UniformRand();
899  targetOscillator = numberOfOscillators-1; //initialization, last oscillator
900  G4double XSsum = 0.;
901  for (size_t i=0;i<numberOfOscillators-1;i++)
902  {
903  XSsum += theXS->GetNormalizedShellCrossSection(i,kineticEnergy);
904  if (XSsum > TST)
905  {
906  targetOscillator = (G4int) i;
907  break;
908  }
909  }
910 
911  if (verboseLevel > 3)
912  {
913  G4cout << "SampleFinalStatePositron: sampled oscillator #" << targetOscillator << "." << G4endl;
914  G4cout << "Ionisation energy: " << (*theTable)[targetOscillator]->GetIonisationEnergy()/eV
915  << " eV " << G4endl;
916  G4cout << "Resonance energy: : " << (*theTable)[targetOscillator]->GetResonanceEnergy()/eV
917  << " eV " << G4endl;
918  }
919 
920  //Constants
921  G4double rb = kineticEnergy + 2.0*electron_mass_c2;
922  G4double gam = 1.0+kineticEnergy/electron_mass_c2;
923  G4double gam2 = gam*gam;
924  G4double beta2 = (gam2-1.0)/gam2;
925  G4double g12 = (gam+1.0)*(gam+1.0);
926  G4double amol = ((gam-1.0)/gam)*((gam-1.0)/gam);
927  //Bhabha coefficients
928  G4double bha1 = amol*(2.0*g12-1.0)/(gam2-1.0);
929  G4double bha2 = amol*(3.0+1.0/g12);
930  G4double bha3 = amol*2.0*gam*(gam-1.0)/g12;
931  G4double bha4 = amol*(gam-1.0)*(gam-1.0)/g12;
932 
933  //
934  //Partial cross section of the active oscillator
935  //
936  G4double resEne = (*theTable)[targetOscillator]->GetResonanceEnergy();
937  G4double invResEne = 1.0/resEne;
938  G4double ionEne = (*theTable)[targetOscillator]->GetIonisationEnergy();
939  G4double cutoffEne = (*theTable)[targetOscillator]->GetCutoffRecoilResonantEnergy();
940 
941  G4double XHDL = 0.;
942  G4double XHDT = 0.;
943  G4double QM = 0.;
944  G4double cps = 0.;
945  G4double cp = 0.;
946 
947  //Distant excitations XS (same as for electrons)
948  if (resEne > cutEnergy && resEne < kineticEnergy)
949  {
950  cps = kineticEnergy*rb;
951  cp = std::sqrt(cps);
952  G4double XHDT0 = std::max(std::log(gam2)-beta2-delta,0.);
953  if (resEne > 1.0e-6*kineticEnergy)
954  {
955  G4double cpp = std::sqrt((kineticEnergy-resEne)*(kineticEnergy-resEne+2.0*electron_mass_c2));
956  QM = std::sqrt((cp-cpp)*(cp-cpp)+electron_mass_c2*electron_mass_c2)-electron_mass_c2;
957  }
958  else
959  {
960  QM = resEne*resEne/(beta2*2.0*electron_mass_c2);
961  QM *= (1.0-QM*0.5/electron_mass_c2);
962  }
963  if (QM < cutoffEne)
964  {
965  XHDL = std::log(cutoffEne*(QM+2.0*electron_mass_c2)/(QM*(cutoffEne+2.0*electron_mass_c2)))
966  *invResEne;
967  XHDT = XHDT0*invResEne;
968  }
969  else
970  {
971  QM = cutoffEne;
972  XHDL = 0.;
973  XHDT = 0.;
974  }
975  }
976  else
977  {
978  QM = cutoffEne;
979  cps = 0.;
980  cp = 0.;
981  XHDL = 0.;
982  XHDT = 0.;
983  }
984  //Close collisions (Bhabha)
985  G4double wmaxc = kineticEnergy;
986  G4double wcl = std::max(cutEnergy,cutoffEne);
987  G4double rcl = wcl/kineticEnergy;
988  G4double XHC = 0.;
989  if (wcl < wmaxc)
990  {
991  G4double rl1 = 1.0-rcl;
992  XHC = ((1.0/rcl-1.0)+bha1*std::log(rcl)+bha2*rl1
993  + (bha3/2.0)*(rcl*rcl-1.0)
994  + (bha4/3.0)*(1.0-rcl*rcl*rcl))/kineticEnergy;
995  }
996 
997  //Total cross section per molecule for the active shell, in cm2
998  G4double XHTOT = XHC + XHDL + XHDT;
999 
1000  //very small cross section, do nothing
1001  if (XHTOT < 1.e-14*barn)
1002  {
1003  kineticEnergy1=kineticEnergy;
1004  cosThetaPrimary=1.0;
1005  energySecondary=0.0;
1006  cosThetaSecondary=1.0;
1007  targetOscillator = numberOfOscillators-1;
1008  return;
1009  }
1010 
1011  //decide which kind of interaction we'll have
1012  TST = XHTOT*G4UniformRand();
1013 
1014  // Hard close collision
1015  G4double TS1 = XHC;
1016  if (TST < TS1)
1017  {
1018  G4double rl1 = 1.0-rcl;
1019  G4double rk=0.;
1020  G4bool loopAgain = false;
1021  do
1022  {
1023  loopAgain = false;
1024  rk = rcl/(1.0-G4UniformRand()*rl1);
1025  G4double phi = 1.0-rk*(bha1-rk*(bha2-rk*(bha3-bha4*rk)));
1026  if (G4UniformRand() > phi)
1027  loopAgain = true;
1028  }while(loopAgain);
1029  //energy and scattering angle (primary electron)
1030  G4double deltaE = rk*kineticEnergy;
1031  kineticEnergy1 = kineticEnergy - deltaE;
1032  cosThetaPrimary = std::sqrt(kineticEnergy1*rb/(kineticEnergy*(rb-deltaE)));
1033  //energy and scattering angle of the delta ray
1034  energySecondary = deltaE - ionEne; //subtract ionisation energy
1035  cosThetaSecondary= std::sqrt(deltaE*rb/(kineticEnergy*(deltaE+2.0*electron_mass_c2)));
1036  if (verboseLevel > 3)
1037  G4cout << "SampleFinalStatePositron: sampled close collision " << G4endl;
1038  return;
1039  }
1040 
1041  //Hard distant longitudinal collisions
1042  TS1 += XHDL;
1043  G4double deltaE = resEne;
1044  kineticEnergy1 = kineticEnergy - deltaE;
1045  if (TST < TS1)
1046  {
1047  G4double QS = QM/(1.0+QM*0.5/electron_mass_c2);
1048  G4double Q = QS/(std::pow((QS/cutoffEne)*(1.0+cutoffEne*0.5/electron_mass_c2),G4UniformRand())
1049  - (QS*0.5/electron_mass_c2));
1050  G4double QTREV = Q*(Q+2.0*electron_mass_c2);
1051  G4double cpps = kineticEnergy1*(kineticEnergy1+2.0*electron_mass_c2);
1052  cosThetaPrimary = (cpps+cps-QTREV)/(2.0*cp*std::sqrt(cpps));
1053  if (cosThetaPrimary > 1.)
1054  cosThetaPrimary = 1.0;
1055  //energy and emission angle of the delta ray
1056  energySecondary = deltaE - ionEne;
1057  cosThetaSecondary = 0.5*(deltaE*(kineticEnergy+rb-deltaE)+QTREV)/std::sqrt(cps*QTREV);
1058  if (cosThetaSecondary > 1.0)
1059  cosThetaSecondary = 1.0;
1060  if (verboseLevel > 3)
1061  G4cout << "SampleFinalStatePositron: sampled distant longitudinal collision " << G4endl;
1062  return;
1063  }
1064 
1065  //Hard distant transverse collisions
1066  cosThetaPrimary = 1.0;
1067  //energy and emission angle of the delta ray
1068  energySecondary = deltaE - ionEne;
1069  cosThetaSecondary = 0.5;
1070 
1071  if (verboseLevel > 3)
1072  G4cout << "SampleFinalStatePositron: sampled distant transverse collision " << G4endl;
1073 
1074  return;
1075 }
1076 
1077 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...
1078 
1079 void G4PenelopeIonisationModel::SetParticle(const G4ParticleDefinition* p)
1080 {
1081  if(!fParticle) {
1082  fParticle = p;
1083  }
1084 }
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4PenelopeOscillatorTable * GetOscillatorTableIonisation(const G4Material *)
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
static G4LossTableManager * Instance()
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4double GetSoftStoppingPower(G4double energy) const
Returns the total stopping power due to soft collisions.
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
void SetVerboseLevel(G4int vl)
Setter for the verbosity level.
Int_t index
Definition: macro.C:9
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:176
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:380
static G4Electron * Definition()
Definition: G4Electron.cc:49
G4double GetHardCrossSection(G4double energy) const
Returns hard cross section at the given energy.
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4ParticleDefinition * GetDefinition() const
const G4PenelopeCrossSection * GetCrossSectionTableForCouple(const G4ParticleDefinition *, const G4Material *, const G4double cut) const
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *)
G4double BindingEnergy() const
int G4int
Definition: G4Types.hh:78
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:158
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
Float_t mat
Definition: plot.C:40
string material
Definition: eplot.py:19
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
double energy
Definition: plottest35.C:25
Float_t Z
Definition: plot.C:39
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
TString part[npart]
Definition: Style.C:32
G4double GetTotalCrossSection(G4double energy) const
Returns total cross section at the given energy.
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4ParticleChangeForLoss * fParticleChange
float electron_mass_c2
Definition: hepunit.py:274
static G4PenelopeOscillatorManager * GetOscillatorManager()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *)
G4int G4Mutex
Definition: G4Threading.hh:156
static G4Positron * Positron()
Definition: G4Positron.cc:94
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
std::vector< G4PenelopeOscillator * > G4PenelopeOscillatorTable
G4double GetDensityCorrection(const G4Material *, const G4double energy) const
Returns the density coeection for the material at the given energy.
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
void BuildXSTable(const G4Material *, G4double cut, const G4ParticleDefinition *, G4bool isMaster=true)
This can be inkoved only by the master.
double G4double
Definition: G4Types.hh:76
G4double GetNormalizedShellCrossSection(size_t shellID, G4double energy) const
Returns the hard cross section for the given shell (normalized to 1)
const G4ParticleDefinition * fParticle
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double GetAtomsPerMolecule(const G4Material *)
Returns the total number of atoms per molecule.
G4double bindingEnergy(G4int A, G4int Z)
G4PenelopeIonisationModel(const G4ParticleDefinition *p=0, const G4String &processName="PenIoni")
static G4AtomicTransitionManager * Instance()
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double, G4double, G4double, G4double, G4double)
static G4Gamma * Definition()
Definition: G4Gamma.cc:49
const G4Material * GetMaterial() const