Geant4 10.7.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
76 G4VAtomDeexcitation("UAtomDeexcitation"),
77 minGammaEnergy(DBL_MAX),
78 minElectronEnergy(DBL_MAX)
79{
80 anaPIXEshellCS = nullptr;
81 PIXEshellCS = nullptr;
82 ePIXEshellCS = nullptr;
84 theElectron = G4Electron::Electron();
85 thePositron = G4Positron::Positron();
86 transitionManager = G4AtomicTransitionManager::Instance();
87}
88
90{
91 delete anaPIXEshellCS;
92 delete PIXEshellCS;
93 delete ePIXEshellCS;
94}
95
97{
98 if(!IsFluoActive()) { return; }
99
100 transitionManager->Initialise();
101
102 if(!IsPIXEActive()) { return; }
103
104 if(!anaPIXEshellCS) {
105 anaPIXEshellCS = new G4teoCrossSection("ECPSSR_Analytical");
106 }
107 G4cout << G4endl;
108 G4cout << "### === G4UAtomicDeexcitation::InitialiseForNewRun()" << G4endl;
109
111 G4String namePIXExsModel = param->PIXECrossSectionModel();
112 G4String namePIXExsElectronModel = param->PIXEElectronCrossSectionModel();
113
114 //G4cout << namePIXExsModel << " " << namePIXExsElectronModel << G4endl;
115
116 // Check if old cross section for p/ion should be deleted
117 if(PIXEshellCS && namePIXExsModel != PIXEshellCS->GetName())
118 {
119 delete PIXEshellCS;
120 PIXEshellCS = nullptr;
121 }
122
123 // Instantiate new proton/ion cross section
124 if(!PIXEshellCS) {
125 if (namePIXExsModel == "ECPSSR_FormFactor")
126 {
127 PIXEshellCS = new G4teoCrossSection(namePIXExsModel);
128 }
129 else if(namePIXExsModel == "Empirical")
130 {
131 PIXEshellCS = new G4empCrossSection(namePIXExsModel);
132 }
133 }
134 //G4cout << "PIXE is initialised" << G4endl;
135
136 // Check if old cross section for e+- should be deleted
137 if(ePIXEshellCS && namePIXExsElectronModel != ePIXEshellCS->GetName())
138 {
139 delete ePIXEshellCS;
140 ePIXEshellCS = nullptr;
141 }
142
143 // Instantiate new e+- cross section
144 if(!ePIXEshellCS)
145 {
146 if(namePIXExsElectronModel == "Empirical")
147 {
148 ePIXEshellCS = new G4empCrossSection("Empirical");
149 }
150 else if(namePIXExsElectronModel == "ECPSSR_Analytical")
151 {
152 ePIXEshellCS = new G4teoCrossSection("ECPSSR_Analytical");
153 }
154 else if (namePIXExsElectronModel == "Penelope")
155 {
156 ePIXEshellCS = new G4PenelopeIonisationCrossSection();
157 }
158 else
159 {
160 ePIXEshellCS = new G4LivermoreIonisationCrossSection();
161 }
162 }
163 //G4cout << "ePIXE is initialised" << G4endl;
164}
165
167{}
168
169const G4AtomicShell*
171{
172 return transitionManager->Shell(Z, size_t(shell));
173}
174
176 std::vector<G4DynamicParticle*>* vectorOfParticles,
177 const G4AtomicShell* atomicShell,
178 G4int Z,
179 G4double gammaCut,
180 G4double eCut)
181{
182
183 // Defined initial conditions
184 G4int givenShellId = atomicShell->ShellId();
185 //G4cout << "generating particles for vacancy in shellId: "
186 // << givenShellId << G4endl; // debug
187 minGammaEnergy = gammaCut;
188 minElectronEnergy = eCut;
189
190 // generation secondaries
191 G4DynamicParticle* aParticle=0;
192 G4int provShellId = 0;
193
194//ORIGINAL METHOD BY ALFONSO MANTERO
196{
197//----------------------------
198 G4int counter = 0;
199
200 // let's check that 5<Z<100
201
202 if (Z>5 && Z<100) {
203
204 // The aim of this loop is to generate more than one fluorecence photon
205 // from the same ionizing event
206 do
207 {
208 if (counter == 0)
209 // First call to GenerateParticles(...):
210 // givenShellId is given by the process
211 {
212 provShellId = SelectTypeOfTransition(Z, givenShellId);
213
214 if ( provShellId >0)
215 {
216 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
217 //if (aParticle != 0) {
218 // G4cout << "****FLUO!_1**** "
219 // << aParticle->GetParticleDefinition()->GetParticleType()
220 // << " " << aParticle->GetKineticEnergy()/keV << G4endl ;}
221 }
222 else if ( provShellId == -1)
223 {
224 // G4cout << "Try to generate Auger 1" << G4endl;
225 aParticle = GenerateAuger(Z, givenShellId);
226 // if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;}
227 }
228 else
229 {
230 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()",
231 // "de0002",JustWarning, "Energy deposited locally");
232 }
233 }
234 else
235 // Following calls to GenerateParticles(...):
236 // newShellId is given by GenerateFluorescence(...)
237 {
238 provShellId = SelectTypeOfTransition(Z,newShellId);
239 if (provShellId >0)
240 {
241 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
242 //if (aParticle != 0) { G4cout << "****FLUO!_2****" << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl;} //debug
243 }
244 else if ( provShellId == -1)
245 {
246 // G4cout << "Try to generate Auger 2" << G4endl; //debug
247 aParticle = GenerateAuger(Z, newShellId);
248 // if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
249 }
250 else
251 {
252 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
253 }
254 }
255 counter++;
256 if (aParticle != 0)
257 {
258 vectorOfParticles->push_back(aParticle);
259 //G4cout << "Deexcitation Occurred!" << G4endl; //debug
260 }
261 else {provShellId = -2;}
262 }
263 while (provShellId > -2);
264 }
265 else
266 {
267 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
268 }
269
270 //G4cout << "---------FATTO!---------" << G4endl; //debug
271
272} // Auger cascade is not active
273
274//END OF ORIGINAL METHOD BY ALFONSO MANTERO
275//----------------------
276
277// NEW METHOD
278// Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
279
281{
282//----------------------
283
284 vacancyArray.push_back(givenShellId);
285
286 // let's check that 5<Z<100
287 if (Z<6 || Z>99){
288 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
289 return;
290 }
291
292 // as long as there is vacancy to be filled by either fluo or auger, stay in the loop.
293 while(!vacancyArray.empty()){
294
295// prepare to process the last element, and then delete it from the vector.
296 givenShellId = vacancyArray[0];
297 provShellId = SelectTypeOfTransition(Z,givenShellId);
298
299 //G4cout<<"\n------ Atom Transition with Z: "<<Z<<"\tbetween current:"
300 // <<givenShellId<<" & target:"<<provShellId<<G4endl;
301 if(provShellId>0){
302 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
303// if (aParticle != 0) {
304// G4cout << "****FLUO!_1**** "
305// << aParticle->GetParticleDefinition()->GetParticleType()
306// << " " << aParticle->GetKineticEnergy()/keV << G4endl ;}
307 }
308 else if(provShellId == -1){
309 aParticle = GenerateAuger(Z, givenShellId);
310// if (aParticle != 0) { G4cout << "****AUGER!****" <<
311// aParticle->GetParticleDefinition()->GetParticleType()
312// << " " << aParticle->GetKineticEnergy()/keV << G4endl ; }
313// else G4cout<<G4endl;
314 }
315 //else
316 // G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
317
318// if a particle is created, put it in the vector of new particles
319 if(aParticle!=0)
320 vectorOfParticles->push_back(aParticle);
321 else{;}
322// one vacancy has been processed. Erase it.
323 vacancyArray.erase(vacancyArray.begin());
324 }
325
326
327//----------------------
328//End of Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
329
330} // Auger cascade is active
331
332//ENDSI
333}
334
337 const G4ParticleDefinition* pdef,
338 G4int Z,
339 G4AtomicShellEnumerator shellEnum,
340 G4double kineticEnergy,
341 const G4Material* mat)
342{
343 // we must put a control on the shell that are passed:
344 // some shells should not pass (line "0" or "2")
345 //G4cout << pdef->GetParticleName() << " Z= " << Z << " Shell= " << shellEnum
346 // << " E= " << kineticEnergy << G4endl;
347
348 // check atomic number
349 G4double xsec = 0.0;
350 if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
351 G4int idx = G4int(shellEnum);
352 if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
353
354 //G4cout << pdef->GetParticleName() << " Z= " << Z << " " << PIXEshellCS
355 // << " " << ePIXEshellCS << G4endl;
356 //
357 if(pdef == theElectron || pdef == thePositron) {
358 xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
359 return xsec;
360 }
361
362 G4double mass = pdef->GetPDGMass();
363 G4double escaled = kineticEnergy;
364 G4double q2 = 0.0;
365
366 // scaling to protons
367 if ((pdef->GetParticleName() != "proton" && pdef->GetParticleName() != "alpha" ) )
368 {
369 mass = proton_mass_c2;
370 escaled = kineticEnergy*mass/(pdef->GetPDGMass());
371
372 if(mat) {
373 q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
374 } else {
375 G4double q = pdef->GetPDGCharge()/eplus;
376 q2 = q*q;
377 }
378 }
379
380 if(PIXEshellCS) { xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat); }
381 if(xsec < 1e-100) {
382
383 xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
384
385 }
386
387 if (q2) {xsec *= q2;}
388
389 return xsec;
390}
391
393{
394 minGammaEnergy = cut;
395}
396
398{
399 minElectronEnergy = cut;
400}
401
404 const G4ParticleDefinition* p,
405 G4int Z,
407 G4double kinE,
408 const G4Material* mat)
409{
410 return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
411}
412
413G4int G4UAtomicDeexcitation::SelectTypeOfTransition(G4int Z, G4int shellId)
414{
415 if (shellId <=0 ) {
416 //G4Exception("G4UAtomicDeexcitation::SelecttypeOfTransition()","de0002",
417 // JustWarning, "Energy deposited locally");
418 return 0;
419 }
420 //G4bool fluoTransitionFoundFlag = false;
421
422 G4int provShellId = -1;
423 G4int shellNum = 0;
424 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
425
426 const G4FluoTransition* refShell =
427 transitionManager->ReachableShell(Z,maxNumOfShells-1);
428
429 // This loop gives shellNum the value of the index of shellId
430 // in the vector storing the list of the shells reachable through
431 // a radiative transition
432 if ( shellId <= refShell->FinalShellId())
433 {
434 while (shellId != transitionManager->ReachableShell(Z,shellNum)->FinalShellId())
435 {
436 if(shellNum ==maxNumOfShells-1)
437 {
438 break;
439 }
440 shellNum++;
441 }
442 G4int transProb = 0; //AM change 29/6/07 was 1
443
444 G4double partialProb = G4UniformRand();
445 G4double partSum = 0;
446 const G4FluoTransition* aShell = transitionManager->ReachableShell(Z,shellNum);
447 G4int trSize = (aShell->TransitionProbabilities()).size();
448
449 // Loop over the shells wich can provide an electron for a
450 // radiative transition towards shellId:
451 // in every loop the partial sum of the first transProb shells
452 // is calculated and compared with a random number [0,1].
453 // If the partial sum is greater, the shell whose index is transProb
454 // is chosen as the starting shell for a radiative transition
455 // and its identity is returned
456 // Else, terminateded the loop, -1 is returned
457 while(transProb < trSize){
458
459 partSum += aShell->TransitionProbability(transProb);
460
461 if(partialProb <= partSum)
462 {
463 provShellId = aShell->OriginatingShellId(transProb);
464 //fluoTransitionFoundFlag = true;
465
466 break;
467 }
468 transProb++;
469 }
470
471 // here provShellId is the right one or is -1.
472 // if -1, the control is passed to the Auger generation part of the package
473 }
474 else
475 {
476 provShellId = -1;
477 }
478 //G4cout << "FlagTransition= " << provShellId << " ecut(MeV)= " << minElectronEnergy
479 // << " gcut(MeV)= " << minGammaEnergy << G4endl;
480 return provShellId;
481}
482
484G4UAtomicDeexcitation::GenerateFluorescence(G4int Z, G4int shellId,
485 G4int provShellId )
486{
487 if (shellId <=0 )
488 {
489 //G4Exception("G4UAtomicDeexcitation::GenerateFluorescence()","de0002",JustWarning, "Energy deposited locally");
490 return 0;
491 }
492
493
494 //isotropic angular distribution for the outcoming photon
495 G4double newcosTh = 1.-2.*G4UniformRand();
496 G4double newsinTh = std::sqrt((1.-newcosTh)*(1. + newcosTh));
497 G4double newPhi = twopi*G4UniformRand();
498
499 G4double xDir = newsinTh*std::sin(newPhi);
500 G4double yDir = newsinTh*std::cos(newPhi);
501 G4double zDir = newcosTh;
502
503 G4ThreeVector newGammaDirection(xDir,yDir,zDir);
504
505 G4int shellNum = 0;
506 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
507
508 // find the index of the shell named shellId
509 while (shellId != transitionManager->
510 ReachableShell(Z,shellNum)->FinalShellId())
511 {
512 if(shellNum == maxNumOfShells-1)
513 {
514 break;
515 }
516 shellNum++;
517 }
518 // number of shell from wich an electron can reach shellId
519 size_t transitionSize = transitionManager->
520 ReachableShell(Z,shellNum)->OriginatingShellIds().size();
521
522 size_t index = 0;
523
524 // find the index of the shell named provShellId in the vector
525 // storing the shells from which shellId can be reached
526 while (provShellId != transitionManager->
527 ReachableShell(Z,shellNum)->OriginatingShellId(index))
528 {
529 if(index == transitionSize-1)
530 {
531 break;
532 }
533 index++;
534 }
535 // energy of the gamma leaving provShellId for shellId
536 G4double transitionEnergy = transitionManager->
537 ReachableShell(Z,shellNum)->TransitionEnergy(index);
538
539 if (transitionEnergy < minGammaEnergy) return 0;
540
541 // This is the shell where the new vacancy is: it is the same
542 // shell where the electron came from
543 newShellId = transitionManager->
544 ReachableShell(Z,shellNum)->OriginatingShellId(index);
545
546
548 newGammaDirection,
549 transitionEnergy);
550 //SI
551 //Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
552 if (IsAugerCascadeActive()) vacancyArray.push_back(newShellId);
553 //ENDSI
554
555 return newPart;
556}
557
558G4DynamicParticle* G4UAtomicDeexcitation::GenerateAuger(G4int Z, G4int shellId)
559{
560 if(!IsAugerActive()) {
561 // G4cout << "auger inactive!" << G4endl; //debug
562 return 0;
563 }
564
565 if (shellId <=0 ) {
566 //G4Exception("G4UAtomicDeexcitation::GenerateAuger()","de0002",
567 // JustWarning, "Energy deposited locally");
568 return 0;
569 }
570
571 // G4int provShellId = -1;
572 G4int maxNumOfShells = transitionManager->NumberOfReachableAugerShells(Z);
573
574 const G4AugerTransition* refAugerTransition =
575 transitionManager->ReachableAugerShell(Z,maxNumOfShells-1);
576
577 // This loop gives to shellNum the value of the index of shellId
578 // in the vector storing the list of the vacancies in the variuos shells
579 // that can originate a NON-radiative transition
580
581 G4int shellNum = 0;
582
583 if ( shellId <= refAugerTransition->FinalShellId() )
584 // "FinalShellId" is final from the point of view of the electron
585 // who makes the transition,
586 // being the Id of the shell in which there is a vacancy
587 {
588 G4int pippo = transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId();
589 if (shellId != pippo ) {
590 do {
591 shellNum++;
592 if(shellNum == maxNumOfShells)
593 {
594 // G4cout << "No Auger transition found" << G4endl; //debug
595 return 0;
596 }
597 }
598 while (shellId != (transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId()) );
599 }
600
601 // Now we have that shellnum is the shellIndex of the shell named ShellId
602 // G4cout << " the index of the shell is: "<<shellNum<<G4endl;
603 // But we have now to select two shells: one for the transition,
604 // and another for the auger emission.
605
606 G4int transitionLoopShellIndex = 0;
607 G4double partSum = 0;
608 const G4AugerTransition* anAugerTransition =
609 transitionManager->ReachableAugerShell(Z,shellNum);
610
611 //G4cout << " corresponding to the ID: "
612 //<< anAugerTransition->FinalShellId()<< G4endl;
613
614 G4int transitionSize =
615 (anAugerTransition->TransitionOriginatingShellIds())->size();
616 while (transitionLoopShellIndex < transitionSize) {
617
618 std::vector<G4int>::const_iterator pos =
619 anAugerTransition->TransitionOriginatingShellIds()->begin();
620
621 G4int transitionLoopShellId = *(pos+transitionLoopShellIndex);
622 G4int numberOfPossibleAuger =
623 (anAugerTransition->AugerTransitionProbabilities(transitionLoopShellId))->size();
624 G4int augerIndex = 0;
625 // G4int partSum2 = 0;
626
627 if (augerIndex < numberOfPossibleAuger) {
628 do
629 {
630 G4double thisProb = anAugerTransition->AugerTransitionProbability(augerIndex,
631 transitionLoopShellId);
632 partSum += thisProb;
633 augerIndex++;
634
635 } while (augerIndex < numberOfPossibleAuger);
636 }
637 transitionLoopShellIndex++;
638 }
639
640 // Now we have the entire probability of an auger transition for the vacancy
641 // located in shellNum (index of shellId)
642
643 // AM *********************** F I X E D **************************** AM
644 // Here we duplicate the previous loop, this time looking to the sum of the probabilities
645 // to be under the random number shoot by G4 UniformRdandom. This could have been done in the
646 // previuos loop, while integrating the probabilities. There is a bug that will be fixed
647 // 5 minutes from now: a line:
648 // G4int numberOfPossibleAuger = (anAugerTransition->
649 // AugerTransitionProbabilities(transitionLoopShellId))->size();
650 // to be inserted.
651 // AM *********************** F I X E D **************************** AM
652
653 // Remains to get the same result with a single loop.
654
655 // AM *********************** F I X E D **************************** AM
656 // Another Bug: in EADL Auger Transition are normalized to all the transitions deriving from
657 // a vacancy in one shell, but not all of these are present in data tables. So if a transition
658 // doesn't occur in the main one a local energy deposition must occur, instead of (like now)
659 // generating the last transition present in EADL data.
660 // AM *********************** F I X E D **************************** AM
661
662 G4double totalVacancyAugerProbability = partSum;
663
664 //And now we start to select the right auger transition and emission
665 G4int transitionRandomShellIndex = 0;
666 G4int transitionRandomShellId = 1;
667 G4int augerIndex = 0;
668 partSum = 0;
669 G4double partialProb = G4UniformRand();
670 // G4int augerOriginatingShellId = 0;
671
672 G4int numberOfPossibleAuger = 0;
673
674 G4bool foundFlag = false;
675
676 while (transitionRandomShellIndex < transitionSize) {
677
678 std::vector<G4int>::const_iterator pos =
679 anAugerTransition->TransitionOriginatingShellIds()->begin();
680
681 transitionRandomShellId = *(pos+transitionRandomShellIndex);
682
683 augerIndex = 0;
684 numberOfPossibleAuger = (anAugerTransition->
685 AugerTransitionProbabilities(transitionRandomShellId))->size();
686
687 while (augerIndex < numberOfPossibleAuger) {
688 G4double thisProb =anAugerTransition->AugerTransitionProbability(augerIndex,
689 transitionRandomShellId);
690
691 partSum += thisProb;
692
693 if (partSum >= (partialProb*totalVacancyAugerProbability) ) { // was /
694 foundFlag = true;
695 break;
696 }
697 augerIndex++;
698 }
699 if (partSum >= (partialProb*totalVacancyAugerProbability) ) {break;} // was /
700 transitionRandomShellIndex++;
701 }
702
703 // Now we have the index of the shell from wich comes the auger electron (augerIndex),
704 // and the id of the shell, from which the transition e- come (transitionRandomShellid)
705 // If no Transition has been found, 0 is returned.
706
707 if (!foundFlag) {
708 // G4cout << "Auger not found (foundflag = false) " << G4endl; //debug
709 return 0;
710 }
711
712 // Isotropic angular distribution for the outcoming e-
713 G4double newcosTh = 1.-2.*G4UniformRand();
714 G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
715 G4double newPhi = twopi*G4UniformRand();
716
717 G4double xDir = newsinTh*std::sin(newPhi);
718 G4double yDir = newsinTh*std::cos(newPhi);
719 G4double zDir = newcosTh;
720
721 G4ThreeVector newElectronDirection(xDir,yDir,zDir);
722
723 // energy of the auger electron emitted
724
725 G4double transitionEnergy =
726 anAugerTransition->AugerTransitionEnergy(augerIndex, transitionRandomShellId);
727 /*
728 G4cout << "AUger TransitionId " << anAugerTransition->FinalShellId() << G4endl;
729 G4cout << "augerIndex: " << augerIndex << G4endl;
730 G4cout << "transitionShellId: " << transitionRandomShellId << G4endl;
731 */
732
733 if (transitionEnergy < minElectronEnergy) {
734 // G4cout << "Problem! (transitionEnergy < minElectronEnergy)" << G4endl; // debug
735 // G4cout << "minElectronEnergy(KeV): " << minElectronEnergy/keV << G4endl; // debug
736 // G4cout << "transitionEnergy(KeV): " << transitionEnergy/keV << G4endl; // debug
737 return 0;
738 }
739
740 // This is the shell where the new vacancy is: it is the same
741 // shell where the electron came from
742 newShellId = transitionRandomShellId;
743
744 //SI
745 //Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
747 {
748 vacancyArray.push_back(newShellId);
749 vacancyArray.push_back(anAugerTransition->AugerOriginatingShellId(augerIndex,transitionRandomShellId));
750 }
751 //ENDSI
752
754 newElectronDirection,
755 transitionEnergy);
756 }
757 else
758 {
759 // G4cout << "G4UAtomicDeexcitation: no auger transition found" << G4endl ;
760 // G4cout << "( shellId <= refAugerTransition->FinalShellId() )" << G4endl;
761 return 0;
762 }
763}
G4AtomicShellEnumerator
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#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
static G4AtomicTransitionManager * Instance()
G4int NumberOfReachableAugerShells(G4int Z) const
G4int AugerOriginatingShellId(G4int index, G4int startShellId) 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:93
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
static G4EmParameters * Instance()
const G4String & PIXECrossSectionModel()
const G4String & PIXEElectronCrossSectionModel()
const G4DataVector & TransitionProbabilities() const
G4int OriginatingShellId(G4int index) const
G4double TransitionProbability(G4int index) const
G4int FinalShellId() const
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)
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
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