Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 (G4double ekin, G4int Z, G4int A)
 
G4double ElementCrossSection (G4double ekin, 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
 
G4String name
 

Detailed Description

Definition at line 59 of file G4GammaNuclearXS.hh.

Constructor & Destructor Documentation

◆ G4GammaNuclearXS() [1/2]

G4GammaNuclearXS::G4GammaNuclearXS ( )

Definition at line 71 of file G4GammaNuclearXS.cc.

73{
74 verboseLevel = 0;
75 if (verboseLevel > 0) {
76 G4cout << "G4GammaNuclearXS::G4GammaNuclearXS Initialise for Z < "
77 << MAXZGAMMAXS << G4endl;
78 }
79 ggXsection =
81 if (ggXsection == nullptr)
82 ggXsection = new G4PhotoNuclearCrossSection();
84
85 // full data set is uploaded once
86 if (nullptr == data) {
87 data = new G4ElementData(MAXZGAMMAXS);
88 data->SetName("gNuclear");
89 for (G4int Z=1; Z<MAXZGAMMAXS; ++Z) {
90 Initialise(Z);
91 }
92 }
93}
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()
void SetName(const G4String &nam)
static const char * Default_Name()
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~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 243 of file G4GammaNuclearXS.cc.

244{
245 if(verboseLevel > 1) {
246 G4cout << "G4GammaNuclearXS::BuildPhysicsTable for "
247 << p.GetParticleName() << G4endl;
248 }
249 if(p.GetParticleName() != "gamma") {
251 ed << p.GetParticleName() << " is a wrong particle type -"
252 << " only gamma is allowed";
253 G4Exception("G4GammaNuclearXS::BuildPhysicsTable(..)","had012",
254 FatalException, ed, "");
255 return;
256 }
257
258 // prepare isotope selection
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
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
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

Definition at line 67 of file G4GammaNuclearXS.hh.

67{ return "GammaNuclearXS"; }

Referenced by G4ElectroVDNuclearModel::G4ElectroVDNuclearModel().

◆ ElementCrossSection()

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

Definition at line 159 of file G4GammaNuclearXS.cc.

160{
161 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
162 return GetElementCrossSection(&theGamma, ZZ);
163}
CLHEP::Hep3Vector G4ThreeVector
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4GammaNuclearXS.cc.

121{
122 // check cache
123 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
124 const G4double ekin = aParticle->GetKineticEnergy();
125 if(Z == fZ && ekin == fEkin) { return fXS; }
126 fZ = Z;
127 fEkin = ekin;
128
129 auto pv = data->GetElementData(Z);
130 if(pv == nullptr || 1 == Z) {
131 fXS = ggXsection->GetElementCrossSection(aParticle, Z, mat);
132 return fXS;
133 }
134 const G4double emax = pv->GetMaxEnergy();
135
136 // low energy based on data
137 if(ekin <= emax) {
138 fXS = pv->Value(ekin);
139 // high energy CHIPS parameterisation
140 } else if(ekin >= eTransitionBound) {
141 fXS = ggXsection->GetElementCrossSection(aParticle, Z, mat);
142 // linear interpolation
143 } else {
144 const G4double rxs = xs150[Z];
145 const G4double lxs = pv->Value(emax);
146 fXS = lxs + (ekin - emax)*(rxs - lxs)/(eTransitionBound - emax);
147 }
148
149#ifdef G4VERBOSE
150 if(verboseLevel > 1) {
151 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
152 << ", nElmXS(b)= " << fXS/CLHEP::barn
153 << G4endl;
154 }
155#endif
156 return fXS;
157}
double G4double
Definition G4Types.hh:83
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()

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 179 of file G4GammaNuclearXS.cc.

183{
184 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
185 // cross section per element
186 G4double xs = GetElementCrossSection(aParticle, Z, mat);
187 const G4double ekin = aParticle->GetKineticEnergy();
188
189 if (Z > 2) {
190 xs *= A/aeff[Z];
191 } else {
192 G4int AA = A - amin[Z];
193 if(ekin >= 10.*CLHEP::GeV && AA >=0 && AA <=2) {
194 xs *= coeff[Z][AA];
195 } else {
196 xs = ggXsection->GetIsoCrossSection(aParticle, Z, A);
197 }
198 }
199
200#ifdef G4VERBOSE
201 if(verboseLevel > 1) {
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}
const G4double A[17]
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()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 105 of file G4GammaNuclearXS.cc.

107{
108 return true;
109}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 111 of file G4GammaNuclearXS.cc.

114{
115 return true;
116}

◆ IsoCrossSection()

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

Definition at line 173 of file G4GammaNuclearXS.cc.

174{
175 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
176 return GetIsoCrossSection(&theGamma, Z, A);
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()

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

Definition at line 165 of file G4GammaNuclearXS.cc.

166{
167 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
168 auto pv = data->GetElementData(Z);
169 return pv->Value(ekin);
170}
G4double Value(const G4double energy, std::size_t &lastidx) const

◆ 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 210 of file G4GammaNuclearXS.cc.

212{
213 std::size_t nIso = anElement->GetNumberOfIsotopes();
214 const G4Isotope* iso = anElement->GetIsotope(0);
215
216 if(1 == nIso) { return iso; }
217
218 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
219 G4double sum = 0.0;
220 G4int Z = anElement->GetZasInt();
221
222 // use isotope cross sections
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 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
228 // << " abund= " << abundVector[j] << G4endl;
229 sum += abundVector[j]*
230 IsoCrossSection(kinEnergy, Z, anElement->GetIsotope((G4int)j)->GetN());
231 temp[j] = sum;
232 }
233 sum *= G4UniformRand();
234 for (std::size_t j = 0; j<nIso; ++j) {
235 if(temp[j] >= sum) {
236 iso = anElement->GetIsotope((G4int)j);
237 break;
238 }
239 }
240 return iso;
241}
#define G4UniformRand()
Definition Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
G4int GetZasInt() const
Definition G4Element.hh:120
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
G4int GetN() const
Definition G4Isotope.hh:83

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