Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UAtomicDeexcitation.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// -------------------------------------------------------------------
28//
29// Geant4 Class file
30//
31// Authors: Alfonso Mantero ([email protected])
32//
33// Created 22 April 2010 from old G4UAtomicDeexcitation class
34//
35// Modified:
36// ---------
37// 20 Oct 2011 Alf modified to take into account ECPSSR form Form Factor
38// 03 Nov 2011 Alf Extended Empirical and Form Factor ionisation XS models
39// out thei ranges with Analytical one.
40// 07 Nov 2011 Alf Restored original ioniation XS for alphas,
41// letting scaled ones for other ions.
42// 20 Mar 2012 LP Register G4PenelopeIonisationCrossSection
43//
44// -------------------------------------------------------------------
45//
46// Class description:
47// Implementation of atomic deexcitation
48//
49// -------------------------------------------------------------------
50
53#include "G4SystemOfUnits.hh"
54#include "Randomize.hh"
55#include "G4Gamma.hh"
57#include "G4FluoTransition.hh"
58#include "G4Electron.hh"
59#include "G4Positron.hh"
60#include "G4Proton.hh"
61#include "G4Alpha.hh"
62
63#include "G4teoCrossSection.hh"
64#include "G4empCrossSection.hh"
67#include "G4EmCorrections.hh"
68#include "G4LossTableManager.hh"
69#include "G4EmParameters.hh"
70#include "G4Material.hh"
71#include "G4AtomicShells.hh"
72
73using namespace std;
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
78 G4VAtomDeexcitation("UAtomDeexcitation"),
79 minGammaEnergy(DBL_MAX),
80 minElectronEnergy(DBL_MAX),
81 newShellId(-1)
82{
83 anaPIXEshellCS = nullptr;
84 PIXEshellCS = nullptr;
85 ePIXEshellCS = nullptr;
87 theElectron = G4Electron::Electron();
88 thePositron = G4Positron::Positron();
89 transitionManager = G4AtomicTransitionManager::Instance();
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
95{
96 delete anaPIXEshellCS;
97 delete PIXEshellCS;
98 delete ePIXEshellCS;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102
104{
105 if(!IsFluoActive()) { return; }
106 transitionManager->Initialise();
107 if(!IsPIXEActive()) { return; }
108
109 if(!anaPIXEshellCS) {
110 anaPIXEshellCS = new G4teoCrossSection("ECPSSR_Analytical");
111 }
112 G4cout << G4endl;
113 G4cout << "### === G4UAtomicDeexcitation::InitialiseForNewRun()" << G4endl;
114
116 G4String namePIXExsModel = param->PIXECrossSectionModel();
117 G4String namePIXExsElectronModel = param->PIXEElectronCrossSectionModel();
118
119 // Check if old cross section for p/ion should be deleted
120 if(PIXEshellCS && namePIXExsModel != PIXEshellCS->GetName())
121 {
122 delete PIXEshellCS;
123 PIXEshellCS = nullptr;
124 }
125
126 // Instantiate new proton/ion cross section
127 if(!PIXEshellCS) {
128 if (namePIXExsModel == "ECPSSR_FormFactor")
129 {
130 PIXEshellCS = new G4teoCrossSection(namePIXExsModel);
131 }
132 else if(namePIXExsModel == "ECPSSR_ANSTO")
133 {
134 PIXEshellCS = new G4teoCrossSection(namePIXExsModel);
135 }
136 else if(namePIXExsModel == "Empirical")
137 {
138 PIXEshellCS = new G4empCrossSection(namePIXExsModel);
139 }
140 }
141
142 // Check if old cross section for e+- should be deleted
143 if(ePIXEshellCS && namePIXExsElectronModel != ePIXEshellCS->GetName())
144 {
145 delete ePIXEshellCS;
146 ePIXEshellCS = nullptr;
147 }
148
149 // Instantiate new e+- cross section
150 if(nullptr == ePIXEshellCS)
151 {
152 if(namePIXExsElectronModel == "Empirical")
153 {
154 ePIXEshellCS = new G4empCrossSection("Empirical");
155 }
156 else if(namePIXExsElectronModel == "ECPSSR_Analytical")
157 {
158 ePIXEshellCS = new G4teoCrossSection("ECPSSR_Analytical");
159 }
160 else if (namePIXExsElectronModel == "Penelope")
161 {
162 ePIXEshellCS = new G4PenelopeIonisationCrossSection();
163 }
164 else
165 {
166 ePIXEshellCS = new G4LivermoreIonisationCrossSection();
167 }
168 }
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172
175
176//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177
178const G4AtomicShell*
180{
181 return transitionManager->Shell(Z, (std::size_t)shell);
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185
187 std::vector<G4DynamicParticle*>* vectorOfParticles,
188 const G4AtomicShell* atomicShell,
189 G4int Z,
190 G4double gammaCut,
191 G4double eCut)
192{
193 // Defined initial conditions
194 G4int givenShellId = atomicShell->ShellId();
195 minGammaEnergy = gammaCut;
196 minElectronEnergy = eCut;
197 vacancyArray.clear();
198
199 // generation secondaries
200 G4DynamicParticle* aParticle=0;
201 G4int provShellId = 0;
202
203 //ORIGINAL METHOD BY ALFONSO MANTERO
205 {
206 //----------------------------
207 G4int counter = 0;
208
209 // limits of the EPDL data
210 if (Z>5 && Z<105) {
211
212 // The aim of this loop is to generate more than one fluorecence photon
213 // from the same ionizing event
214 do
215 {
216 if (counter == 0)
217 // First call to GenerateParticles(...):
218 // givenShellId is given by the process
219 {
220 provShellId = SelectTypeOfTransition(Z, givenShellId);
221
222 if (provShellId >0)
223 {
224 aParticle =
225 GenerateFluorescence(Z, givenShellId, provShellId);
226 }
227 else if (provShellId == -1)
228 {
229 aParticle = GenerateAuger(Z, givenShellId);
230 }
231 }
232 else
233 // Following calls to GenerateParticles(...):
234 // newShellId is given by GenerateFluorescence(...)
235 {
236 provShellId = SelectTypeOfTransition(Z,newShellId);
237 if (provShellId >0)
238 {
239 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
240 }
241 else if ( provShellId == -1)
242 {
243 aParticle = GenerateAuger(Z, newShellId);
244 }
245 }
246 ++counter;
247 if (aParticle != 0)
248 {
249 vectorOfParticles->push_back(aParticle);
250 }
251 else {provShellId = -2;}
252 }
253 while (provShellId > -2);
254 }
255 } // Auger cascade is not active
256
257 //END OF ORIGINAL METHOD BY ALFONSO MANTERO
258 //----------------------
259
260 // NEW METHOD
261 // Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
263 {
264 //----------------------
265 vacancyArray.push_back(givenShellId);
266
267 // let's check that 5<Z<100
268 if (Z<6 || Z>104){
269 return;
270 }
271
272 // as long as there is vacancy to be filled by either fluo or auger, stay in the loop.
273 while(!vacancyArray.empty()){
274 // prepare to process the last element, and then delete it from the vector.
275 givenShellId = vacancyArray[0];
276 provShellId = SelectTypeOfTransition(Z,givenShellId);
277
278 //G4cout<<"\n------ Atom Transition with Z: "<<Z<<"\tbetween current:"
279 // <<givenShellId<<" & target:"<<provShellId<<G4endl;
280 if(provShellId>0){
281 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
282 }
283 else if(provShellId == -1){
284 aParticle = GenerateAuger(Z, givenShellId);
285 }
286 // if a particle is created, put it in the vector of new particles
287 if(aParticle!=0)
288 vectorOfParticles->push_back(aParticle);
289
290 // one vacancy has been processed. Erase it.
291 vacancyArray.erase(vacancyArray.begin());
292 }
293 //----------------------
294 //End of Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
295
296 } // Auger cascade is active
297}
298
299//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
300
303 const G4ParticleDefinition* pdef,
304 G4int Z,
305 G4AtomicShellEnumerator shellEnum,
306 G4double kineticEnergy,
307 const G4Material* mat)
308{
309 // we must put a control on the shell that are passed:
310 // some shells should not pass (line "0" or "2")
311
312 // check atomic number
313 G4double xsec = 0.0;
314 if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
315 G4int idx = G4int(shellEnum);
316 if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
317
318 if(pdef == theElectron || pdef == thePositron) {
319 xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
320 return xsec;
321 }
322
323 G4double mass = pdef->GetPDGMass();
324 G4double escaled = kineticEnergy;
325 G4double q2 = 0.0;
326
327 // scaling to protons for all particles excluding protons and alpha
328 G4int pdg = pdef->GetPDGEncoding();
329 if (pdg != 2212 && pdg != 1000020040)
330 {
331 mass = proton_mass_c2;
332 escaled = kineticEnergy*mass/(pdef->GetPDGMass());
333
334 if(mat) {
335 q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
336 } else {
337 G4double q = pdef->GetPDGCharge()/eplus;
338 q2 = q*q;
339 }
340 }
341
342 if(PIXEshellCS) {
343 xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
344 }
345 if(xsec < 1e-100) {
346 xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
347 }
348
349 if (q2) {xsec *= q2;}
350
351 return xsec;
352}
353
354//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
355
357{
358 minGammaEnergy = cut;
359}
360
361//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
362
364{
365 minElectronEnergy = cut;
366}
367
368//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
369
379
380//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
381
382G4int G4UAtomicDeexcitation::SelectTypeOfTransition(G4int Z, G4int shellId)
383{
384 if (shellId <=0 ) {
385 return 0;
386 }
387
388 G4int provShellId = -1;
389 G4int shellNum = 0;
390 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
391
392 const G4FluoTransition* refShell =
393 transitionManager->ReachableShell(Z,maxNumOfShells-1);
394
395 // This loop gives shellNum the value of the index of shellId
396 // in the vector storing the list of the shells reachable through
397 // a radiative transition
398 if ( shellId <= refShell->FinalShellId())
399 {
400 while (shellId != transitionManager->ReachableShell(Z,shellNum)->FinalShellId())
401 {
402 if(shellNum ==maxNumOfShells-1)
403 {
404 break;
405 }
406 shellNum++;
407 }
408 G4int transProb = 0; //AM change 29/6/07 was 1
409
410 G4double partialProb = G4UniformRand();
411 G4double partSum = 0;
412 const G4FluoTransition* aShell = transitionManager->ReachableShell(Z,shellNum);
413 G4int trSize = (G4int)(aShell->TransitionProbabilities()).size();
414
415 // Loop over the shells wich can provide an electron for a
416 // radiative transition towards shellId:
417 // in every loop the partial sum of the first transProb shells
418 // is calculated and compared with a random number [0,1].
419 // If the partial sum is greater, the shell whose index is transProb
420 // is chosen as the starting shell for a radiative transition
421 // and its identity is returned
422 // Else, terminateded the loop, -1 is returned
423 while(transProb < trSize){
424 partSum += aShell->TransitionProbability(transProb);
425
426 if(partialProb <= partSum)
427 {
428 provShellId = aShell->OriginatingShellId(transProb);
429 break;
430 }
431 ++transProb;
432 }
433 // here provShellId is the right one or is -1.
434 // if -1, the control is passed to the Auger generation part of the package
435 }
436 else
437 {
438 provShellId = -1;
439 }
440 return provShellId;
441}
442
443//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
444
446G4UAtomicDeexcitation::GenerateFluorescence(G4int Z, G4int shellId,
447 G4int provShellId )
448{
449 if (shellId <=0 )
450 {
451 return nullptr;
452 }
453
454 //isotropic angular distribution for the outcoming photon
455 G4double newcosTh = 1.-2.*G4UniformRand();
456 G4double newsinTh = std::sqrt((1.-newcosTh)*(1. + newcosTh));
457 G4double newPhi = twopi*G4UniformRand();
458
459 G4double xDir = newsinTh*std::sin(newPhi);
460 G4double yDir = newsinTh*std::cos(newPhi);
461 G4double zDir = newcosTh;
462
463 G4ThreeVector newGammaDirection(xDir,yDir,zDir);
464
465 G4int shellNum = 0;
466 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
467
468 // find the index of the shell named shellId
469 while (shellId != transitionManager->
470 ReachableShell(Z,shellNum)->FinalShellId())
471 {
472 if(shellNum == maxNumOfShells-1)
473 {
474 break;
475 }
476 ++shellNum;
477 }
478 // number of shell from wich an electron can reach shellId
479 G4int transitionSize = (G4int)transitionManager->
480 ReachableShell(Z,shellNum)->OriginatingShellIds().size();
481
482 G4int index = 0;
483
484 // find the index of the shell named provShellId in the vector
485 // storing the shells from which shellId can be reached
486 while (provShellId != transitionManager->
487 ReachableShell(Z,shellNum)->OriginatingShellId(index))
488 {
489 if(index == transitionSize-1)
490 {
491 break;
492 }
493 ++index;
494 }
495 // energy of the gamma leaving provShellId for shellId
496 G4double transitionEnergy = transitionManager->
497 ReachableShell(Z,shellNum)->TransitionEnergy(index);
498
499 if (transitionEnergy < minGammaEnergy) return nullptr;
500
501 // This is the shell where the new vacancy is: it is the same
502 // shell where the electron came from
503 newShellId = transitionManager->
504 ReachableShell(Z,shellNum)->OriginatingShellId(index);
505
506 G4DynamicParticle* newPart = new G4DynamicParticle(G4Gamma::Gamma(),
507 newGammaDirection,
508 transitionEnergy);
509
510 //Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
511 if (IsAugerCascadeActive()) vacancyArray.push_back(newShellId);
512
513 return newPart;
514}
515
516//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
517
518G4DynamicParticle* G4UAtomicDeexcitation::GenerateAuger(G4int Z, G4int shellId)
519{
520 if(!IsAugerActive()) {
521 // G4cout << "auger inactive!" << G4endl; //debug
522 return nullptr;
523 }
524
525 if (shellId <=0 ) {
526 //G4Exception("G4UAtomicDeexcitation::GenerateAuger()","de0002",
527 // JustWarning, "Energy deposited locally");
528 return nullptr;
529 }
530
531 G4int maxNumOfShells = transitionManager->NumberOfReachableAugerShells(Z);
532
533 const G4AugerTransition* refAugerTransition =
534 transitionManager->ReachableAugerShell(Z,maxNumOfShells-1);
535
536 // This loop gives to shellNum the value of the index of shellId
537 // in the vector storing the list of the vacancies in the variuos shells
538 // that can originate a NON-radiative transition
539 G4int shellNum = 0;
540
541 if ( shellId <= refAugerTransition->FinalShellId() )
542 // "FinalShellId" is final from the point of view of the electron
543 // who makes the transition,
544 // being the Id of the shell in which there is a vacancy
545 {
546 G4int pippo = transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId();
547 if (shellId != pippo ) {
548 do {
549 ++shellNum;
550 if(shellNum == maxNumOfShells)
551 {
552 // G4cout << "No Auger transition found" << G4endl; //debug
553 return 0;
554 }
555 }
556 while (shellId != (transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId()) );
557 }
558
559 // Now we have that shellnum is the shellIndex of the shell named ShellId
560 // G4cout << " the index of the shell is: "<<shellNum<<G4endl;
561 // But we have now to select two shells: one for the transition,
562 // and another for the auger emission.
563 G4int transitionLoopShellIndex = 0;
564 G4double partSum = 0;
565 const G4AugerTransition* anAugerTransition =
566 transitionManager->ReachableAugerShell(Z,shellNum);
567
568 G4int transitionSize = (G4int)
569 (anAugerTransition->TransitionOriginatingShellIds())->size();
570 while (transitionLoopShellIndex < transitionSize) {
571
572 std::vector<G4int>::const_iterator pos =
573 anAugerTransition->TransitionOriginatingShellIds()->cbegin();
574
575 G4int transitionLoopShellId = *(pos+transitionLoopShellIndex);
576 G4int numberOfPossibleAuger = (G4int)
577 (anAugerTransition->AugerTransitionProbabilities(transitionLoopShellId))->size();
578 G4int augerIndex = 0;
579
580 if (augerIndex < numberOfPossibleAuger) {
581 do
582 {
583 G4double thisProb = anAugerTransition->AugerTransitionProbability(augerIndex,
584 transitionLoopShellId);
585 partSum += thisProb;
586 augerIndex++;
587
588 } while (augerIndex < numberOfPossibleAuger);
589 }
590 ++transitionLoopShellIndex;
591 }
592
593 G4double totalVacancyAugerProbability = partSum;
594
595 //And now we start to select the right auger transition and emission
596 G4int transitionRandomShellIndex = 0;
597 G4int transitionRandomShellId = 1;
598 G4int augerIndex = 0;
599 partSum = 0;
600 G4double partialProb = G4UniformRand();
601
602 G4int numberOfPossibleAuger = 0;
603 G4bool foundFlag = false;
604
605 while (transitionRandomShellIndex < transitionSize) {
606
607 std::vector<G4int>::const_iterator pos =
608 anAugerTransition->TransitionOriginatingShellIds()->begin();
609
610 transitionRandomShellId = *(pos+transitionRandomShellIndex);
611
612 augerIndex = 0;
613 numberOfPossibleAuger = (G4int)(anAugerTransition->
614 AugerTransitionProbabilities(transitionRandomShellId))->size();
615
616 while (augerIndex < numberOfPossibleAuger) {
617 G4double thisProb =anAugerTransition->AugerTransitionProbability(augerIndex,
618 transitionRandomShellId);
619
620 partSum += thisProb;
621
622 if (partSum >= (partialProb*totalVacancyAugerProbability) ) { // was /
623 foundFlag = true;
624 break;
625 }
626 augerIndex++;
627 }
628 if (partSum >= (partialProb*totalVacancyAugerProbability) ) {break;} // was /
629 ++transitionRandomShellIndex;
630 }
631
632 // Now we have the index of the shell from wich comes the auger electron (augerIndex),
633 // and the id of the shell, from which the transition e- come (transitionRandomShellid)
634 // If no Transition has been found, 0 is returned.
635 if (!foundFlag) {
636 return nullptr;
637 }
638
639 // Isotropic angular distribution for the outcoming e-
640 G4double newcosTh = 1.-2.*G4UniformRand();
641 G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
642 G4double newPhi = twopi*G4UniformRand();
643
644 G4double xDir = newsinTh*std::sin(newPhi);
645 G4double yDir = newsinTh*std::cos(newPhi);
646 G4double zDir = newcosTh;
647
648 G4ThreeVector newElectronDirection(xDir,yDir,zDir);
649
650 // energy of the auger electron emitted
651 G4double transitionEnergy =
652 anAugerTransition->AugerTransitionEnergy(augerIndex, transitionRandomShellId);
653
654 if (transitionEnergy < minElectronEnergy) {
655 return nullptr;
656 }
657
658 // This is the shell where the new vacancy is: it is the same
659 // shell where the electron came from
660 newShellId = transitionRandomShellId;
661
662 //Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
664 {
665 vacancyArray.push_back(newShellId);
666 vacancyArray.push_back(anAugerTransition->AugerOriginatingShellId(augerIndex,transitionRandomShellId));
667 }
668
669 return new G4DynamicParticle(G4Electron::Electron(),
670 newElectronDirection,
671 transitionEnergy);
672 }
673 else
674 {
675 return nullptr;
676 }
677}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4int ShellId() const
static G4int GetNumberOfShells(G4int Z)
static G4AtomicTransitionManager * Instance()
G4int AugerOriginatingShellId(G4int index, G4int startShellId) const
const G4DataVector * AugerTransitionProbabilities(G4int startShellId) const
G4double AugerTransitionEnergy(G4int index, G4int startShellId) const
const std::vector< G4int > * TransitionOriginatingShellIds() const
Returns the ids of the shells from wich an electron cuuld fill the vacancy in finalShellId.
G4double AugerTransitionProbability(G4int index, G4int startShellId) const
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
const G4String & PIXECrossSectionModel()
const G4String & PIXEElectronCrossSectionModel()
const G4DataVector & TransitionProbabilities() const
Return the probabilities of the transitions.
G4int OriginatingShellId(G4int index) const
Given the index of the originating shells returns its identity.
G4double TransitionProbability(G4int index) const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetCutForSecondaryPhotons(G4double cut)
Set threshold energy for fluorescence.
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut) override
generation of deexcitation for given atom, shell vacancy and cuts
G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr) override
access or compute PIXE cross section
void InitialiseForNewRun() override
initialisation methods
void SetCutForAugerElectrons(G4double cut)
Set threshold energy for Auger electron production.
void InitialiseForExtraAtom(G4int Z) override
const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell) override
G4double ComputeShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr) override
access or compute PIXE cross section
G4VAtomDeexcitation(const G4String &modname="Deexcitation")
G4bool IsAugerCascadeActive() const
#define DBL_MAX
Definition templates.hh:62