Geant4  10.01.p02
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 
228 const G4AtomicShell*
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: "
245  // << givenShellId << G4endl; // debug
246  minGammaEnergy = gammaCut;
247  minElectronEnergy = eCut;
248 
249  // V.I. short-cut
250  // if(!IsAugerActive()) { minElectronEnergy = DBL_MAX; }
251 
252  // generation secondaries
253  G4DynamicParticle* aParticle=0;
254  G4int provShellId = 0;
255  G4int counter = 0;
256 
257  // let's check that 5<Z<100
258 
259  if (Z>5 && Z<100) {
260 
261  // The aim of this loop is to generate more than one fluorecence photon
262  // from the same ionizing event
263  do
264  {
265  if (counter == 0)
266  // First call to GenerateParticles(...):
267  // givenShellId is given by the process
268  {
269  provShellId = SelectTypeOfTransition(Z, givenShellId);
270 
271  if ( provShellId >0)
272  {
273  aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
274  //if (aParticle != 0) {
275  // G4cout << "****FLUO!_1**** "
276  // << aParticle->GetParticleDefinition()->GetParticleType()
277  // << " " << aParticle->GetKineticEnergy()/keV << G4endl ;}
278  }
279  else if ( provShellId == -1)
280  {
281  // G4cout << "Try to generate Auger 1" << G4endl;
282  aParticle = GenerateAuger(Z, givenShellId);
283  // if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;}
284  }
285  else
286  {
287  G4Exception("G4UAtomicDeexcitation::GenerateParticles()",
288  "de0002",JustWarning, "Energy deposited locally");
289  }
290  }
291  else
292  // Following calls to GenerateParticles(...):
293  // newShellId is given by GenerateFluorescence(...)
294  {
295  provShellId = SelectTypeOfTransition(Z,newShellId);
296  if (provShellId >0)
297  {
298  aParticle = GenerateFluorescence(Z,newShellId,provShellId);
299  //if (aParticle != 0) { G4cout << "****FLUO!_2****" << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl;} //debug
300  }
301  else if ( provShellId == -1)
302  {
303  // G4cout << "Try to generate Auger 2" << G4endl; //debug
304  aParticle = GenerateAuger(Z, newShellId);
305  // if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
306  }
307  else
308  {
309  G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
310  }
311  }
312  counter++;
313  if (aParticle != 0)
314  {
315  vectorOfParticles->push_back(aParticle);
316  //G4cout << "Deexcitation Occurred!" << G4endl; //debug
317  }
318  else {provShellId = -2;}
319  }
320  while (provShellId > -2);
321  }
322  else
323  {
324  G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
325  }
326 
327  //G4cout << "---------FATTO!---------" << G4endl; //debug
328 
329 }
330 
331 G4double
333  const G4ParticleDefinition* pdef,
334  G4int Z,
335  G4AtomicShellEnumerator shellEnum,
336  G4double kineticEnergy,
337  const G4Material* mat)
338 {
339 
340  // we must put a control on the shell that are passed:
341  // some shells should not pass (line "0" or "2")
342 
343  // check atomic number
344  G4double xsec = 0.0;
345  if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
346  G4int idx = G4int(shellEnum);
347  if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
348 
349  //
350  if(pdef == theElectron || pdef == thePositron) {
351  xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
352  return xsec;
353  }
354 
355  G4double mass = pdef->GetPDGMass();
356  G4double escaled = kineticEnergy;
357  G4double q2 = 0.0;
358 
359  // scaling to protons
360  if ((pdef->GetParticleName() != "proton" && pdef->GetParticleName() != "alpha" ) )
361  {
362  mass = proton_mass_c2;
363  escaled = kineticEnergy*mass/(pdef->GetPDGMass());
364 
365  if(mat) {
366  q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
367  } else {
368  G4double q = pdef->GetPDGCharge()/eplus;
369  q2 = q*q;
370  }
371  }
372 
373  if(PIXEshellCS) { xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat); }
374  if(xsec < 1e-100) {
375 
376  xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
377 
378  }
379 
380  if (q2) {xsec *= q2;}
381 
382  return xsec;
383 }
384 
386 {
387  minGammaEnergy = cut;
388 }
389 
391 {
392  minElectronEnergy = cut;
393 }
394 
395 G4double
397  const G4ParticleDefinition* p,
398  G4int Z,
400  G4double kinE,
401  const G4Material* mat)
402 {
403  return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
404 }
405 
407 {
408  if (shellId <=0 ) {
409  G4Exception("G4UAtomicDeexcitation::SelecttypeOfTransition()","de0002",JustWarning, "Energy deposited locally");
410  return 0;
411  }
412  //G4bool fluoTransitionFoundFlag = false;
413 
414  G4int provShellId = -1;
415  G4int shellNum = 0;
416  G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
417 
418  const G4FluoTransition* refShell = transitionManager->ReachableShell(Z,maxNumOfShells-1);
419 
420  // This loop gives shellNum the value of the index of shellId
421  // in the vector storing the list of the shells reachable through
422  // a radiative transition
423  if ( shellId <= refShell->FinalShellId())
424  {
425  while (shellId != transitionManager->ReachableShell(Z,shellNum)->FinalShellId())
426  {
427  if(shellNum ==maxNumOfShells-1)
428  {
429  break;
430  }
431  shellNum++;
432  }
433  G4int transProb = 0; //AM change 29/6/07 was 1
434 
435  G4double partialProb = G4UniformRand();
436  G4double partSum = 0;
437  const G4FluoTransition* aShell = transitionManager->ReachableShell(Z,shellNum);
438  G4int trSize = (aShell->TransitionProbabilities()).size();
439 
440  // Loop over the shells wich can provide an electron for a
441  // radiative transition towards shellId:
442  // in every loop the partial sum of the first transProb shells
443  // is calculated and compared with a random number [0,1].
444  // If the partial sum is greater, the shell whose index is transProb
445  // is chosen as the starting shell for a radiative transition
446  // and its identity is returned
447  // Else, terminateded the loop, -1 is returned
448  while(transProb < trSize){
449 
450  partSum += aShell->TransitionProbability(transProb);
451 
452  if(partialProb <= partSum)
453  {
454  provShellId = aShell->OriginatingShellId(transProb);
455  //fluoTransitionFoundFlag = true;
456 
457  break;
458  }
459  transProb++;
460  }
461 
462  // here provShellId is the right one or is -1.
463  // if -1, the control is passed to the Auger generation part of the package
464  }
465  else
466  {
467  provShellId = -1;
468  }
469  //G4cout << "FlagTransition= " << provShellId << " ecut(MeV)= " << minElectronEnergy
470  // << " gcut(MeV)= " << minGammaEnergy << G4endl;
471  return provShellId;
472 }
473 
476  G4int provShellId )
477 {
478 
479  //if(!IsDeexActive()) { return 0; }
480 
481  if (shellId <=0 )
482 
483  {
484  G4Exception("G4UAtomicDeexcitation::GenerateFluorescence()","de0002",JustWarning, "Energy deposited locally");
485  return 0;
486  }
487 
488 
489  //isotropic angular distribution for the outcoming photon
490  G4double newcosTh = 1.-2.*G4UniformRand();
491  G4double newsinTh = std::sqrt((1.-newcosTh)*(1. + newcosTh));
492  G4double newPhi = twopi*G4UniformRand();
493 
494  G4double xDir = newsinTh*std::sin(newPhi);
495  G4double yDir = newsinTh*std::cos(newPhi);
496  G4double zDir = newcosTh;
497 
498  G4ThreeVector newGammaDirection(xDir,yDir,zDir);
499 
500  G4int shellNum = 0;
501  G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
502 
503  // find the index of the shell named shellId
504  while (shellId != transitionManager->
505  ReachableShell(Z,shellNum)->FinalShellId())
506  {
507  if(shellNum == maxNumOfShells-1)
508  {
509  break;
510  }
511  shellNum++;
512  }
513  // number of shell from wich an electron can reach shellId
514  size_t transitionSize = transitionManager->
515  ReachableShell(Z,shellNum)->OriginatingShellIds().size();
516 
517  size_t index = 0;
518 
519  // find the index of the shell named provShellId in the vector
520  // storing the shells from which shellId can be reached
521  while (provShellId != transitionManager->
522  ReachableShell(Z,shellNum)->OriginatingShellId(index))
523  {
524  if(index == transitionSize-1)
525  {
526  break;
527  }
528  index++;
529  }
530  // energy of the gamma leaving provShellId for shellId
531  G4double transitionEnergy = transitionManager->
532  ReachableShell(Z,shellNum)->TransitionEnergy(index);
533 
534  if (transitionEnergy < minGammaEnergy) return 0;
535 
536  // This is the shell where the new vacancy is: it is the same
537  // shell where the electron came from
539  ReachableShell(Z,shellNum)->OriginatingShellId(index);
540 
541 
543  newGammaDirection,
544  transitionEnergy);
545  return newPart;
546 }
547 
549 {
550  if(!IsAugerActive()) {
551  // G4cout << "auger inactive!" << G4endl; //debug
552  return 0;
553  }
554 
555  if (shellId <=0 ) {
556 
557 
558  G4Exception("G4UAtomicDeexcitation::GenerateAuger()","de0002",JustWarning, "Energy deposited locally");
559 
560  return 0;
561 
562  }
563 
564  // G4int provShellId = -1;
566 
567  const G4AugerTransition* refAugerTransition =
568  transitionManager->ReachableAugerShell(Z,maxNumOfShells-1);
569 
570  // This loop gives to shellNum the value of the index of shellId
571  // in the vector storing the list of the vacancies in the variuos shells
572  // that can originate a NON-radiative transition
573 
574  G4int shellNum = 0;
575 
576  if ( shellId <= refAugerTransition->FinalShellId() )
577  //"FinalShellId" is final from the point of view of the elctron who makes the transition,
578  // being the Id of the shell in which there is a vacancy
579  {
581  if (shellId != pippo ) {
582  do {
583  shellNum++;
584  if(shellNum == maxNumOfShells)
585  {
586  // G4cout << "No Auger transition found" << G4endl; //debug
587  return 0;
588  }
589  }
590  while (shellId != (transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId()) ) ;
591  }
592 
593 
594  // Now we have that shellnum is the shellIndex of the shell named ShellId
595 
596  // G4cout << " the index of the shell is: "<<shellNum<<G4endl;
597 
598  // But we have now to select two shells: one for the transition,
599  // and another for the auger emission.
600 
601  G4int transitionLoopShellIndex = 0;
602  G4double partSum = 0;
603  const G4AugerTransition* anAugerTransition =
605 
606  // G4cout << " corresponding to the ID: "<< anAugerTransition->FinalShellId() << G4endl;
607 
608 
609  G4int transitionSize =
610  (anAugerTransition->TransitionOriginatingShellIds())->size();
611  while (transitionLoopShellIndex < transitionSize) {
612 
613  std::vector<G4int>::const_iterator pos =
614  anAugerTransition->TransitionOriginatingShellIds()->begin();
615 
616  G4int transitionLoopShellId = *(pos+transitionLoopShellIndex);
617  G4int numberOfPossibleAuger =
618  (anAugerTransition->AugerTransitionProbabilities(transitionLoopShellId))->size();
619  G4int augerIndex = 0;
620  // G4int partSum2 = 0;
621 
622 
623  if (augerIndex < numberOfPossibleAuger) {
624 
625  do
626  {
627  G4double thisProb = anAugerTransition->AugerTransitionProbability(augerIndex,
628  transitionLoopShellId);
629  partSum += thisProb;
630  augerIndex++;
631 
632  } while (augerIndex < numberOfPossibleAuger);
633  }
634  transitionLoopShellIndex++;
635  }
636 
637 
638 
639  // Now we have the entire probability of an auger transition for the vacancy
640  // located in shellNum (index of shellId)
641 
642  // AM *********************** F I X E D **************************** AM
643  // Here we duplicate the previous loop, this time looking to the sum of the probabilities
644  // to be under the random number shoot by G4 UniformRdandom. This could have been done in the
645  // previuos loop, while integrating the probabilities. There is a bug that will be fixed
646  // 5 minutes from now: a line:
647  // G4int numberOfPossibleAuger = (anAugerTransition->
648  // AugerTransitionProbabilities(transitionLoopShellId))->size();
649  // to be inserted.
650  // AM *********************** F I X E D **************************** AM
651 
652  // Remains to get the same result with a single loop.
653 
654  // AM *********************** F I X E D **************************** AM
655  // Another Bug: in EADL Auger Transition are normalized to all the transitions deriving from
656  // a vacancy in one shell, but not all of these are present in data tables. So if a transition
657  // doesn't occur in the main one a local energy deposition must occur, instead of (like now)
658  // generating the last transition present in EADL data.
659  // AM *********************** F I X E D **************************** AM
660 
661 
662  G4double totalVacancyAugerProbability = partSum;
663 
664 
665  //And now we start to select the right auger transition and emission
666  G4int transitionRandomShellIndex = 0;
667  G4int transitionRandomShellId = 1;
668  G4int augerIndex = 0;
669  partSum = 0;
670  G4double partialProb = G4UniformRand();
671  // G4int augerOriginatingShellId = 0;
672 
673  G4int numberOfPossibleAuger = 0;
674 
675  G4bool foundFlag = false;
676 
677  while (transitionRandomShellIndex < transitionSize) {
678 
679  std::vector<G4int>::const_iterator pos =
680  anAugerTransition->TransitionOriginatingShellIds()->begin();
681 
682  transitionRandomShellId = *(pos+transitionRandomShellIndex);
683 
684  augerIndex = 0;
685  numberOfPossibleAuger = (anAugerTransition->
686  AugerTransitionProbabilities(transitionRandomShellId))->size();
687 
688  while (augerIndex < numberOfPossibleAuger) {
689  G4double thisProb =anAugerTransition->AugerTransitionProbability(augerIndex,
690  transitionRandomShellId);
691 
692  partSum += thisProb;
693 
694  if (partSum >= (partialProb*totalVacancyAugerProbability) ) { // was /
695  foundFlag = true;
696  break;
697  }
698  augerIndex++;
699  }
700  if (partSum >= (partialProb*totalVacancyAugerProbability) ) {break;} // was /
701  transitionRandomShellIndex++;
702  }
703 
704  // Now we have the index of the shell from wich comes the auger electron (augerIndex),
705  // and the id of the shell, from which the transition e- come (transitionRandomShellid)
706  // If no Transition has been found, 0 is returned.
707 
708  if (!foundFlag) {
709  // G4cout << "Auger not found (foundflag = false) " << G4endl; //debug
710  return 0;
711 
712 }
713 
714  // Isotropic angular distribution for the outcoming e-
715  G4double newcosTh = 1.-2.*G4UniformRand();
716  G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
717  G4double newPhi = twopi*G4UniformRand();
718 
719  G4double xDir = newsinTh*std::sin(newPhi);
720  G4double yDir = newsinTh*std::cos(newPhi);
721  G4double zDir = newcosTh;
722 
723  G4ThreeVector newElectronDirection(xDir,yDir,zDir);
724 
725  // energy of the auger electron emitted
726 
727 
728  G4double transitionEnergy = anAugerTransition->AugerTransitionEnergy(augerIndex, transitionRandomShellId);
729  /*
730  G4cout << "AUger TransitionId " << anAugerTransition->FinalShellId() << G4endl;
731  G4cout << "augerIndex: " << augerIndex << G4endl;
732  G4cout << "transitionShellId: " << transitionRandomShellId << G4endl;
733  */
734 
735  if (transitionEnergy < minElectronEnergy) {
736  // G4cout << "Problem! (transitionEnergy < minElectronEnergy)" << G4endl; // debug
737  // G4cout << "minElectronEnergy(KeV): " << minElectronEnergy/keV << G4endl; // debug
738  // G4cout << "transitionEnergy(KeV): " << transitionEnergy/keV << G4endl; // debug
739  return 0;
740  }
741 
742  // This is the shell where the new vacancy is: it is the same
743  // shell where the electron came from
744  newShellId = transitionRandomShellId;
745 
747  newElectronDirection,
748  transitionEnergy);
749  }
750  else
751  {
752  // G4cout << "G4UAtomicDeexcitation: no auger transition found" << G4endl ;
753  // G4cout << "( shellId <= refAugerTransition->FinalShellId() )" << G4endl;
754  return 0;
755  }
756 }
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
G4AtomicTransitionManager * transitionManager
G4int FinalShellId() const
G4int OriginatingShellId(G4int index) const
const G4String & PIXEElectronCrossSectionModel() const
virtual void InitialiseForNewRun()
void SetCutForSecondaryPhotons(G4double cut)
#define G4UniformRand()
Definition: Randomize.hh:93
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
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)