Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id: G4UAtomicDeexcitation.cc,v 1.11
27// GEANT4 tag $Name: not supported by cvs2svn $
28//
29// -------------------------------------------------------------------
30//
31// Geant4 Class file
32//
33// Authors: Alfonso Mantero (Alfonso.Mantero@ge.infn.it)
34//
35// Created 22 April 2010 from old G4UAtomicDeexcitation class
36//
37// Modified:
38// ---------
39// 20 Oct 2011 Alf modified to take into account ECPSSR form Form Factor
40// 03 Nov 2011 Alf Extended Empirical and Form Factor ionisation XS models
41// out thei ranges with Analytical one.
42// 07 Nov 2011 Alf Restored original ioniation XS for alphas,
43// letting scaled ones for other ions.
44// 20 Mar 2012 LP Register G4PenelopeIonisationCrossSection
45//
46// -------------------------------------------------------------------
47//
48// Class description:
49// Implementation of atomic deexcitation
50//
51// -------------------------------------------------------------------
52
55#include "G4SystemOfUnits.hh"
56#include "Randomize.hh"
57#include "G4Gamma.hh"
59#include "G4FluoTransition.hh"
60#include "G4Electron.hh"
61#include "G4Positron.hh"
62#include "G4Proton.hh"
63#include "G4Alpha.hh"
64
65#include "G4teoCrossSection.hh"
66#include "G4empCrossSection.hh"
69#include "G4EmCorrections.hh"
70#include "G4LossTableManager.hh"
71#include "G4Material.hh"
72#include "G4AtomicShells.hh"
73
74using namespace std;
75
77 G4VAtomDeexcitation("UAtomDeexcitation"),
78 minGammaEnergy(DBL_MAX),
79 minElectronEnergy(DBL_MAX),
80 emcorr(0)
81{
82 PIXEshellCS = 0;
83 ePIXEshellCS = 0;
85 theElectron = G4Electron::Electron();
86 thePositron = G4Positron::Positron();
87 transitionManager = 0;
88 anaPIXEshellCS = 0;
89}
90
92{
93 delete PIXEshellCS;
94 delete anaPIXEshellCS;
95 delete ePIXEshellCS;
96}
97
99{
100 if(!IsFluoActive()) { return; }
101 transitionManager = G4AtomicTransitionManager::Instance();
102 if(IsPIXEActive()) {
103 G4cout << G4endl;
104 G4cout << "### === G4UAtomicDeexcitation::InitialiseForNewRun()" << G4endl;
105 anaPIXEshellCS = new G4teoCrossSection("Analytical");
106
107 }
108 else {return;}
109 // initializing PIXE x-section name
110 //
111 if (PIXECrossSectionModel() == "" ||
112 PIXECrossSectionModel() == "Empirical" ||
113 PIXECrossSectionModel() == "empirical")
114 {
115 SetPIXECrossSectionModel("Empirical");
116 }
117 else if (PIXECrossSectionModel() == "ECPSSR_Analytical" ||
118 PIXECrossSectionModel() == "Analytical" ||
119 PIXECrossSectionModel() == "analytical")
120 {
121 SetPIXECrossSectionModel("Analytical");
122 }
123 else if (PIXECrossSectionModel() == "ECPSSR_FormFactor" ||
124 PIXECrossSectionModel() == "ECPSSR_Tabulated" ||
125 PIXECrossSectionModel() == "Analytical_Tabulated")
126 {
127 SetPIXECrossSectionModel("ECPSSR_FormFactor");
128 }
129 else
130 {
131 G4cout << "### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
132 << G4endl;
133 G4cout << " PIXE cross section name " << PIXECrossSectionModel()
134 << " is unknown, Analytical cross section will be used" << G4endl;
135 SetPIXECrossSectionModel("Analytical");
136 }
137
138 // Check if old model should be deleted
139 if(PIXEshellCS)
140 {
141 if(PIXECrossSectionModel() != PIXEshellCS->GetName())
142 {
143 delete PIXEshellCS;
144 PIXEshellCS = 0;
145 }
146 }
147
148 // Instantiate empirical model
149 if(!PIXEshellCS) {
150 if (PIXECrossSectionModel() == "Empirical")
151 {
152 PIXEshellCS = new G4empCrossSection("Empirical");
153 }
154
155 if (PIXECrossSectionModel() == "ECPSSR_FormFactor")
156 {
157 PIXEshellCS = new G4teoCrossSection("ECPSSR_FormFactor");
158 }
159 }
160
161 // Electron cross section
162 // initializing PIXE x-section name
163 //
164 if (PIXEElectronCrossSectionModel() == "" ||
165 PIXEElectronCrossSectionModel() == "Livermore")
166 {
168 }
169 else if (PIXEElectronCrossSectionModel() == "ProtonAnalytical" ||
170 PIXEElectronCrossSectionModel() == "Analytical" ||
171 PIXEElectronCrossSectionModel() == "analytical")
172 {
173 SetPIXEElectronCrossSectionModel("ProtonAnalytical");
174 }
175 else if (PIXEElectronCrossSectionModel() == "ProtonEmpirical" ||
176 PIXEElectronCrossSectionModel() == "Empirical" ||
177 PIXEElectronCrossSectionModel() == "empirical")
178 {
179 SetPIXEElectronCrossSectionModel("ProtonEmpirical");
180 }
181 else if (PIXEElectronCrossSectionModel() == "Penelope")
183 else
184 {
185 G4cout << "### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
186 << G4endl;
187 G4cout << " PIXE e- cross section name " << PIXEElectronCrossSectionModel()
188 << " is unknown, PIXE is disabled" << G4endl;
190 }
191
192 // Check if old model should be deleted
193 if(ePIXEshellCS)
194 {
195 if(PIXEElectronCrossSectionModel() != ePIXEshellCS->GetName())
196 {
197 delete ePIXEshellCS;
198 ePIXEshellCS = 0;
199 }
200 }
201
202 // Instantiate empirical model
203 if(!ePIXEshellCS)
204 {
205 if(PIXEElectronCrossSectionModel() == "Empirical")
206 {
207 ePIXEshellCS = new G4empCrossSection("Empirical");
208 }
209
210 else if(PIXEElectronCrossSectionModel() == "Analytical")
211 {
212 ePIXEshellCS = new G4teoCrossSection("Analytical");
213 }
214
215 else if(PIXEElectronCrossSectionModel() == "Livermore")
216 {
217 ePIXEshellCS = new G4LivermoreIonisationCrossSection();
218 }
219 else if (PIXEElectronCrossSectionModel() == "Penelope")
220 {
221 ePIXEshellCS = new G4PenelopeIonisationCrossSection();
222 }
223 }
224}
225
227{}
228
230{
231 return transitionManager->Shell(Z, size_t(shell));
232}
233
235 std::vector<G4DynamicParticle*>* vectorOfParticles,
236 const G4AtomicShell* atomicShell,
237 G4int Z,
238 G4double gammaCut,
239 G4double eCut)
240{
241
242 // Defined initial conditions
243 G4int givenShellId = atomicShell->ShellId();
244 // G4cout << "generating particles for vacancy in shellId: " << givenShellId << G4endl;
245 minGammaEnergy = gammaCut;
246 minElectronEnergy = eCut;
247
248 // V.I. short-cut
249 if(!IsAugerActive()) { minElectronEnergy = DBL_MAX; }
250
251 // generation secondaries
252 G4DynamicParticle* aParticle=0;
253 G4int provShellId = 0;
254 G4int counter = 0;
255
256 // let's check that 5<Z<100
257
258 if (Z>5 && Z<100) {
259
260 // The aim of this loop is to generate more than one fluorecence photon
261 // from the same ionizing event
262 do
263 {
264 if (counter == 0)
265 // First call to GenerateParticles(...):
266 // givenShellId is given by the process
267 {
268 provShellId = SelectTypeOfTransition(Z, givenShellId);
269
270 if ( provShellId >0)
271 {
272 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
273 //if (aParticle != 0) { G4cout << "****FLUO!_1**** " << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl ;} //debug
274 }
275 else if ( provShellId == -1)
276 {
277 aParticle = GenerateAuger(Z, givenShellId);
278 //if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
279 }
280 else
281 {
282 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
283 }
284 }
285 else
286 // Following calls to GenerateParticles(...):
287 // newShellId is given by GenerateFluorescence(...)
288 {
289 provShellId = SelectTypeOfTransition(Z,newShellId);
290 if (provShellId >0)
291 {
292 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
293 //if (aParticle != 0) { G4cout << "****FLUO!_2****" << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl;} //debug
294 }
295 else if ( provShellId == -1)
296 {
297 aParticle = GenerateAuger(Z, newShellId);
298 //if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
299 }
300 else
301 {
302 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
303 }
304 }
305 counter++;
306 if (aParticle != 0)
307 {
308 vectorOfParticles->push_back(aParticle);
309 //G4cout << "Deexcitation Occurred!" << G4endl; //debug
310 }
311 else {provShellId = -2;}
312 }
313 while (provShellId > -2);
314 }
315 else
316 {
317 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
318 }
319
320 // G4cout << "---------FATTO!---------" << G4endl; //debug
321
322}
323
326 const G4ParticleDefinition* pdef,
327 G4int Z,
328 G4AtomicShellEnumerator shellEnum,
329 G4double kineticEnergy,
330 const G4Material* mat)
331{
332
333 // we must put a control on the shell that are passed:
334 // some shells should not pass (line "0" or "2")
335
336 // check atomic number
337 G4double xsec = 0.0;
338 if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
339 G4int idx = G4int(shellEnum);
340 if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
341
342 //
343 if(pdef == theElectron || pdef == thePositron) {
344 xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
345 return xsec;
346 }
347
348 G4double mass = pdef->GetPDGMass();
349 G4double escaled = kineticEnergy;
350 G4double q2 = 0.0;
351
352 // scaling to protons
353 if ((pdef->GetParticleName() != "proton" && pdef->GetParticleName() != "alpha" ) )
354 {
355 mass = proton_mass_c2;
356 escaled = kineticEnergy*mass/(pdef->GetPDGMass());
357
358 if(mat) {
359 q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
360 } else {
361 G4double q = pdef->GetPDGCharge()/eplus;
362 q2 = q*q;
363 }
364 }
365
366 if(PIXEshellCS) { xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat); }
367 if(xsec < 1e-100) {
368
369 xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
370
371 }
372
373 if (q2) {xsec *= q2;}
374
375 return xsec;
376}
377
379{
380 minGammaEnergy = cut;
381}
382
384{
385 minElectronEnergy = cut;
386}
387
390 const G4ParticleDefinition* p,
391 G4int Z,
393 G4double kinE,
394 const G4Material* mat)
395{
396 return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
397}
398
399G4int G4UAtomicDeexcitation::SelectTypeOfTransition(G4int Z, G4int shellId)
400{
401 if (shellId <=0 ) {
402 G4Exception("G4UAtomicDeexcitation::SelecttypeOfTransition()","de0002",JustWarning, "Energy deposited locally");
403 return 0;
404 }
405 //G4bool fluoTransitionFoundFlag = false;
406
407 G4int provShellId = -1;
408 G4int shellNum = 0;
409 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
410
411 const G4FluoTransition* refShell = 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
443 partSum += aShell->TransitionProbability(transProb);
444
445 if(partialProb <= partSum)
446 {
447 provShellId = aShell->OriginatingShellId(transProb);
448 //fluoTransitionFoundFlag = true;
449
450 break;
451 }
452 transProb++;
453 }
454
455 // here provShellId is the right one or is -1.
456 // if -1, the control is passed to the Auger generation part of the package
457 }
458 else
459 {
460 provShellId = -1;
461 }
462 //G4cout << "FlagTransition= " << provShellId << " ecut(MeV)= " << minElectronEnergy
463 // << " gcut(MeV)= " << minGammaEnergy << G4endl;
464 return provShellId;
465}
466
468G4UAtomicDeexcitation::GenerateFluorescence(G4int Z, G4int shellId,
469 G4int provShellId )
470{
471
472 //if(!IsDeexActive()) { return 0; }
473
474 if (shellId <=0 )
475
476 {
477 G4Exception("G4UAtomicDeexcitation::GenerateFluorescence()","de0002",JustWarning, "Energy deposited locally");
478 return 0;
479 }
480
481
482 //isotropic angular distribution for the outcoming photon
483 G4double newcosTh = 1.-2.*G4UniformRand();
484 G4double newsinTh = std::sqrt((1.-newcosTh)*(1. + newcosTh));
485 G4double newPhi = twopi*G4UniformRand();
486
487 G4double xDir = newsinTh*std::sin(newPhi);
488 G4double yDir = newsinTh*std::cos(newPhi);
489 G4double zDir = newcosTh;
490
491 G4ThreeVector newGammaDirection(xDir,yDir,zDir);
492
493 G4int shellNum = 0;
494 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
495
496 // find the index of the shell named shellId
497 while (shellId != transitionManager->
498 ReachableShell(Z,shellNum)->FinalShellId())
499 {
500 if(shellNum == maxNumOfShells-1)
501 {
502 break;
503 }
504 shellNum++;
505 }
506 // number of shell from wich an electron can reach shellId
507 size_t transitionSize = transitionManager->
508 ReachableShell(Z,shellNum)->OriginatingShellIds().size();
509
510 size_t index = 0;
511
512 // find the index of the shell named provShellId in the vector
513 // storing the shells from which shellId can be reached
514 while (provShellId != transitionManager->
515 ReachableShell(Z,shellNum)->OriginatingShellId(index))
516 {
517 if(index == transitionSize-1)
518 {
519 break;
520 }
521 index++;
522 }
523 // energy of the gamma leaving provShellId for shellId
524 G4double transitionEnergy = transitionManager->
525 ReachableShell(Z,shellNum)->TransitionEnergy(index);
526
527 if (transitionEnergy < minGammaEnergy) return 0;
528
529 // This is the shell where the new vacancy is: it is the same
530 // shell where the electron came from
531 newShellId = transitionManager->
532 ReachableShell(Z,shellNum)->OriginatingShellId(index);
533
534
536 newGammaDirection,
537 transitionEnergy);
538 return newPart;
539}
540
541G4DynamicParticle* G4UAtomicDeexcitation::GenerateAuger(G4int Z, G4int shellId)
542{
543 if(!IsAugerActive()) { return 0; }
544
545 if (shellId <=0 ) {
546 {
547
548 G4Exception("G4UAtomicDeexcitation::GenerateAuger()","de0002",JustWarning, "Energy deposited locally");
549
550 return 0;
551 }
552 }
553 // G4int provShellId = -1;
554 G4int maxNumOfShells = transitionManager->NumberOfReachableAugerShells(Z);
555
556 const G4AugerTransition* refAugerTransition =
557 transitionManager->ReachableAugerShell(Z,maxNumOfShells-1);
558
559 // This loop gives to shellNum the value of the index of shellId
560 // in the vector storing the list of the vacancies in the variuos shells
561 // that can originate a NON-radiative transition
562
563 G4int shellNum = 0;
564
565 if ( shellId <= refAugerTransition->FinalShellId() )
566 //"FinalShellId" is final from the point of view of the elctron who makes the transition,
567 // being the Id of the shell in which there is a vacancy
568 {
569 G4int pippo = transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId();
570 if (shellId != pippo ) {
571 do {
572 shellNum++;
573 if(shellNum == maxNumOfShells)
574 {
575 //G4Exception("G4UAtomicDeexcitation: No Auger transition found");
576 return 0;
577 }
578 }
579 while (shellId != (transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId()) ) ;
580 }
581
582
583 // Now we have that shellnum is the shellIndex of the shell named ShellId
584
585 // G4cout << " the index of the shell is: "<<shellNum<<G4endl;
586
587 // But we have now to select two shells: one for the transition,
588 // and another for the auger emission.
589
590 G4int transitionLoopShellIndex = 0;
591 G4double partSum = 0;
592 const G4AugerTransition* anAugerTransition =
593 transitionManager->ReachableAugerShell(Z,shellNum);
594
595 // G4cout << " corresponding to the ID: "<< anAugerTransition->FinalShellId() << G4endl;
596
597
598 G4int transitionSize =
599 (anAugerTransition->TransitionOriginatingShellIds())->size();
600 while (transitionLoopShellIndex < transitionSize) {
601
602 std::vector<G4int>::const_iterator pos =
603 anAugerTransition->TransitionOriginatingShellIds()->begin();
604
605 G4int transitionLoopShellId = *(pos+transitionLoopShellIndex);
606 G4int numberOfPossibleAuger =
607 (anAugerTransition->AugerTransitionProbabilities(transitionLoopShellId))->size();
608 G4int augerIndex = 0;
609 // G4int partSum2 = 0;
610
611
612 if (augerIndex < numberOfPossibleAuger) {
613
614 do
615 {
616 G4double thisProb = anAugerTransition->AugerTransitionProbability(augerIndex,
617 transitionLoopShellId);
618 partSum += thisProb;
619 augerIndex++;
620
621 } while (augerIndex < numberOfPossibleAuger);
622 }
623 transitionLoopShellIndex++;
624 }
625
626
627
628 // Now we have the entire probability of an auger transition for the vacancy
629 // located in shellNum (index of shellId)
630
631 // AM *********************** F I X E D **************************** AM
632 // Here we duplicate the previous loop, this time looking to the sum of the probabilities
633 // to be under the random number shoot by G4 UniformRdandom. This could have been done in the
634 // previuos loop, while integrating the probabilities. There is a bug that will be fixed
635 // 5 minutes from now: a line:
636 // G4int numberOfPossibleAuger = (anAugerTransition->
637 // AugerTransitionProbabilities(transitionLoopShellId))->size();
638 // to be inserted.
639 // AM *********************** F I X E D **************************** AM
640
641 // Remains to get the same result with a single loop.
642
643 // AM *********************** F I X E D **************************** AM
644 // Another Bug: in EADL Auger Transition are normalized to all the transitions deriving from
645 // a vacancy in one shell, but not all of these are present in data tables. So if a transition
646 // doesn't occur in the main one a local energy deposition must occur, instead of (like now)
647 // generating the last transition present in EADL data.
648 // AM *********************** F I X E D **************************** AM
649
650
651 G4double totalVacancyAugerProbability = partSum;
652
653
654 //And now we start to select the right auger transition and emission
655 G4int transitionRandomShellIndex = 0;
656 G4int transitionRandomShellId = 1;
657 G4int augerIndex = 0;
658 partSum = 0;
659 G4double partialProb = G4UniformRand();
660 // G4int augerOriginatingShellId = 0;
661
662 G4int numberOfPossibleAuger = 0;
663
664 G4bool foundFlag = false;
665
666 while (transitionRandomShellIndex < transitionSize) {
667
668 std::vector<G4int>::const_iterator pos =
669 anAugerTransition->TransitionOriginatingShellIds()->begin();
670
671 transitionRandomShellId = *(pos+transitionRandomShellIndex);
672
673 augerIndex = 0;
674 numberOfPossibleAuger = (anAugerTransition->
675 AugerTransitionProbabilities(transitionRandomShellId))->size();
676
677 while (augerIndex < numberOfPossibleAuger) {
678 G4double thisProb =anAugerTransition->AugerTransitionProbability(augerIndex,
679 transitionRandomShellId);
680
681 partSum += thisProb;
682
683 if (partSum >= (partialProb*totalVacancyAugerProbability) ) { // was /
684 foundFlag = true;
685 break;
686 }
687 augerIndex++;
688 }
689 if (partSum >= (partialProb*totalVacancyAugerProbability) ) {break;} // was /
690 transitionRandomShellIndex++;
691 }
692
693 // Now we have the index of the shell from wich comes the auger electron (augerIndex),
694 // and the id of the shell, from which the transition e- come (transitionRandomShellid)
695 // If no Transition has been found, 0 is returned.
696
697 if (!foundFlag) {return 0;}
698
699 // Isotropic angular distribution for the outcoming e-
700 G4double newcosTh = 1.-2.*G4UniformRand();
701 G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
702 G4double newPhi = twopi*G4UniformRand();
703
704 G4double xDir = newsinTh*std::sin(newPhi);
705 G4double yDir = newsinTh*std::cos(newPhi);
706 G4double zDir = newcosTh;
707
708 G4ThreeVector newElectronDirection(xDir,yDir,zDir);
709
710 // energy of the auger electron emitted
711
712
713 G4double transitionEnergy = anAugerTransition->AugerTransitionEnergy(augerIndex, transitionRandomShellId);
714 /*
715 G4cout << "AUger TransitionId " << anAugerTransition->FinalShellId() << G4endl;
716 G4cout << "augerIndex: " << augerIndex << G4endl;
717 G4cout << "transitionShellId: " << transitionRandomShellId << G4endl;
718 */
719
720 if (transitionEnergy < minElectronEnergy) {
721 //G4cout << "Problem!" << G4endl; // debug
722 return 0;
723 }
724
725 // This is the shell where the new vacancy is: it is the same
726 // shell where the electron came from
727 newShellId = transitionRandomShellId;
728
730 newElectronDirection,
731 transitionEnergy);
732 }
733 else
734 {
735 //G4Exception("G4UAtomicDeexcitation: no auger transition found");
736 return 0;
737 }
738}
G4AtomicShellEnumerator
@ JustWarning
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
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
static G4AtomicTransitionManager * Instance()
G4int NumberOfReachableAugerShells(G4int Z) const
G4int FinalShellId() const
const G4DataVector * AugerTransitionProbabilities(G4int startShellId) const
G4double AugerTransitionEnergy(G4int index, G4int startShellId) const
const std::vector< G4int > * TransitionOriginatingShellIds() const
G4double AugerTransitionProbability(G4int index, G4int startShellId) const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4DataVector & TransitionProbabilities() const
G4int OriginatingShellId(G4int index) const
G4double TransitionProbability(G4int index) const
G4int FinalShellId() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Positron * Positron()
Definition: G4Positron.cc:94
void SetCutForSecondaryPhotons(G4double cut)
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)
virtual void InitialiseForExtraAtom(G4int Z)
virtual G4double ComputeShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)
virtual void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut)
virtual void InitialiseForNewRun()
void SetCutForAugerElectrons(G4double cut)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)
G4bool IsAugerActive() const
const G4String & PIXECrossSectionModel() const
void SetPIXECrossSectionModel(const G4String &)
const G4String & PIXEElectronCrossSectionModel() const
void SetPIXEElectronCrossSectionModel(const G4String &)
virtual G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)=0
const G4String & GetName() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define DBL_MAX
Definition: templates.hh:83