Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4AugerTransition Class Reference

#include <G4AugerTransition.hh>

Public Member Functions

 G4AugerTransition (G4int finalShell, std::vector< G4int > transIds, const std::map< G4int, std::vector< G4int >, std::less< G4int > > *idMap, const std::map< G4int, G4DataVector, std::less< G4int > > *energyMap, const std::map< G4int, G4DataVector, std::less< G4int > > *probabilityMap)
 ~G4AugerTransition ()
const std::vector< G4int > * AugerOriginatingShellIds (G4int startShellId) const
const std::vector< G4int > * TransitionOriginatingShellIds () const
 Returns the ids of the shells from wich an electron cuuld fill the vacancy in finalShellId.
const G4DataVectorAugerTransitionEnergies (G4int startShellId) const
const G4DataVectorAugerTransitionProbabilities (G4int startShellId) const
G4int FinalShellId () const
 returns the id of the shell in wich the transition electron arrives
G4int AugerOriginatingShellId (G4int index, G4int startShellId) const
G4double AugerTransitionEnergy (G4int index, G4int startShellId) const
G4double AugerTransitionProbability (G4int index, G4int startShellId) const
G4int TransitionOriginatingShellId (G4int index) const
 Returns the id of the shell form wich the transition electron come from.

Detailed Description

Definition at line 51 of file G4AugerTransition.hh.

Constructor & Destructor Documentation

◆ G4AugerTransition()

G4AugerTransition::G4AugerTransition ( G4int finalShell,
std::vector< G4int > transIds,
const std::map< G4int, std::vector< G4int >, std::less< G4int > > * idMap,
const std::map< G4int, G4DataVector, std::less< G4int > > * energyMap,
const std::map< G4int, G4DataVector, std::less< G4int > > * probabilityMap )
explicit

Definition at line 44 of file G4AugerTransition.cc.

48{
49 finalShellId = finalShell;
50 augerOriginatingShellIdsMap = *idMap;
51 augerTransitionEnergiesMap = *energyMap;
52 augerTransitionProbabilitiesMap = *probabilityMap;
53 transitionOriginatingShellIds = std::move(transIds);
54}

◆ ~G4AugerTransition()

G4AugerTransition::~G4AugerTransition ( )

Definition at line 58 of file G4AugerTransition.cc.

59{;}

Member Function Documentation

◆ AugerOriginatingShellId()

G4int G4AugerTransition::AugerOriginatingShellId ( G4int index,
G4int startShellId ) const

Returns the id of the shell from wich come the auger electron , given the shell from wich the transition electron comes from and the index number.

Definition at line 132 of file G4AugerTransition.cc.

133{
134 const std::vector<G4int>* ids = AugerOriginatingShellIds(startShellId);
135 // G4int i =
136 std::vector<G4int>::const_iterator pos = ids->begin();
137 G4int n = 0;
138 n = *(pos+index);
139 return n;
140}
int G4int
Definition G4Types.hh:85
const std::vector< G4int > * AugerOriginatingShellIds(G4int startShellId) const

◆ AugerOriginatingShellIds()

const std::vector< G4int > * G4AugerTransition::AugerOriginatingShellIds ( G4int startShellId) const

All the data stored and provided by this class are relative to a given vacancy, whose identity is provided by the FinalShellId() method, in an atom of a given material Returns the ids of the shells from wich an auger electron culd came from, given the shell from wich the transition electron comes from.

Definition at line 64 of file G4AugerTransition.cc.

65{
66 auto shellId = augerOriginatingShellIdsMap.find(startShellId);
67 if (shellId == augerOriginatingShellIdsMap.end())
68 {
69 G4Exception("G4AugerTransition::AugerOriginatingShellIds()",
70 "em2199",JustWarning,"Error: no Auger ID found");
71 return nullptr;
72 }
73 const std::vector<G4int>* dataSet = &(*shellId).second;
74 if (dataSet->empty())
75 G4Exception("G4AugerTransition::AugerOriginatingShellIds()",
76 "em2198",JustWarning,"Error: no Auger ID found");
77 return dataSet;
78}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

