Geant4 9.6.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 55 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 45 of file G4AugerTransition.cc.

49{
50 finalShellId = finalShell;
51 augerOriginatingShellIdsMap = *idMap;
52 augerTransitionEnergiesMap = *energyMap;
53 augerTransitionProbabilitiesMap = *probabilityMap;
54 transitionOriginatingShellIds = transIds;
55
56
57}

◆ ~G4AugerTransition()

G4AugerTransition::~G4AugerTransition ( )

Definition at line 59 of file G4AugerTransition.cc.

60{
61
62}

Member Function Documentation

◆ AugerOriginatingShellId()

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

Definition at line 148 of file G4AugerTransition.cc.

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

◆ AugerOriginatingShellIds()

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

Definition at line 67 of file G4AugerTransition.cc.

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

Referenced by AugerOriginatingShellId().

◆ AugerTransitionEnergies()

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

Definition at line 96 of file G4AugerTransition.cc.

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

Referenced by AugerTransitionEnergy().

◆ AugerTransitionEnergy()

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

Definition at line 161 of file G4AugerTransition.cc.

162{
163 const G4DataVector* energies = AugerTransitionEnergies(startShellId);
164 G4double energy = 0;
165 if (index < (G4int) energies->size()) {
166 G4DataVector::const_iterator pos = energies->begin();
167 energy = *(pos+index);
168 }
169 return energy;
170}
double G4double
Definition: G4Types.hh:64
const G4DataVector * AugerTransitionEnergies(G4int startShellId) const

◆ AugerTransitionProbabilities()

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

Definition at line 115 of file G4AugerTransition.cc.

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

Referenced by AugerTransitionProbability().

◆ AugerTransitionProbability()

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

Definition at line 175 of file G4AugerTransition.cc.

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

◆ FinalShellId()

G4int G4AugerTransition::FinalShellId ( ) const

Definition at line 140 of file G4AugerTransition.cc.

141{
142 return finalShellId;
143}

◆ TransitionOriginatingShellId()

G4int G4AugerTransition::TransitionOriginatingShellId ( G4int  index) const

Definition at line 188 of file G4AugerTransition.cc.

189{
190 return transitionOriginatingShellIds[index];
191}

◆ TransitionOriginatingShellIds()

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

Definition at line 86 of file G4AugerTransition.cc.

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

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