#include <G4GammaNuclearXS.hh>
|
| G4GammaNuclearXS () |
|
| ~G4GammaNuclearXS () override=default |
|
G4bool | IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final |
|
G4bool | IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *mat) final |
|
G4double | GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final |
|
G4double | GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final |
|
const G4Isotope * | SelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final |
|
void | BuildPhysicsTable (const G4ParticleDefinition &) final |
|
G4double | IsoCrossSection (G4double ekin, G4int Z, G4int A) |
|
G4double | ElementCrossSection (G4double ekin, G4int Z) |
|
G4double | LowEnergyCrossSection (G4double ekin, G4int Z) |
|
void | CrossSectionDescription (std::ostream &) const final |
|
G4GammaNuclearXS & | operator= (const G4GammaNuclearXS &right)=delete |
|
| G4GammaNuclearXS (const G4GammaNuclearXS &)=delete |
|
| G4VCrossSectionDataSet (const G4String &nam="") |
|
virtual | ~G4VCrossSectionDataSet () |
|
G4double | GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr) |
|
G4double | ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr) |
|
virtual G4double | ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr) |
|
virtual G4double | ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) |
|
virtual void | DumpPhysicsTable (const G4ParticleDefinition &) |
|
virtual void | SetVerboseLevel (G4int value) |
|
G4double | GetMinKinEnergy () const |
|
void | SetMinKinEnergy (G4double value) |
|
G4double | GetMaxKinEnergy () const |
|
void | SetMaxKinEnergy (G4double value) |
|
bool | ForAllAtomsAndEnergies () const |
|
void | SetForAllAtomsAndEnergies (G4bool val) |
|
const G4String & | GetName () const |
|
void | SetName (const G4String &nam) |
|
G4VCrossSectionDataSet & | operator= (const G4VCrossSectionDataSet &right)=delete |
|
| G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete |
|
Definition at line 59 of file G4GammaNuclearXS.hh.
◆ G4GammaNuclearXS() [1/2]
G4GammaNuclearXS::G4GammaNuclearXS |
( |
| ) |
|
Definition at line 71 of file G4GammaNuclearXS.cc.
73{
76 G4cout <<
"G4GammaNuclearXS::G4GammaNuclearXS Initialise for Z < "
78 }
79 ggXsection =
81 if (ggXsection == nullptr)
84
85
86 if (nullptr == data) {
89 for (
G4int Z=1; Z<MAXZGAMMAXS; ++Z) {
90 Initialise(Z);
91 }
92 }
93}
G4GLOB_DLL std::ostream G4cout
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=false)
static G4CrossSectionDataSetRegistry * Instance()
void SetName(const G4String &nam)
static const char * Default_Name()
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)
◆ ~G4GammaNuclearXS()
G4GammaNuclearXS::~G4GammaNuclearXS |
( |
| ) |
|
|
overridedefault |
◆ G4GammaNuclearXS() [2/2]
◆ BuildPhysicsTable()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 243 of file G4GammaNuclearXS.cc.
244{
246 G4cout <<
"G4GammaNuclearXS::BuildPhysicsTable for "
248 }
252 << " only gamma is allowed";
253 G4Exception(
"G4GammaNuclearXS::BuildPhysicsTable(..)",
"had012",
255 return;
256 }
257
258
260 std::size_t nIso = temp.size();
261 for ( auto & elm : *table ) {
262 std::size_t
n = elm->GetNumberOfIsotopes();
263 if (n > nIso) { nIso =
n; }
264 }
265 temp.resize(nIso, 0.0);
266}
std::vector< G4Element * > G4ElementTable
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
static G4ElementTable * GetElementTable()
const G4String & GetParticleName() const
◆ CrossSectionDescription()
void G4GammaNuclearXS::CrossSectionDescription |
( |
std::ostream & | outFile | ) |
const |
|
finalvirtual |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 95 of file G4GammaNuclearXS.cc.
96{
97 outFile << "G4GammaNuclearXS calculates the gamma nuclear\n"
98 << "cross-section for GDR energy region on nuclei using "
99 << "data from the high precision\n"
100 << "IAEA photonuclear database (2019). Then liniear connection\n"
101 <<
"implemented with previous CHIPS photonuclear model." <<
G4endl;
102}
◆ Default_Name()
static const char * G4GammaNuclearXS::Default_Name |
( |
| ) |
|
|
inlinestatic |
◆ ElementCrossSection()
Definition at line 159 of file G4GammaNuclearXS.cc.
160{
163}
CLHEP::Hep3Vector G4ThreeVector
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
◆ GetElementCrossSection()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 119 of file G4GammaNuclearXS.cc.
121{
122
123 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
125 if(Z == fZ && ekin == fEkin) { return fXS; }
126 fZ = Z;
127 fEkin = ekin;
128
130 if(pv == nullptr || 1 == Z) {
132 return fXS;
133 }
135
136
137 if(ekin <= emax) {
138 fXS = pv->Value(ekin);
139
140 } else if(ekin >= eTransitionBound) {
142
143 } else {
145 const G4double lxs = pv->Value(emax);
146 fXS = lxs + (ekin - emax)*(rxs - lxs)/(eTransitionBound - emax);
147 }
148
149#ifdef G4VERBOSE
151 G4cout <<
"Z= " << Z <<
" Ekin(MeV)= " << ekin/CLHEP::MeV
152 << ", nElmXS(b)= " << fXS/CLHEP::barn
154 }
155#endif
156 return fXS;
157}
G4double GetKineticEnergy() const
G4PhysicsVector * GetElementData(G4int Z) const
G4double GetMaxEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
Referenced by ElementCrossSection(), and GetIsoCrossSection().
◆ GetIsoCrossSection()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 179 of file G4GammaNuclearXS.cc.
183{
184 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
185
188
189 if (Z > 2) {
191 } else {
193 if(ekin >= 10.*CLHEP::GeV && AA >=0 && AA <=2) {
194 xs *= coeff[Z][AA];
195 } else {
197 }
198 }
199
200#ifdef G4VERBOSE
202 G4cout <<
"G4GammaNuclearXS::IsoXS: Z= " << Z <<
" A= " <<
A
203 << " Ekin(MeV)= " << ekin/CLHEP::MeV
204 <<
", ElmXS(b)= " << xs/CLHEP::barn <<
G4endl;
205 }
206#endif
207 return xs;
208}
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
Referenced by IsoCrossSection().
◆ IsElementApplicable()
◆ IsIsoApplicable()
◆ IsoCrossSection()
Definition at line 173 of file G4GammaNuclearXS.cc.
174{
177}
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
Referenced by SelectIsotope().
◆ LowEnergyCrossSection()
Definition at line 165 of file G4GammaNuclearXS.cc.
166{
167 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
169 return pv->
Value(ekin);
170}
G4double Value(const G4double energy, std::size_t &lastidx) const
◆ operator=()
◆ SelectIsotope()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 210 of file G4GammaNuclearXS.cc.
212{
215
216 if(1 == nIso) { return iso; }
217
221
222
223 std::size_t
nn = temp.size();
224 if(nn < nIso) { temp.resize(nIso, 0.); }
225
226 for (std::size_t j=0; j<nIso; ++j) {
227
228
229 sum += abundVector[j]*
231 temp[j] = sum;
232 }
234 for (std::size_t j = 0; j<nIso; ++j) {
235 if(temp[j] >= sum) {
237 break;
238 }
239 }
240 return iso;
241}
G4double * GetRelativeAbundanceVector() const
const G4Isotope * GetIsotope(G4int iso) const
size_t GetNumberOfIsotopes() const
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
The documentation for this class was generated from the following files: