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