Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4teoCrossSection.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28//
29// History:
30// -----------
31// 21 Apr 2009 ALF 1st implementation
32// 29 Apr 2009 ALF Updated Desing for Integration
33// 15 Mar 2011 ALF introduced the usage of G4AtomicShellEnumerator
34// 20 Oct 2011 ALF updated to take into account ECPSSR Form Factor
35// 09 Mar 2012 LP update methods
36// 09 Mar 2012 ALF update for M-shells Simulation
37// 10 Nov 2021 S. Guatelli & S. Bakr Added ECPSSR form factor documented
38// in Bakr et al, NIM B, vol. 436, pp:285-291, 2018 and
39// called here ECPSSR_ANSTO
40
41#include "globals.hh"
42#include "G4teoCrossSection.hh"
43#include "G4Proton.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
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}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
88{
89 delete ecpssrShellK;
90 delete ecpssrShellLi;
91 delete ecpssrShellMi;
92}
93
94//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
95
97 G4double incidentEnergy,
98 G4double mass,
100 const G4Material*)
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}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121
123 G4double incidentEnergy,
124 G4double mass,
125 const G4Material*)
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}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174
176 G4double incidentEnergy,
177 G4double mass,
178 G4double deltaEnergy,
179 const G4Material*)
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}
191
192//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193
195 totalCS = val;
196}
197
198
199
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
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
void SetTotalCS(G4double) override
G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat) override
std::vector< G4double > GetCrossSection(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=nullptr) override
std::vector< G4double > Probabilities(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=nullptr) override
G4teoCrossSection(const G4String &name)