Geant4 11.1.1
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 = 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 126 of file G4AugerTransition.cc.

127{
128 const std::vector<G4int>* ids = AugerOriginatingShellIds(startShellId);
129 // G4int i =
130 std::vector<G4int>::const_iterator pos = ids->begin();
131 G4int n = 0;
132 n = *(pos+index);
133 return n;
134}
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
68 const std::vector<G4int>* dataSet = &(*shellId).second;
69 if (dataSet->empty())
70 G4cout << "Error: no auger Id found"<< G4endl;
71 return dataSet;
72}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

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 85 of file G4AugerTransition.cc.

86{
87 auto shellId = augerTransitionEnergiesMap.find(startShellId);
88
89 if (shellId == augerTransitionEnergiesMap.end() )
90 {
91 G4Exception("G4AugerTransition::AugerTransitionEnergies()","de0002",JustWarning,
92 "corresponding map element not found, energy deposited locally");
93 return 0;
94 }
95 const G4DataVector* dataSet = &(*shellId).second;
96 return dataSet;
97}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59

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 139 of file G4AugerTransition.cc.

140{
141 const G4DataVector* energies = AugerTransitionEnergies(startShellId);
142 G4double energy = 0;
143 if (index < (G4int) energies->size()) {
144 G4DataVector::const_iterator pos = energies->begin();
145 energy = *(pos+index);
146 }
147 return energy;
148}
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 102 of file G4AugerTransition.cc.

103{
104 auto shellId = augerTransitionProbabilitiesMap.find(startShellId);
105 if (shellId == augerTransitionProbabilitiesMap.end() )
106 {
107
108 G4Exception("G4AugerTransition::AugerTransitionProbabilities()","de0002",
109 JustWarning,"corresponding map element not found, energy deposited locally");
110 return 0;
111 }
112
113 const G4DataVector* dataSet = &(*shellId).second;
114 return dataSet;
115}

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 153 of file G4AugerTransition.cc.

154{
155
156 const G4DataVector *probabilities = AugerTransitionProbabilities(startShellId);
157 G4DataVector::const_iterator pos = probabilities->begin();
158
159 G4double probability = 0;
160 probability = *(pos+index);
161
162 return probability;
163}
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 118 of file G4AugerTransition.cc.

119{
120 return finalShellId;
121}

◆ TransitionOriginatingShellId()

G4int G4AugerTransition::TransitionOriginatingShellId ( G4int  index) const

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

Definition at line 166 of file G4AugerTransition.cc.

167{
168 return transitionOriginatingShellIds[index];
169}

◆ 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 76 of file G4AugerTransition.cc.

77{
78 const std::vector<G4int>* dataSet = &transitionOriginatingShellIds;
79 return dataSet;
80}

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