Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UAtomicDeexcitation.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// -------------------------------------------------------------------
28//
29// Geant4 Header G4UAtomicDeexcitation
30//
31// Authors: Alfonso Mantero ([email protected])
32//
33// Created 22 April 2010 from old G4AtomicDeexcitation class
34//
35// Modified:
36// ---------
37//
38//
39// -------------------------------------------------------------------
40//
41// Class description:
42// Implementation of atomic deexcitation
43//
44// -------------------------------------------------------------------
45
46#ifndef G4UAtomicDeexcitation_h
47#define G4UAtomicDeexcitation_h 1
48
50#include "G4AtomicShell.hh"
51#include "globals.hh"
52#include "G4DynamicParticle.hh"
53#include <vector>
54
57class G4EmCorrections;
58class G4Material;
59
61{
62public:
63 explicit G4UAtomicDeexcitation();
64 virtual ~G4UAtomicDeexcitation();
65
66 //=================================================================
67 // methods that are requested to be implemented by the interface
68 //=================================================================
69 /// initialisation methods
70 void InitialiseForNewRun() override;
71 void InitialiseForExtraAtom(G4int Z) override;
72
73 /// Set threshold energy for fluorescence
75
76 /// Set threshold energy for Auger electron production
78
79
80 /// Get atomic shell by shell index, used by discrete processes
81 /// (for example, photoelectric), when shell vacancy sampled by the model
83 G4AtomicShellEnumerator shell) override;
84
85 /// generation of deexcitation for given atom, shell vacancy and cuts
86 void GenerateParticles(std::vector<G4DynamicParticle*>* secVect,
87 const G4AtomicShell*,
88 G4int Z,
89 G4double gammaCut,
90 G4double eCut) override;
91
92 /// access or compute PIXE cross section
94 G4int Z,
96 G4double kinE,
97 const G4Material* mat = nullptr) override;
98
99 /// access or compute PIXE cross section
101 G4int Z,
103 G4double kinE,
104 const G4Material* mat = nullptr) override;
105
108
109private:
110 /// Decides wether a radiative transition is possible and, if it is,
111 /// returns the identity of the starting shell for the transition
112 G4int SelectTypeOfTransition(G4int Z, G4int shellId);
113
114 /// Generates a particle from a radiative transition and returns it
115 G4DynamicParticle* GenerateFluorescence(G4int Z, G4int shellId,
116 G4int provShellId);
117
118 /// Generates a particle from a non-radiative transition and returns it
119 G4DynamicParticle* GenerateAuger(G4int Z, G4int shellId);
120
121 ///Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
122 ///Generates auger electron cascade.
123 G4DynamicParticle* GenerateAuger(G4int Z, G4int shellId, G4int& newAugerShellId);
124 G4AtomicTransitionManager* transitionManager;
125
126 /// Data member for the calculation of the proton and alpha ionisation XS
127 G4VhShellCrossSection* PIXEshellCS;
128 G4VhShellCrossSection* anaPIXEshellCS;
129 G4VhShellCrossSection* ePIXEshellCS;
130 G4EmCorrections* emcorr;
131
132 const G4ParticleDefinition* theElectron;
133 const G4ParticleDefinition* thePositron;
134
135 //Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
136 //Data member to keep track of cascading vacancies.
137 std::vector<int> vacancyArray;
138
139 /// Data member which stores the shells to be filled by
140 /// the radiative transition
141 G4double minGammaEnergy;
142 G4double minElectronEnergy;
143 G4int newShellId;
144};
145
146#endif
147
148
149
150
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void SetCutForSecondaryPhotons(G4double cut)
Set threshold energy for fluorescence.
G4UAtomicDeexcitation(G4UAtomicDeexcitation &)=delete
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut) override
generation of deexcitation for given atom, shell vacancy and cuts
G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr) override
access or compute PIXE cross section
void InitialiseForNewRun() override
initialisation methods
G4UAtomicDeexcitation & operator=(const G4UAtomicDeexcitation &right)=delete
void SetCutForAugerElectrons(G4double cut)
Set threshold energy for Auger electron production.
void InitialiseForExtraAtom(G4int Z) override
const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell) override
G4double ComputeShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr) override
access or compute PIXE cross section