Referenced by AugerOriginatingShellId().

◆ AugerTransitionEnergies()

const G4DataVector * G4AugerTransition::AugerTransitionEnergies ( G4int startShellId) const

Returns the energiess of the possible auger electrons, given th shell from wich the transition electron comes from.

Definition at line 91 of file G4AugerTransition.cc.

92{
93 auto shellId = augerTransitionEnergiesMap.find(startShellId);
94
95 if (shellId == augerTransitionEnergiesMap.end() )
96 {
97 G4Exception("G4AugerTransition::AugerTransitionEnergies()","de0002",JustWarning,
98 "corresponding map element not found, energy deposited locally");
99 return 0;
100 }
101 const G4DataVector* dataSet = &(*shellId).second;
102 return dataSet;
103}

Referenced by AugerTransitionEnergy().

◆ AugerTransitionEnergy()

G4double G4AugerTransition::AugerTransitionEnergy ( G4int index,
G4int startShellId ) const

Returns the energy of the auger electron, given the shell from wich the transition electron comes from and the index number.

Definition at line 145 of file G4AugerTransition.cc.

146{
147 const G4DataVector* energies = AugerTransitionEnergies(startShellId);
148 G4double energy = 0;
149 if (index < (G4int) energies->size()) {
150 G4DataVector::const_iterator pos = energies->begin();
151 energy = *(pos+index);
152 }
153 return energy;
154}
double G4double
Definition G4Types.hh:83
const G4DataVector * AugerTransitionEnergies(G4int startShellId) const
G4double energy(const ThreeVector &p, const G4double m)

◆ AugerTransitionProbabilities()

const G4DataVector * G4AugerTransition::AugerTransitionProbabilities ( G4int startShellId) const

Returns the emission probabilities of the auger electrons, given th shell from wich the transition electron comes from.

Definition at line 108 of file G4AugerTransition.cc.

109{
110 auto shellId = augerTransitionProbabilitiesMap.find(startShellId);
111 if (shellId == augerTransitionProbabilitiesMap.end() )
112 {
113
114 G4Exception("G4AugerTransition::AugerTransitionProbabilities()","de0002",
115 JustWarning,"corresponding map element not found, energy deposited locally");
116 return 0;
117 }
118
119 const G4DataVector* dataSet = &(*shellId).second;
120 return dataSet;
121}

Referenced by AugerTransitionProbability().

◆ AugerTransitionProbability()

G4double G4AugerTransition::AugerTransitionProbability ( G4int index,
G4int startShellId ) const

Returns the probability of the auger emission, given the shell from wich the transition electron comes from and the index number.

Definition at line 159 of file G4AugerTransition.cc.

160{
161
162 const G4DataVector *probabilities = AugerTransitionProbabilities(startShellId);
163 G4DataVector::const_iterator pos = probabilities->begin();
164
165 G4double probability = 0;
166 probability = *(pos+index);
167
168 return probability;
169}
const G4DataVector * AugerTransitionProbabilities(G4int startShellId) const

◆ FinalShellId()

G4int G4AugerTransition::FinalShellId ( ) const

returns the id of the shell in wich the transition electron arrives

Definition at line 124 of file G4AugerTransition.cc.

125{
126 return finalShellId;
127}

◆ TransitionOriginatingShellId()

G4int G4AugerTransition::TransitionOriginatingShellId ( G4int index) const

Returns the id of the shell form wich the transition electron come from.

Definition at line 172 of file G4AugerTransition.cc.

173{
174 return transitionOriginatingShellIds[index];
175}

◆ TransitionOriginatingShellIds()

const std::vector< G4int > * G4AugerTransition::TransitionOriginatingShellIds ( ) const

Returns the ids of the shells from wich an electron cuuld fill the vacancy in finalShellId.

Definition at line 82 of file G4AugerTransition.cc.

83{
84 const std::vector<G4int>* dataSet = &transitionOriginatingShellIds;
85 return dataSet;
86}

The documentation for this class was generated from the following files: