Geant4 11.2.2
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=nullptr) override
 
G4double CrossSection (G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat) override
 
std::vector< G4doubleProbabilities (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=nullptr) override
 
void SetTotalCS (G4double) override
 
 G4teoCrossSection (const G4teoCrossSection &)=delete
 
G4teoCrossSectionoperator= (const G4teoCrossSection &right)=delete
 
- 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)
 
const G4StringGetName () const
 
 G4VhShellCrossSection (const G4VhShellCrossSection &)=delete
 
G4VhShellCrossSectionoperator= (const G4VhShellCrossSection &right)=delete
 

Detailed Description

Definition at line 54 of file G4teoCrossSection.hh.

Constructor & Destructor Documentation

◆ G4teoCrossSection() [1/2]

G4teoCrossSection::G4teoCrossSection ( const G4String & name)
explicit

Definition at line 54 of file G4teoCrossSection.cc.

55 :G4VhShellCrossSection(nam),totalCS(0.0)
56{
57 ecpssrShellMi = nullptr;
58 if (nam == "ECPSSR_Analytical")
59 {
60 ecpssrShellK = new G4ecpssrBaseKxsModel();
61 ecpssrShellLi = new G4ecpssrBaseLixsModel();
62 }
63 else if (nam == "ECPSSR_FormFactor")
64 {
65 ecpssrShellK = new G4ecpssrFormFactorKxsModel();
66 ecpssrShellLi = new G4ecpssrFormFactorLixsModel();
67 ecpssrShellMi = new G4ecpssrFormFactorMixsModel();
68 }
69 else if (nam == "ECPSSR_ANSTO")
70 {
71 ecpssrShellK = new G4ANSTOecpssrKxsModel();
72 ecpssrShellLi = new G4ANSTOecpssrLixsModel();
73 ecpssrShellMi = new G4ANSTOecpssrMixsModel();
74 }
75 else
76 {
77 G4cout << "G4teoCrossSection::G4teoCrossSection: ERROR "
78 << " in cross section name ECPSSR_Analytical is used"
79 << G4endl;
80 ecpssrShellK = new G4ecpssrBaseKxsModel();
81 ecpssrShellLi = new G4ecpssrBaseLixsModel();
82 }
83}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VhShellCrossSection(const G4String &xname="")

◆ ~G4teoCrossSection()

G4teoCrossSection::~G4teoCrossSection ( )
virtual

Definition at line 87 of file G4teoCrossSection.cc.

88{
89 delete ecpssrShellK;
90 delete ecpssrShellLi;
91 delete ecpssrShellMi;
92}

◆ G4teoCrossSection() [2/2]

G4teoCrossSection::G4teoCrossSection ( const G4teoCrossSection & )
delete

Member Function Documentation

◆ CrossSection()

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

Implements G4VhShellCrossSection.

Definition at line 122 of file G4teoCrossSection.cc.

126{
127 G4double res = 0.0;
128 if(shell > 3 && !ecpssrShellMi) {
129 return res;
130 }
131 else if(shell > 8) {
132 return res;
133 }
134 else if(fKShell == shell)
135 {
136 res = ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy);
137 }
138 else if(fL1Shell == shell)
139 {
140 res = ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy);
141 }
142 else if(fL2Shell == shell)
143 {
144 res = ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy);
145 }
146 else if(fL3Shell == shell)
147 {
148 res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
149 }
150 else if(fM1Shell == shell)
151 {
152 res = ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy);
153 }
154 else if(fM2Shell == shell)
155 {
156 res = ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy);
157 }
158 else if(fM3Shell == shell)
159 {
160 res = ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy);
161 }
162 else if(fM4Shell == shell)
163 {
164 res = ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy);
165 }
166 else if(fM5Shell == shell)
167 {
168 res = ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy);
169 }
170 return res;
171}
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 = nullptr )
overridevirtual

Implements G4VhShellCrossSection.

Definition at line 96 of file G4teoCrossSection.cc.

101{
102 std::vector<G4double> crossSections;
103
104 crossSections.push_back( ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy) );
105
106 crossSections.push_back( ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy) );
107 crossSections.push_back( ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy) );
108 crossSections.push_back( ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy) );
109
110 if (ecpssrShellMi) {
111 crossSections.push_back( ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy) );
112 crossSections.push_back( ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy) );
113 crossSections.push_back( ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy) );
114 crossSections.push_back( ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy) );
115 crossSections.push_back( ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy) );
116 }
117 return crossSections;
118}

Referenced by Probabilities().

◆ operator=()

G4teoCrossSection & G4teoCrossSection::operator= ( const G4teoCrossSection & right)
delete

◆ Probabilities()

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

Implements G4VhShellCrossSection.

Definition at line 175 of file G4teoCrossSection.cc.

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

◆ SetTotalCS()

void G4teoCrossSection::SetTotalCS ( G4double val)
overridevirtual

Reimplemented from G4VhShellCrossSection.

Definition at line 194 of file G4teoCrossSection.cc.

194 {
195 totalCS = val;
196}

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