Geant4 11.1.1
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
174{}
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
198 // generation secondaries
199 G4DynamicParticle* aParticle=0;
200 G4int provShellId = 0;
201
202 //ORIGINAL METHOD BY ALFONSO MANTERO
204 {
205 //----------------------------
206 G4int counter = 0;
207
208 // limits of the EPDL data
209 if (Z>5 && Z<105) {
210
211 // The aim of this loop is to generate more than one fluorecence photon
212 // from the same ionizing event
213 do
214 {
215 if (counter == 0)
216 // First call to GenerateParticles(...):
217 // givenShellId is given by the process
218 {
219 provShellId = SelectTypeOfTransition(Z, givenShellId);
220
221 if (provShellId >0)
222 {
223 aParticle =
224 GenerateFluorescence(Z, givenShellId, provShellId);
225 }
226 else if (provShellId == -1)
227 {
228 aParticle = GenerateAuger(Z, givenShellId);
229 }
230 }
231 else
232 // Following calls to GenerateParticles(...):
233 // newShellId is given by GenerateFluorescence(...)
234 {
235 provShellId = SelectTypeOfTransition(Z,newShellId);
236 if (provShellId >0)
237 {
238 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
239 }
240 else if ( provShellId == -1)
241 {
242 aParticle = GenerateAuger(Z, newShellId);
243 }
244 }
245 ++counter;
246 if (aParticle != 0)
247 {
248 vectorOfParticles->push_back(aParticle);
249 }
250 else {provShellId = -2;}
251 }
252 while (provShellId > -2);
253 }
254 } // Auger cascade is not active
255
256 //END OF ORIGINAL METHOD BY ALFONSO MANTERO
257 //----------------------
258
259 // NEW METHOD
260 // Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
262 {
263 //----------------------
264 vacancyArray.push_back(givenShellId);
265
266 // let's check that 5<Z<100
267 if (Z<6 || Z>104){
268 return;
269 }
270
271 // as long as there is vacancy to be filled by either fluo or auger, stay in the loop.
272 while(!vacancyArray.empty()){
273 // prepare to process the last element, and then delete it from the vector.
274 givenShellId = vacancyArray[0];
275 provShellId = SelectTypeOfTransition(Z,givenShellId);
276
277 //G4cout<<"\n------ Atom Transition with Z: "<<Z<<"\tbetween current:"
278 // <<givenShellId<<" & target:"<<provShellId<<G4endl;
279 if(provShellId>0){
280 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
281 }
282 else if(provShellId == -1){
283 aParticle = GenerateAuger(Z, givenShellId);
284 }
285 // if a particle is created, put it in the vector of new particles
286 if(aParticle!=0)
287 vectorOfParticles->push_back(aParticle);
288
289 // one vacancy has been processed. Erase it.
290 vacancyArray.erase(vacancyArray.begin());
291 }
292 //----------------------
293 //End of Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
294
295 } // Auger cascade is active
296}
297
298//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
299
302 const G4ParticleDefinition* pdef,
303 G4int Z,
304 G4AtomicShellEnumerator shellEnum,
305 G4double kineticEnergy,
306 const G4Material* mat)
307{
308 // we must put a control on the shell that are passed:
309 // some shells should not pass (line "0" or "2")
310
311 // check atomic number
312 G4double xsec = 0.0;
313 if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
314 G4int idx = G4int(shellEnum);
315 if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
316
317 if(pdef == theElectron || pdef == thePositron) {
318 xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
319 return xsec;
320 }
321
322 G4double mass = pdef->GetPDGMass();
323 G4double escaled = kineticEnergy;
324 G4double q2 = 0.0;
325
326 // scaling to protons for all particles excluding protons and alpha
327 G4int pdg = pdef->GetPDGEncoding();
328 if (pdg != 2212 && pdg != 1000020040)
329 {
330 mass = proton_mass_c2;
331 escaled = kineticEnergy*mass/(pdef->GetPDGMass());
332
333 if(mat) {
334 q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
335 } else {
336 G4double q = pdef->GetPDGCharge()/eplus;
337 q2 = q*q;
338 }
339 }
340
341 if(PIXEshellCS) {
342 xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
343 }
344 if(xsec < 1e-100) {
345 xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
346 }
347
348 if (q2) {xsec *= q2;}
349
350 return xsec;
351}
352
353//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
354
356{
357 minGammaEnergy = cut;
358}
359
360//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
361
363{
364 minElectronEnergy = cut;
365}
366
367//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
368
370 const G4ParticleDefinition* p,
371 G4int Z,
373 G4double kinE,
374 const G4Material* mat)
375{
376 return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
377}
378
379//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
380
381G4int G4UAtomicDeexcitation::SelectTypeOfTransition(G4int Z, G4int shellId)
382{
383 if (shellId <=0 ) {
384 return 0;
385 }
386
387 G4int provShellId = -1;
388 G4int shellNum = 0;
389 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
390
391 const G4FluoTransition* refShell =
392 transitionManager->ReachableShell(Z,maxNumOfShells-1);
393
394 // This loop gives shellNum the value of the index of shellId
395 // in the vector storing the list of the shells reachable through
396 // a radiative transition
397 if ( shellId <= refShell->FinalShellId())
398 {
399 while (shellId != transitionManager->ReachableShell(Z,shellNum)->FinalShellId())
400 {
401 if(shellNum ==maxNumOfShells-1)
402 {
403 break;
404 }
405 shellNum++;
406 }
407 G4int transProb = 0; //AM change 29/6/07 was 1
408
409 G4double partialProb = G4UniformRand();
410 G4double partSum = 0;
411 const G4FluoTransition* aShell = transitionManager->ReachableShell(Z,shellNum);
412 G4int trSize = (G4int)(aShell->TransitionProbabilities()).size();
413
414 // Loop over the shells wich can provide an electron for a
415 // radiative transition towards shellId:
416 // in every loop the partial sum of the first transProb shells
417 // is calculated and compared with a random number [0,1].
418 // If the partial sum is greater, the shell whose index is transProb
419 // is chosen as the starting shell for a radiative transition
420 // and its identity is returned
421 // Else, terminateded the loop, -1 is returned
422 while(transProb < trSize){
423 partSum += aShell->TransitionProbability(transProb);
424
425 if(partialProb <= partSum)
426 {
427 provShellId = aShell->OriginatingShellId(transProb);
428 break;
429 }
430 ++transProb;
431 }
432 // here provShellId is the right one or is -1.
433 // if -1, the control is passed to the Auger generation part of the package
434 }
435 else
436 {
437 provShellId = -1;
438 }
439 return provShellId;
440}
441
442//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
443
445G4UAtomicDeexcitation::GenerateFluorescence(G4int Z, G4int shellId,
446 G4int provShellId )
447{
448 if (shellId <=0 )
449 {
450 return nullptr;
451 }
452
453 //isotropic angular distribution for the outcoming photon
454 G4double newcosTh = 1.-2.*G4UniformRand();
455 G4double newsinTh = std::sqrt((1.-newcosTh)*(1. + newcosTh));
456 G4double newPhi = twopi*G4UniformRand();
457
458 G4double xDir = newsinTh*std::sin(newPhi);
459 G4double yDir = newsinTh*std::cos(newPhi);
460 G4double zDir = newcosTh;
461
462 G4ThreeVector newGammaDirection(xDir,yDir,zDir);
463
464 G4int shellNum = 0;
465 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
466
467 // find the index of the shell named shellId
468 while (shellId != transitionManager->
469 ReachableShell(Z,shellNum)->FinalShellId())
470 {
471 if(shellNum == maxNumOfShells-1)
472 {
473 break;
474 }
475 ++shellNum;
476 }
477 // number of shell from wich an electron can reach shellId
478 G4int transitionSize = (G4int)transitionManager->
479 ReachableShell(Z,shellNum)->OriginatingShellIds().size();
480
481 G4int index = 0;
482
483 // find the index of the shell named provShellId in the vector
484 // storing the shells from which shellId can be reached
485 while (provShellId != transitionManager->
486 ReachableShell(Z,shellNum)->OriginatingShellId(index))
487 {
488 if(index == transitionSize-1)
489 {
490 break;
491 }
492 ++index;
493 }
494 // energy of the gamma leaving provShellId for shellId
495 G4double transitionEnergy = transitionManager->
496 ReachableShell(Z,shellNum)->TransitionEnergy(index);
497
498 if (transitionEnergy < minGammaEnergy) return nullptr;
499
500 // This is the shell where the new vacancy is: it is the same
501 // shell where the electron came from
502 newShellId = transitionManager->
503 ReachableShell(Z,shellNum)->OriginatingShellId(index);
504
506 newGammaDirection,
507 transitionEnergy);
508
509 //Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
510 if (IsAugerCascadeActive()) vacancyArray.push_back(newShellId);
511
512 return newPart;
513}
514
515//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
516
517G4DynamicParticle* G4UAtomicDeexcitation::GenerateAuger(G4int Z, G4int shellId)
518{
519 if(!IsAugerActive()) {
520 // G4cout << "auger inactive!" << G4endl; //debug
521 return nullptr;
522 }
523
524 if (shellId <=0 ) {
525 //G4Exception("G4UAtomicDeexcitation::GenerateAuger()","de0002",
526 // JustWarning, "Energy deposited locally");
527 return nullptr;
528 }
529
530 G4int maxNumOfShells = transitionManager->NumberOfReachableAugerShells(Z);
531
532 const G4AugerTransition* refAugerTransition =
533 transitionManager->ReachableAugerShell(Z,maxNumOfShells-1);
534
535 // This loop gives to shellNum the value of the index of shellId
536 // in the vector storing the list of the vacancies in the variuos shells
537 // that can originate a NON-radiative transition
538 G4int shellNum = 0;
539
540 if ( shellId <= refAugerTransition->FinalShellId() )
541 // "FinalShellId" is final from the point of view of the electron
542 // who makes the transition,
543 // being the Id of the shell in which there is a vacancy
544 {
545 G4int pippo = transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId();
546 if (shellId != pippo ) {
547 do {
548 ++shellNum;
549 if(shellNum == maxNumOfShells)
550 {
551 // G4cout << "No Auger transition found" << G4endl; //debug
552 return 0;
553 }
554 }
555 while (shellId != (transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId()) );
556 }
557
558 // Now we have that shellnum is the shellIndex of the shell named ShellId
559 // G4cout << " the index of the shell is: "<<shellNum<<G4endl;
560 // But we have now to select two shells: one for the transition,
561 // and another for the auger emission.
562 G4int transitionLoopShellIndex = 0;
563 G4double partSum = 0;
564 const G4AugerTransition* anAugerTransition =
565 transitionManager->ReachableAugerShell(Z,shellNum);
566
567 G4int transitionSize = (G4int)
568 (anAugerTransition->TransitionOriginatingShellIds())->size();
569 while (transitionLoopShellIndex < transitionSize) {
570
571 std::vector<G4int>::const_iterator pos =
572 anAugerTransition->TransitionOriginatingShellIds()->cbegin();
573
574 G4int transitionLoopShellId = *(pos+transitionLoopShellIndex);
575 G4int numberOfPossibleAuger = (G4int)
576 (anAugerTransition->AugerTransitionProbabilities(transitionLoopShellId))->size();
577 G4int augerIndex = 0;
578
579 if (augerIndex < numberOfPossibleAuger) {
580 do
581 {
582 G4double thisProb = anAugerTransition->AugerTransitionProbability(augerIndex,
583 transitionLoopShellId);
584 partSum += thisProb;
585 augerIndex++;
586
587 } while (augerIndex < numberOfPossibleAuger);
588 }
589 ++transitionLoopShellIndex;
590 }
591
592 G4double totalVacancyAugerProbability = partSum;
593
594 //And now we start to select the right auger transition and emission
595 G4int transitionRandomShellIndex = 0;
596 G4int transitionRandomShellId = 1;
597 G4int augerIndex = 0;
598 partSum = 0;
599 G4double partialProb = G4UniformRand();
600
601 G4int numberOfPossibleAuger = 0;
602 G4bool foundFlag = false;
603
604 while (transitionRandomShellIndex < transitionSize) {
605
606 std::vector<G4int>::const_iterator pos =
607 anAugerTransition->TransitionOriginatingShellIds()->begin();
608
609 transitionRandomShellId = *(pos+transitionRandomShellIndex);
610
611 augerIndex = 0;
612 numberOfPossibleAuger = (G4int)(anAugerTransition->
613 AugerTransitionProbabilities(transitionRandomShellId))->size();
614
615 while (augerIndex < numberOfPossibleAuger) {
616 G4double thisProb =anAugerTransition->AugerTransitionProbability(augerIndex,
617 transitionRandomShellId);
618
619 partSum += thisProb;
620
621 if (partSum >= (partialProb*totalVacancyAugerProbability) ) { // was /
622 foundFlag = true;
623 break;
624 }
625 augerIndex++;
626 }
627 if (partSum >= (partialProb*totalVacancyAugerProbability) ) {break;} // was /
628 ++transitionRandomShellIndex;
629 }
630
631 // Now we have the index of the shell from wich comes the auger electron (augerIndex),
632 // and the id of the shell, from which the transition e- come (transitionRandomShellid)
633 // If no Transition has been found, 0 is returned.
634 if (!foundFlag) {
635 return nullptr;
636 }
637
638 // Isotropic angular distribution for the outcoming e-
639 G4double newcosTh = 1.-2.*G4UniformRand();
640 G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
641 G4double newPhi = twopi*G4UniformRand();
642
643 G4double xDir = newsinTh*std::sin(newPhi);
644 G4double yDir = newsinTh*std::cos(newPhi);
645 G4double zDir = newcosTh;
646
647 G4ThreeVector newElectronDirection(xDir,yDir,zDir);
648
649 // energy of the auger electron emitted
650 G4double transitionEnergy =
651 anAugerTransition->AugerTransitionEnergy(augerIndex, transitionRandomShellId);
652
653 if (transitionEnergy < minElectronEnergy) {
654 return nullptr;
655 }
656
657 // This is the shell where the new vacancy is: it is the same
658 // shell where the electron came from
659 newShellId = transitionRandomShellId;
660
661 //Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
663 {
664 vacancyArray.push_back(newShellId);
665 vacancyArray.push_back(anAugerTransition->AugerOriginatingShellId(augerIndex,transitionRandomShellId));
666 }
667
669 newElectronDirection,
670 transitionEnergy);
671 }
672 else
673 {
674 return nullptr;
675 }
676}
G4AtomicShellEnumerator
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4int ShellId() const
static G4int GetNumberOfShells(G4int Z)
G4int NumberOfReachableShells(G4int Z) const
const G4AugerTransition * ReachableAugerShell(G4int Z, G4int shellIndex) const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
const G4FluoTransition * ReachableShell(G4int Z, size_t shellIndex) const
void Initialise()
needs to be called once from other code before start of run
static G4AtomicTransitionManager * Instance()
G4int NumberOfReachableAugerShells(G4int Z) const
G4int AugerOriginatingShellId(G4int index, G4int startShellId) const
G4int FinalShellId() const
returns the id of the shell in wich the transition electron arrives
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:93
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
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
G4int FinalShellId() const
Return the identity if the vacancy.
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
G4double GetPDGCharge() const
static G4Positron * Positron()
Definition: G4Positron.cc:93
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
G4bool IsAugerActive() const
G4bool IsAugerCascadeActive() const
virtual G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)=0
const G4String & GetName() const
#define DBL_MAX
Definition: templates.hh:62