Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaNuclearXS.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// GEANT4 Class file
29//
30//
31// File name: G4GammaNuclearXS
32//
33// Authors V.Ivantchenko, Geant4, 20 October 2020
34// B.Kutsenko, BINP/NSU, 10 August 2021
35//
36// Modifications:
37//
38
39#include "G4GammaNuclearXS.hh"
40#include "G4DynamicParticle.hh"
41#include "G4ThreeVector.hh"
42#include "G4ElementTable.hh"
43#include "G4Material.hh"
44#include "G4Element.hh"
45#include "G4ElementData.hh"
46#include "G4PhysicsVector.hh"
52#include "Randomize.hh"
53#include "G4SystemOfUnits.hh"
54#include "G4Gamma.hh"
55#include "G4IsotopeList.hh"
56#include "G4AutoLock.hh"
57
58#include <fstream>
59#include <sstream>
60#include <vector>
61
62G4ElementData* G4GammaNuclearXS::data = nullptr;
63
64G4double G4GammaNuclearXS::coeff[3][3];
65G4double G4GammaNuclearXS::xs150[] = {0.0};
66const G4double G4GammaNuclearXS::eTransitionBound = 150.*CLHEP::MeV;
67const G4int G4GammaNuclearXS::freeVectorException[] = {
684, 6, 7, 8, 27, 39, 45, 65, 67, 69, 73};
69G4String G4GammaNuclearXS::gDataDirectory = "";
70
72 : G4VCrossSectionDataSet(Default_Name()), gamma(G4Gamma::Gamma())
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}
94
95void G4GammaNuclearXS::CrossSectionDescription(std::ostream& outFile) const
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}
103
104G4bool
106 G4int, const G4Material*)
107{
108 return true;
109}
110
112 G4int, G4int,
113 const G4Element*, const G4Material*)
114{
115 return true;
116}
117
120 G4int ZZ, const G4Material* mat)
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}
158
160{
161 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
162 return GetElementCrossSection(&theGamma, ZZ);
163}
164
166{
167 const G4int Z = (ZZ < MAXZGAMMAXS) ? ZZ : MAXZGAMMAXS - 1;
168 auto pv = data->GetElementData(Z);
169 return pv->Value(ekin);
170}
171
174{
175 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
176 return GetIsoCrossSection(&theGamma, Z, A);
177}
178
180 const G4DynamicParticle* aParticle,
181 G4int ZZ, G4int A,
182 const G4Isotope*, const G4Element*, const G4Material* mat)
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}
209
211 const G4Element* anElement, G4double kinEnergy, G4double)
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}
242
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}
267
268const G4String& G4GammaNuclearXS::FindDirectoryPath()
269{
270 // build the complete string identifying the file with the data set
271 if(gDataDirectory.empty()) {
272 std::ostringstream ost;
273 ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/gamma/inel";
274 gDataDirectory = ost.str();
275 }
276 return gDataDirectory;
277}
278
279void G4GammaNuclearXS::Initialise(G4int Z)
280{
281 // upload data from file
282 std::ostringstream ost;
283 ost << FindDirectoryPath() << Z ;
284 G4PhysicsVector* v = RetrieveVector(ost, true, Z);
285
286 data->InitialiseForElement(Z, v);
287 /*
288 G4cout << "G4GammaNuclearXS::Initialise for Z= " << Z
289 << " A= " << Amean << " Amin= " << amin[Z]
290 << " Amax= " << amax[Z] << G4endl;
291 */
292 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), eTransitionBound);
293 xs150[Z] = ggXsection->GetElementCrossSection(&theGamma, Z, 0);
294
295 // compute corrections for low Z data
296 if(Z <= 2){
297 theGamma.SetKineticEnergy(10*CLHEP::GeV);
298 if(amax[Z] > amin[Z]) {
299 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
300 G4int AA = A - amin[Z];
301 if(AA >= 0 && AA <= 2) {
302 G4double sig1 = ggXsection->GetIsoCrossSection(&theGamma, Z, A);
303 G4double sig2 = ggXsection->GetElementCrossSection(&theGamma, Z, 0);
304 if(sig2 > 0.) { coeff[Z][AA] = (sig1/sig2); }
305 else { coeff[Z][AA] = 1.; }
306 }
307 }
308 }
309 }
310}
311
313G4GammaNuclearXS::RetrieveVector(std::ostringstream& ost, G4bool warn, G4int Z)
314{
315 G4PhysicsVector* v = nullptr;
316
317 std::ifstream filein(ost.str().c_str());
318 if (!filein.is_open()) {
319 if(warn) {
321 ed << "Data file <" << ost.str().c_str()
322 << "> is not opened!";
323 G4Exception("G4GammaNuclearXS::RetrieveVector(..)","had014",
324 FatalException, ed, "Check G4PARTICLEXSDATA");
325 }
326 } else {
327 if(verboseLevel > 1) {
328 G4cout << "File " << ost.str()
329 << " is opened by G4GammaNuclearXS" << G4endl;
330 }
331 // retrieve data from DB
332 if(std::find(std::begin(freeVectorException), std::end(freeVectorException), Z)
333 == std::end(freeVectorException)) {
334 v = new G4PhysicsLinearVector(false);
335 } else {
336 v = new G4PhysicsFreeVector(false);
337 }
338 if(!v->Retrieve(filein, true)) {
340 ed << "Data file <" << ost.str().c_str()
341 << "> is not retrieved!";
342 G4Exception("G4GammaNuclearXS::RetrieveVector(..)","had015",
343 FatalException, ed, "Check G4PARTICLEXSDATA");
344 }
345 }
346 return v;
347}
348
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=false)
static G4CrossSectionDataSetRegistry * Instance()
G4double GetKineticEnergy() const
G4PhysicsVector * GetElementData(G4int Z) const
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
void SetName(const G4String &nam)
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
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 GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double LowEnergyCrossSection(G4double ekin, G4int Z)
G4double ElementCrossSection(G4double ekin, G4int Z)
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
void CrossSectionDescription(std::ostream &) const final
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
static G4HadronicParameters * Instance()
const G4String & GetDirPARTICLEXS() const
G4int GetN() const
Definition G4Isotope.hh:83
const G4String & GetParticleName() const
G4double GetMaxEnergy() const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
void SetForAllAtomsAndEnergies(G4bool val)