Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GammaNuclearXS Class Referencefinal

#include <G4GammaNuclearXS.hh>

+ Inheritance diagram for G4GammaNuclearXS:

Public Member Functions

 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 G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
G4double IsoCrossSection (const G4double ekin, const G4int Z, const G4int A)
 
G4double ElementCrossSection (const G4double ekin, const G4int Z)
 
G4double LowEnergyCrossSection (G4double ekin, G4int Z)
 
void CrossSectionDescription (std::ostream &) const final
 
G4GammaNuclearXSoperator= (const G4GammaNuclearXS &right)=delete
 
 G4GammaNuclearXS (const G4GammaNuclearXS &)=delete
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 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 G4StringGetName () const
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel {0}
 
G4String name
 

Detailed Description

Definition at line 60 of file G4GammaNuclearXS.hh.

Constructor & Destructor Documentation

◆ G4GammaNuclearXS() [1/2]

G4GammaNuclearXS::G4GammaNuclearXS ( )

Definition at line 78 of file G4GammaNuclearXS.cc.

80{
81 verboseLevel = 0;
82 if (verboseLevel > 0) {
83 G4cout << "G4GammaNuclearXS::G4GammaNuclearXS Initialise for Z < "
84 << MAXZGAMMAXS << G4endl;
85 }
86 ggXsection = dynamic_cast<G4PhotoNuclearCrossSection*>
88 if (ggXsection == nullptr) {
89 ggXsection = new G4PhotoNuclearCrossSection();
90 }
92
93 // full data set is uploaded once
94 if (nullptr == data) {
95 data = new G4ElementData(MAXZGAMMAXS);
96 data->SetName("gNuclear");
97 for (G4int Z=1; Z<MAXZGAMMAXS; ++Z) {
98 Initialise(Z);
99 }
100 }
101}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=false)
static G4CrossSectionDataSetRegistry * Instance()
static const char * Default_Name()
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

Referenced by G4GammaNuclearXS(), and operator=().

◆ ~G4GammaNuclearXS()

G4GammaNuclearXS::~G4GammaNuclearXS ( )
overridedefault

◆ G4GammaNuclearXS() [2/2]

G4GammaNuclearXS::G4GammaNuclearXS ( const G4GammaNuclearXS & )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4GammaNuclearXS::BuildPhysicsTable ( const G4ParticleDefinition & p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 255 of file G4GammaNuclearXS.cc.

256{
257 if(verboseLevel > 1) {
258 G4cout << "G4GammaNuclearXS::BuildPhysicsTable for "
259 << p.GetParticleName() << G4endl;
260 }
261 if(p.GetParticleName() != "gamma") {
263 ed << p.GetParticleName() << " is a wrong particle type -"
264 << " only gamma is allowed";
265 G4Exception("G4GammaNuclearXS::BuildPhysicsTable(..)","had012",
266 FatalException, ed, "");
267 return;
268 }
269
270 // prepare isotope selection
272 std::size_t nIso = temp.size();
273 for (auto const & elm : *table ) {
274 std::size_t n = elm->GetNumberOfIsotopes();
275 if (n > nIso) { nIso = n; }
276 }
277 temp.resize(nIso, 0.0);
278}
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
const G4String & GetParticleName() const

Referenced by G4InterfaceToXS::G4InterfaceToXS().

◆ CrossSectionDescription()

void G4GammaNuclearXS::CrossSectionDescription ( std::ostream & outFile) const
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 103 of file G4GammaNuclearXS.cc.

104{
105 outFile << "G4GammaNuclearXS calculates the gamma nuclear\n"
106 << "cross-section for GDR energy region on nuclei using "
107 << "data from the high precision\n"
108 << "IAEA photonuclear database (2019). Then liniear connection\n"
109 << "implemented with previous CHIPS photonuclear model." << G4endl;
110}

◆ Default_Name()

static const char * G4GammaNuclearXS::Default_Name ( )
inlinestatic

Definition at line 68 of file G4GammaNuclearXS.hh.

68{ return "GammaNuclearXS"; }

Referenced by G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), and G4GammaNuclearXS().

◆ ElementCrossSection()

G4double G4GammaNuclearXS::ElementCrossSection ( const G4double ekin,
const G4int Z )

Definition at line 134 of file G4GammaNuclearXS.cc.

135{
136 // check cache
137 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
138 if(Z == fZ && ekin == fEkin) { return fXS; }
139 fZ = Z;
140 fEkin = ekin;
141
142 auto pv = data->GetElementData(Z);
143 const G4double limCHIPS1 = 25*CLHEP::MeV;
144 const G4double limCHIPS2 = 16*CLHEP::MeV;
145 if (pv == nullptr || 1 == Z || Z == 40 || Z == 74 ||
146 (Z == 24 && ekin >= limCHIPS1) ||
147 (Z == 39 && ekin >= limCHIPS1) ||
148 (Z == 50 && ekin >= limCHIPS2) ||
149 (Z == 64 && ekin >= limCHIPS2)
150 ) {
151 fXS = ggXsection->ComputeElementXSection(ekin, Z);
152 return fXS;
153 }
154 const G4double emax = pv->GetMaxEnergy();
155
156 // low energy based on data
157 if(ekin <= emax) {
158 fXS = pv->Value(ekin);
159 // high energy CHIPS parameterisation
160 } else if(ekin >= eTransitionBound) {
161 fXS = ggXsection->ComputeElementXSection(ekin, Z);
162 // linear interpolation
163 } else {
164 const G4double rxs = xs150[Z];
165 const G4double lxs = pv->Value(emax);
166 fXS = lxs + (ekin - emax)*(rxs - lxs)/(eTransitionBound - emax);
167 }
168
169#ifdef G4VERBOSE
170 if(verboseLevel > 1) {
171 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
172 << ", nElmXS(b)= " << fXS/CLHEP::barn
173 << G4endl;
174 }
175#endif
176 return fXS;
177}
double G4double
Definition G4Types.hh:83

