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

#include <G4teoCrossSection.hh>

+ Inheritance diagram for G4teoCrossSection:

Public Member Functions

 G4teoCrossSection (const G4String &name)
 
virtual ~G4teoCrossSection ()
 
std::vector< G4doubleGetCrossSection (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=0)
 
G4double CrossSection (G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)
 
std::vector< G4doubleProbabilities (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=0)
 
void SetTotalCS (G4double)
 
- Public Member Functions inherited from G4VhShellCrossSection
 G4VhShellCrossSection (const G4String &xname="")
 
virtual ~G4VhShellCrossSection ()
 
G4int SelectRandomShell (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy, const G4Material *mat)
 
virtual std::vector< G4doubleGetCrossSection (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy, const G4Material *mat)=0
 
virtual G4double CrossSection (G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)=0
 
virtual std::vector< G4doubleProbabilities (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy, const G4Material *mat)=0
 
virtual void SetTotalCS (G4double)
 
const G4StringGetName () const
 

Detailed Description

Definition at line 55 of file G4teoCrossSection.hh.

Constructor & Destructor Documentation

◆ G4teoCrossSection()

G4teoCrossSection::G4teoCrossSection ( const G4String name)

Definition at line 49 of file G4teoCrossSection.cc.

50 :G4VhShellCrossSection(nam),totalCS(0.0)
51{
52 ecpssrShellMi = nullptr;
53 if (nam == "ECPSSR_Analytical")
54 {
55 ecpssrShellK = new G4ecpssrBaseKxsModel();
56 ecpssrShellLi = new G4ecpssrBaseLixsModel();
57 }
58 else if (nam == "ECPSSR_FormFactor")
59 {
60 ecpssrShellK = new G4ecpssrFormFactorKxsModel();
61 ecpssrShellLi = new G4ecpssrFormFactorLixsModel();
62 ecpssrShellMi = new G4ecpssrFormFactorMixsModel();
63 }
64 else
65 {
66 G4cout << "G4teoCrossSection::G4teoCrossSection: ERROR "
67 << " in cross section name ECPSSR_Analytical is used"
68 << G4endl;
69 ecpssrShellK = new G4ecpssrBaseKxsModel();
70 ecpssrShellLi = new G4ecpssrBaseLixsModel();
71 }
72}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ ~G4teoCrossSection()

G4teoCrossSection::~G4teoCrossSection ( )
virtual

Definition at line 74 of file G4teoCrossSection.cc.

75{
76 delete ecpssrShellK;
77 delete ecpssrShellLi;
78 delete ecpssrShellMi;
79}

Member Function Documentation

◆ CrossSection()

G4double G4teoCrossSection::CrossSection ( G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  incidentEnergy,
G4double  mass,
const G4Material mat 
)
virtual

Implements G4VhShellCrossSection.

Definition at line 109 of file G4teoCrossSection.cc.

113{
114 G4double res = 0.0;
115 if(shell > 3 && !ecpssrShellMi) {
116 return res;
117 }
118
119 else if(shell > 8) {
120 return res;
121 }
122
123 else if(fKShell == shell)
124 {
125 res = ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy);
126 }
127
128 else if(fL1Shell == shell)
129 {
130 res = ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy);
131 }
132
133 else if(fL2Shell == shell)
134 {
135 res = ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy);
136 }
137
138 else if(fL3Shell == shell)
139 {
140 res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
141 }
142
143 else if(fM1Shell == shell)
144 {
145 res = ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy);
146 }
147
148 else if(fM2Shell == shell)
149 {
150 res = ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy);
151 }
152
153 else if(fM3Shell == shell)
154 {
155 res = ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy);
156 }
157
158 else if(fM4Shell == shell)
159 {
160 res = ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy);
161 }
162
163 else if(fM5Shell == shell)
164 {
165 res = ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy);
166 }
167 return res;
168}
double G4double
Definition: G4Types.hh:83
virtual G4double CalculateCrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateL2CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateL3CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateL1CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateM5CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateM1CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateM2CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateM4CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0
virtual G4double CalculateM3CrossSection(G4int zTarget, G4double massIncident, G4double energyIncident)=0

◆ GetCrossSection()

std::vector< G4double > G4teoCrossSection::GetCrossSection ( G4int  Z,
G4double  incidentEnergy,
G4double  mass,
G4double  deltaEnergy = 0,
const G4Material mat = 0 
)
virtual

Implements G4VhShellCrossSection.

Definition at line 81 of file G4teoCrossSection.cc.

86{
87 std::vector<G4double> crossSections;
88
89 crossSections.push_back( ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy) );
90
91 crossSections.push_back( ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy) );
92 crossSections.push_back( ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy) );
93 crossSections.push_back( ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy) );
94
95 if (ecpssrShellMi) {
96
97 crossSections.push_back( ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy) );
98 crossSections.push_back( ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy) );
99 crossSections.push_back( ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy) );
100 crossSections.push_back( ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy) );
101 crossSections.push_back( ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy) );
102
103 }
104
105
106 return crossSections;
107}

Referenced by Probabilities().

◆ Probabilities()

std::vector< G4double > G4teoCrossSection::Probabilities ( G4int  Z,
G4double  incidentEnergy,
G4double  mass,
G4double  deltaEnergy = 0,
const G4Material mat = 0 
)
virtual

Implements G4VhShellCrossSection.

Definition at line 170 of file G4teoCrossSection.cc.

175{
176 std::vector<G4double> crossSections =
177 GetCrossSection(Z, incidentEnergy, mass, deltaEnergy);
178
179 for (size_t i=0; i<crossSections.size(); i++ ) {
180
181 if (totalCS) {
182 crossSections[i] = crossSections[i]/totalCS;
183 }
184
185 }
186 return crossSections;
187}
std::vector< G4double > GetCrossSection(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=0)

◆ SetTotalCS()

void G4teoCrossSection::SetTotalCS ( G4double  val)
virtual

Reimplemented from G4VhShellCrossSection.

Definition at line 189 of file G4teoCrossSection.cc.

189 {
190
191 totalCS = val;
192 // G4cout << "totalXS set to: " << val / barn << " barns" << G4endl;
193}

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