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