Referenced by GetElementCrossSection(), and IsoCrossSection().

◆ GetElementCrossSection()

G4double G4GammaNuclearXS::GetElementCrossSection ( const G4DynamicParticle * aParticle,
G4int Z,
const G4Material * mat = nullptr )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 127 of file G4GammaNuclearXS.cc.

129{
130 return ElementCrossSection(aParticle->GetKineticEnergy(), Z);
131}
G4double GetKineticEnergy() const
G4double ElementCrossSection(const G4double ekin, const G4int Z)

◆ GetIsoCrossSection()

G4double G4GammaNuclearXS::GetIsoCrossSection ( const G4DynamicParticle * aParticle,
G4int Z,
G4int A,
const G4Isotope * iso = nullptr,
const G4Element * elm = nullptr,
const G4Material * mat = nullptr )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 186 of file G4GammaNuclearXS.cc.

190{
191 return IsoCrossSection(aParticle->GetKineticEnergy(), Z, A);
192}
const G4double A[17]
G4double IsoCrossSection(const G4double ekin, const G4int Z, const G4int A)

◆ IsElementApplicable()

G4bool G4GammaNuclearXS::IsElementApplicable ( const G4DynamicParticle * ,
G4int Z,
const G4Material *  )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 113 of file G4GammaNuclearXS.cc.

115{
116 return true;
117}

◆ IsIsoApplicable()

G4bool G4GammaNuclearXS::IsIsoApplicable ( const G4DynamicParticle * ,
G4int Z,
G4int A,
const G4Element * ,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4GammaNuclearXS.cc.

122{
123 return true;
124}

◆ IsoCrossSection()

G4double G4GammaNuclearXS::IsoCrossSection ( const G4double ekin,
const G4int Z,
const G4int A )

Definition at line 195 of file G4GammaNuclearXS.cc.

196{
197 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
198 // cross section per element
199 G4double xs = ElementCrossSection(ekin, Z);
200
201 if (Z > 2) {
202 xs *= A/aeff[Z];
203 } else {
204 G4int AA = A - amin[Z];
205 if(ekin >= ehigh && AA >=0 && AA <=2) {
206 xs *= coeff[Z][AA];
207 } else {
208 xs = ggXsection->ComputeIsoXSection(ekin, Z, A);
209 }
210 }
211
212#ifdef G4VERBOSE
213 if(verboseLevel > 1) {
214 G4cout << "G4GammaNuclearXS::IsoXS: Z= " << Z << " A= " << A
215 << " Ekin(MeV)= " << ekin/CLHEP::MeV
216 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
217 }
218#endif
219 return xs;
220}

Referenced by GetIsoCrossSection(), and SelectIsotope().

◆ LowEnergyCrossSection()

G4double G4GammaNuclearXS::LowEnergyCrossSection ( G4double ekin,
G4int Z )

Definition at line 179 of file G4GammaNuclearXS.cc.

180{
181 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
182 auto pv = data->GetElementData(Z);
183 return pv->Value(ekin);
184}

◆ operator=()

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

◆ SelectIsotope()

const G4Isotope * G4GammaNuclearXS::SelectIsotope ( const G4Element * anElement,
G4double kinEnergy,
G4double logE )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 222 of file G4GammaNuclearXS.cc.

224{
225 std::size_t nIso = anElement->GetNumberOfIsotopes();
226 const G4Isotope* iso = anElement->GetIsotope(0);
227
228 if(1 == nIso) { return iso; }
229
230 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
231 G4double sum = 0.0;
232 G4int Z = anElement->GetZasInt();
233
234 // use isotope cross sections
235 std::size_t nn = temp.size();
236 if(nn < nIso) { temp.resize(nIso, 0.); }
237
238 for (std::size_t j=0; j<nIso; ++j) {
239 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
240 // << " abund= " << abundVector[j] << G4endl;
241 sum += abundVector[j]*
242 IsoCrossSection(kinEnergy, Z, anElement->GetIsotope((G4int)j)->GetN());
243 temp[j] = sum;
244 }
245 sum *= G4UniformRand();
246 for (std::size_t j = 0; j<nIso; ++j) {
247 if(temp[j] >= sum) {
248 iso = anElement->GetIsotope((G4int)j);
249 break;
250 }
251 }
252 return iso;
253}
#define G4UniformRand()
Definition Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
G4int GetZasInt() const
Definition G4Element.hh:120
G4int GetN() const
Definition G4Isotope.hh:83

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