Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4AtomicDeexcitation.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// Authors: Elena Guardincerri (Elena.Guardincerri@ge.infn.it)
29// Alfonso Mantero (Alfonso.Mantero@ge.infn.it)
30//
31// History:
32// -----------
33//
34// 16 Sept 2001 First committed to cvs
35// 12 Sep 2003 Bug in auger production fixed
36//
37// -------------------------------------------------------------------
38
40#include "Randomize.hh"
42#include "G4SystemOfUnits.hh"
43#include "G4Gamma.hh"
44#include "G4Electron.hh"
46#include "G4FluoTransition.hh"
47
49 minGammaEnergy(100.*eV),
50 minElectronEnergy(100.*eV),
51 fAuger(false)
52{
53
54 G4cout << " ********************************************************** " << G4endl;
55 G4cout << " * W A R N I N G ! ! ! * " << G4endl;
56 G4cout << " ********************************************************** " << G4endl;
57 G4cout << " * * " << G4endl;
58 G4cout << " * Class G4AtomicDeexcitation is obsolete. It has been * " << G4endl;
59 G4cout << " * discontinued and is going to be removed by next Geant4 * " << G4endl;
60 G4cout << " * release please migrate to G4UAtomDeexcitation. * " << G4endl;
61 G4cout << " * * " << G4endl;
62 G4cout << " ********************************************************** " << G4endl;
63
64 augerVacancyId=0;
65 newShellId=0;
66}
67
69{}
70
71std::vector<G4DynamicParticle*>* G4AtomicDeexcitation::GenerateParticles(G4int Z,G4int givenShellId)
72{
73
74 std::vector<G4DynamicParticle*>* vectorOfParticles;
75 vectorOfParticles = new std::vector<G4DynamicParticle*>;
76
77 G4DynamicParticle* aParticle;
78 G4int provShellId = 0;
79 G4int counter = 0;
80
81 // The aim of this loop is to generate more than one fluorecence photon
82 // from the same ionizing event
83 do
84 {
85 if (counter == 0)
86 // First call to GenerateParticles(...):
87 // givenShellId is given by the process
88 {
89 provShellId = SelectTypeOfTransition(Z, givenShellId);
90
91 if ( provShellId >0)
92 {
93 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
94 }
95 else if ( provShellId == -1)
96 {
97 aParticle = GenerateAuger(Z, givenShellId);
98 }
99 else
100 {
101 G4Exception("G4AtomicDeexcitation::Constructor", "de0002", JustWarning, "Transition selection invalid, energy local deposited");
102 }
103 }
104 else
105 // Following calls to GenerateParticles(...):
106 // newShellId is given by GenerateFluorescence(...)
107 {
108 provShellId = SelectTypeOfTransition(Z,newShellId);
109 if (provShellId >0)
110 {
111 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
112 }
113 else if ( provShellId == -1)
114 {
115 aParticle = GenerateAuger(Z, newShellId);
116 }
117 else
118 {
119 G4Exception("G4AtomicDeexcitation::constructor", "de0002", JustWarning, "Transition selection invalid, energy local deposited" );
120 }
121 }
122 counter++;
123 if (aParticle != 0) {vectorOfParticles->push_back(aParticle);}
124 else {provShellId = -2;}
125 }
126
127 // Look this in a particular way: only one auger emitted! // ????
128 while (provShellId > -2);
129
130 // debug
131 // if (vectorOfParticles->size() > 0) {
132 // G4cout << " DEEXCITATION!" << G4endl;
133 // }
134
135 return vectorOfParticles;
136}
137
138G4int G4AtomicDeexcitation::SelectTypeOfTransition(G4int Z, G4int shellId)
139{
140 if (shellId <=0 )
141 {G4Exception("G4AtomicDeexcitation::SelectTypeOfTransition()","de0002", JustWarning ,"zero or negative shellId");}
142
143 //G4bool fluoTransitionFoundFlag = false;
144
145 const G4AtomicTransitionManager* transitionManager =
147 G4int provShellId = -1;
148 G4int shellNum = 0;
149 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
150
151 const G4FluoTransition* refShell = transitionManager->ReachableShell(Z,maxNumOfShells-1);
152
153 // This loop gives shellNum the value of the index of shellId
154 // in the vector storing the list of the shells reachable through
155 // a radiative transition
156 if ( shellId <= refShell->FinalShellId())
157 {
158 while (shellId != transitionManager->ReachableShell(Z,shellNum)->FinalShellId())
159 {
160 if(shellNum ==maxNumOfShells-1)
161 {
162 break;
163 }
164 shellNum++;
165 }
166 G4int transProb = 0; //AM change 29/6/07 was 1
167
168 G4double partialProb = G4UniformRand();
169 G4double partSum = 0;
170 const G4FluoTransition* aShell = transitionManager->ReachableShell(Z,shellNum);
171 G4int trSize = (aShell->TransitionProbabilities()).size();
172
173 // Loop over the shells wich can provide an electron for a
174 // radiative transition towards shellId:
175 // in every loop the partial sum of the first transProb shells
176 // is calculated and compared with a random number [0,1].
177 // If the partial sum is greater, the shell whose index is transProb
178 // is chosen as the starting shell for a radiative transition
179 // and its identity is returned
180 // Else, terminateded the loop, -1 is returned
181 while(transProb < trSize){
182
183 partSum += aShell->TransitionProbability(transProb);
184
185 if(partialProb <= partSum)
186 {
187 provShellId = aShell->OriginatingShellId(transProb);
188 //fluoTransitionFoundFlag = true;
189
190 break;
191 }
192 transProb++;
193 }
194
195 // here provShellId is the right one or is -1.
196 // if -1, the control is passed to the Auger generation part of the package
197 }
198
199
200
201 else
202 {
203
204 provShellId = -1;
205
206 }
207 return provShellId;
208}
209
210G4DynamicParticle* G4AtomicDeexcitation::GenerateFluorescence(G4int Z,
211 G4int shellId,
212 G4int provShellId )
213{
214
215
217 // G4int provenienceShell = provShellId;
218
219 //isotropic angular distribution for the outcoming photon
220 G4double newcosTh = 1.-2.*G4UniformRand();
221 G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
222 G4double newPhi = twopi*G4UniformRand();
223
224 G4double xDir = newsinTh*std::sin(newPhi);
225 G4double yDir = newsinTh*std::cos(newPhi);
226 G4double zDir = newcosTh;
227
228 G4ThreeVector newGammaDirection(xDir,yDir,zDir);
229
230 G4int shellNum = 0;
231 G4int maxNumOfShells = transitionManager->NumberOfReachableShells(Z);
232
233 // find the index of the shell named shellId
234 while (shellId != transitionManager->
235 ReachableShell(Z,shellNum)->FinalShellId())
236 {
237 if(shellNum == maxNumOfShells-1)
238 {
239 break;
240 }
241 shellNum++;
242 }
243 // number of shell from wich an electron can reach shellId
244 size_t transitionSize = transitionManager->
245 ReachableShell(Z,shellNum)->OriginatingShellIds().size();
246
247 size_t index = 0;
248
249 // find the index of the shell named provShellId in the vector
250 // storing the shells from which shellId can be reached
251 while (provShellId != transitionManager->
252 ReachableShell(Z,shellNum)->OriginatingShellId(index))
253 {
254 if(index == transitionSize-1)
255 {
256 break;
257 }
258 index++;
259 }
260 // energy of the gamma leaving provShellId for shellId
261 G4double transitionEnergy = transitionManager->
262 ReachableShell(Z,shellNum)->TransitionEnergy(index);
263
264 // This is the shell where the new vacancy is: it is the same
265 // shell where the electron came from
266 newShellId = transitionManager->
267 ReachableShell(Z,shellNum)->OriginatingShellId(index);
268
269
271 newGammaDirection,
272 transitionEnergy);
273 return newPart;
274}
275
276G4DynamicParticle* G4AtomicDeexcitation::GenerateAuger(G4int Z, G4int shellId)
277{
278 if(!fAuger) return 0;
279
280
281 const G4AtomicTransitionManager* transitionManager =
283
284
285
286 if (shellId <=0 )
287 {G4Exception("G4AtomicDeexcitation::GenerateAuger()","de0002", JustWarning ,"zero or negative shellId");}
288
289 // G4int provShellId = -1;
290 G4int maxNumOfShells = transitionManager->NumberOfReachableAugerShells(Z);
291
292 const G4AugerTransition* refAugerTransition =
293 transitionManager->ReachableAugerShell(Z,maxNumOfShells-1);
294
295
296 // This loop gives to shellNum the value of the index of shellId
297 // in the vector storing the list of the vacancies in the variuos shells
298 // that can originate a NON-radiative transition
299
300 // ---- MGP ---- Next line commented out to remove compilation warning
301 // G4int p = refAugerTransition->FinalShellId();
302
303 G4int shellNum = 0;
304
305
306 if ( shellId <= refAugerTransition->FinalShellId() )
307 //"FinalShellId" is final from the point of view of the elctron who makes the transition,
308 // being the Id of the shell in which there is a vacancy
309 {
310 G4int pippo = transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId();
311 if (shellId != pippo ) {
312 do {
313 shellNum++;
314 if(shellNum == maxNumOfShells)
315 {
316
317 //G4Exception("G4AtomicDeexcitation: No Auger transition found");
318 return 0;
319 }
320 }
321 while (shellId != (transitionManager->ReachableAugerShell(Z,shellNum)->FinalShellId()) ) ;
322 }
323
324
325 // Now we have that shellnum is the shellIndex of the shell named ShellId
326
327 // G4cout << " the index of the shell is: "<<shellNum<<G4endl;
328
329 // But we have now to select two shells: one for the transition,
330 // and another for the auger emission.
331
332 G4int transitionLoopShellIndex = 0;
333 G4double partSum = 0;
334 const G4AugerTransition* anAugerTransition =
335 transitionManager->ReachableAugerShell(Z,shellNum);
336
337 // G4cout << " corresponding to the ID: "<< anAugerTransition->FinalShellId() << G4endl;
338
339
340 G4int transitionSize =
341 (anAugerTransition->TransitionOriginatingShellIds())->size();
342 while (transitionLoopShellIndex < transitionSize) {
343
344 std::vector<G4int>::const_iterator pos =
345 anAugerTransition->TransitionOriginatingShellIds()->begin();
346
347 G4int transitionLoopShellId = *(pos+transitionLoopShellIndex);
348 G4int numberOfPossibleAuger =
349 (anAugerTransition->AugerTransitionProbabilities(transitionLoopShellId))->size();
350 G4int augerIndex = 0;
351 // G4int partSum2 = 0;
352
353
354 if (augerIndex < numberOfPossibleAuger) {
355
356 do
357 {
358 G4double thisProb = anAugerTransition->AugerTransitionProbability(augerIndex,
359 transitionLoopShellId);
360 partSum += thisProb;
361 augerIndex++;
362
363 } while (augerIndex < numberOfPossibleAuger);
364 }
365 transitionLoopShellIndex++;
366 }
367
368
369
370 // Now we have the entire probability of an auger transition for the vacancy
371 // located in shellNum (index of shellId)
372
373 // AM *********************** F I X E D **************************** AM
374 // Here we duplicate the previous loop, this time looking to the sum of the probabilities
375 // to be under the random number shoot by G4 UniformRdandom. This could have been done in the
376 // previuos loop, while integrating the probabilities. There is a bug that will be fixed
377 // 5 minutes from now: a line:
378 // G4int numberOfPossibleAuger = (anAugerTransition->
379 // AugerTransitionProbabilities(transitionLoopShellId))->size();
380 // to be inserted.
381 // AM *********************** F I X E D **************************** AM
382
383 // Remains to get the same result with a single loop.
384
385 // AM *********************** F I X E D **************************** AM
386 // Another Bug: in EADL Auger Transition are normalized to all the transitions deriving from
387 // a vacancy in one shell, but not all of these are present in data tables. So if a transition
388 // doesn't occur in the main one a local energy deposition must occur, instead of (like now)
389 // generating the last transition present in EADL data.
390 // AM *********************** F I X E D **************************** AM
391
392
393 G4double totalVacancyAugerProbability = partSum;
394
395
396 //And now we start to select the right auger transition and emission
397 G4int transitionRandomShellIndex = 0;
398 G4int transitionRandomShellId = 1;
399 G4int augerIndex = 0;
400 partSum = 0;
401 G4double partialProb = G4UniformRand();
402 // G4int augerOriginatingShellId = 0;
403
404 G4int numberOfPossibleAuger = 0;
405
406 G4bool foundFlag = false;
407
408 while (transitionRandomShellIndex < transitionSize) {
409
410 std::vector<G4int>::const_iterator pos =
411 anAugerTransition->TransitionOriginatingShellIds()->begin();
412
413 transitionRandomShellId = *(pos+transitionRandomShellIndex);
414
415 augerIndex = 0;
416 numberOfPossibleAuger = (anAugerTransition->
417 AugerTransitionProbabilities(transitionRandomShellId))->size();
418
419 while (augerIndex < numberOfPossibleAuger) {
420 G4double thisProb =anAugerTransition->AugerTransitionProbability(augerIndex,
421 transitionRandomShellId);
422
423 partSum += thisProb;
424
425 if (partSum >= (partialProb*totalVacancyAugerProbability) ) { // was /
426 foundFlag = true;
427 break;
428 }
429 augerIndex++;
430 }
431 if (partSum >= (partialProb*totalVacancyAugerProbability) ) {break;} // was /
432 transitionRandomShellIndex++;
433 }
434
435 // Now we have the index of the shell from wich comes the auger electron (augerIndex),
436 // and the id of the shell, from which the transition e- come (transitionRandomShellid)
437 // If no Transition has been found, 0 is returned.
438
439 if (!foundFlag) {return 0;}
440
441 // Isotropic angular distribution for the outcoming e-
442 G4double newcosTh = 1.-2.*G4UniformRand();
443 G4double newsinTh = std::sqrt(1.-newcosTh*newcosTh);
444 G4double newPhi = twopi*G4UniformRand();
445
446 G4double xDir = newsinTh*std::sin(newPhi);
447 G4double yDir = newsinTh*std::cos(newPhi);
448 G4double zDir = newcosTh;
449
450 G4ThreeVector newElectronDirection(xDir,yDir,zDir);
451
452 // energy of the auger electron emitted
453
454
455 G4double transitionEnergy = anAugerTransition->AugerTransitionEnergy(augerIndex, transitionRandomShellId);
456 /*
457 G4cout << "AUger TransitionId " << anAugerTransition->FinalShellId() << G4endl;
458 G4cout << "augerIndex: " << augerIndex << G4endl;
459 G4cout << "transitionShellId: " << transitionRandomShellId << G4endl;
460 */
461
462 // This is the shell where the new vacancy is: it is the same
463 // shell where the electron came from
464 newShellId = transitionRandomShellId;
465
466
468 newElectronDirection,
469 transitionEnergy);
470 return newPart;
471
472 }
473 else
474 {
475 //G4Exception("G4AtomicDeexcitation: no auger transition found");
476 return 0;
477 }
478
479}
480
482{
483 minGammaEnergy = cut;
484}
485
487{
488 minElectronEnergy = cut;
489}
490
492{
493 fAuger = val;
494}
495
496
497
498
499
500
501
@ JustWarning
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
void ActivateAugerElectronProduction(G4bool val)
void SetCutForSecondaryPhotons(G4double cut)
std::vector< G4DynamicParticle * > * GenerateParticles(G4int Z, G4int shellId)
void SetCutForAugerElectrons(G4double cut)
G4int NumberOfReachableShells(G4int Z) const
const G4AugerTransition * ReachableAugerShell(G4int Z, G4int shellIndex) const
const G4FluoTransition * ReachableShell(G4int Z, size_t shellIndex) const
static G4AtomicTransitionManager * Instance()
G4int NumberOfReachableAugerShells(G4int Z) const
G4int FinalShellId() const
const G4DataVector * AugerTransitionProbabilities(G4int startShellId) const
G4double AugerTransitionEnergy(G4int index, G4int startShellId) const
const std::vector< G4int > * TransitionOriginatingShellIds() const
G4double AugerTransitionProbability(G4int index, G4int startShellId) const
static G4Electron * Electron()
Definition: G4Electron.cc:94
const G4DataVector & TransitionProbabilities() const
G4int OriginatingShellId(G4int index) const
G4double TransitionProbability(G4int index) const
G4int FinalShellId() 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