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