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