Geant4 10.7.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
 
const G4DataVectorAugerTransitionEnergies (G4int startShellId) const
 
const G4DataVectorAugerTransitionProbabilities (G4int startShellId) const
 
G4int FinalShellId () const
 
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
 

Detailed Description

Definition at line 54 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 
)

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
55
56}

◆ ~G4AugerTransition()

G4AugerTransition::~G4AugerTransition ( )

Definition at line 58 of file G4AugerTransition.cc.

59{
60
61}

Member Function Documentation

◆ AugerOriginatingShellId()

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

Definition at line 147 of file G4AugerTransition.cc.

148{
149 const std::vector<G4int>* ids = AugerOriginatingShellIds(startShellId);
150 // G4int i =
151 std::vector<G4int>::const_iterator pos = ids->begin();
152 G4int n = 0;
153 n = *(pos+index);
154 return n;
155}
int G4int
Definition: G4Types.hh:85
const std::vector< G4int > * AugerOriginatingShellIds(G4int startShellId) const

◆ AugerOriginatingShellIds()

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

Definition at line 66 of file G4AugerTransition.cc.

67{
68 std::map<G4int,std::vector<G4int>,std::less<G4int> >::const_iterator shellId = augerOriginatingShellIdsMap.find(startShellId);
69
70 const std::vector<G4int>* dataSet = &(*shellId).second;
71 //const std::vector<G4int>* dataOut = 0;
72
73 if (dataSet->size() == 0) {G4cout << "Error: no auger Id found"<< G4endl;}
74 else {
75
76 // dataOut = &dataSet;
77
78 }
79
80 return dataSet;
81}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

Referenced by AugerOriginatingShellId().

◆ AugerTransitionEnergies()

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

Definition at line 95 of file G4AugerTransition.cc.

96{
97 std::map<G4int,G4DataVector,std::less<G4int> >::const_iterator shellId = augerTransitionEnergiesMap.find(startShellId);
98
99 if (shellId == augerTransitionEnergiesMap.end() )
100 {
101 G4Exception("G4AugerTransition::AugerTransitionEnergies()","de0002",JustWarning,"corresponding map element not found, energy deposited locally");
102 return 0;
103 }
104
105 const G4DataVector* dataSet = &(*shellId).second;
106
107
108 return dataSet;
109}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

Referenced by AugerTransitionEnergy().

◆ AugerTransitionEnergy()

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

Definition at line 160 of file G4AugerTransition.cc.

161{
162 const G4DataVector* energies = AugerTransitionEnergies(startShellId);
163 G4double energy = 0;
164 if (index < (G4int) energies->size()) {
165 G4DataVector::const_iterator pos = energies->begin();
166 energy = *(pos+index);
167 }
168 return energy;
169}
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

Definition at line 114 of file G4AugerTransition.cc.

115{
116
117 //debugging
118 //if (startShellId == 1){G4cout <<"OI!!!"<< G4endl;}
119
120 std::map<G4int,G4DataVector,std::less<G4int> >::const_iterator shellId = augerTransitionProbabilitiesMap.find(startShellId);
121
122 if (shellId == augerTransitionProbabilitiesMap.end() )
123 {
124
125 G4Exception("G4AugerTransition::AugerTransitionProbabilities()","de0002",JustWarning,"corresponding map element not found, energy deposited locally");
126 return 0;
127 }
128
129 const G4DataVector* dataSet = &(*shellId).second;
130 // debugging purpose:
131 /* G4cout << "id: " << shellId->first << G4endl;
132 G4cout << "size:" << dataSet->size() << G4endl;
133 for (G4int i = 0; i < dataSet->size(); i++){
134 G4cout << (dataSet[0])[i] << G4endl;
135 }*/
136 return dataSet;
137}

Referenced by AugerTransitionProbability().

◆ AugerTransitionProbability()

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

Definition at line 174 of file G4AugerTransition.cc.

175{
176
177 const G4DataVector *probabilities = AugerTransitionProbabilities(startShellId);
178 G4DataVector::const_iterator pos = probabilities->begin();
179
180 G4double probability = 0;
181 probability = *(pos+index);
182
183 return probability;
184
185}
const G4DataVector * AugerTransitionProbabilities(G4int startShellId) const

◆ FinalShellId()

G4int G4AugerTransition::FinalShellId ( ) const

Definition at line 139 of file G4AugerTransition.cc.

140{
141 return finalShellId;
142}

◆ TransitionOriginatingShellId()

G4int G4AugerTransition::TransitionOriginatingShellId ( G4int  index) const

Definition at line 187 of file G4AugerTransition.cc.

188{
189 return transitionOriginatingShellIds[index];
190}

◆ TransitionOriginatingShellIds()

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

Definition at line 85 of file G4AugerTransition.cc.

86{
87
88 const std::vector<G4int>* dataSet = &transitionOriginatingShellIds;
89 return dataSet;
90}

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