Geant4 11.1.1
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"
48#include "Randomize.hh"
49#include "G4SystemOfUnits.hh"
50#include "G4Gamma.hh"
51#include "G4IsotopeList.hh"
52
53#include <fstream>
54#include <sstream>
55#include <vector>
56
57G4ElementData* G4GammaNuclearXS::data = nullptr;
58
59G4double G4GammaNuclearXS::coeff[3][3];
60G4double G4GammaNuclearXS::xs150[MAXZGAMMAXS] = {0.0};
61G4String G4GammaNuclearXS::gDataDirectory = "";
62
63#ifdef G4MULTITHREADED
64 G4Mutex G4GammaNuclearXS::gNuclearXSMutex = G4MUTEX_INITIALIZER;
65#endif
66
68 : G4VCrossSectionDataSet(Default_Name()),
69 gamma(G4Gamma::Gamma())
70{
71 // verboseLevel = 0;
72 if(verboseLevel > 0){
73 G4cout << "G4GammaNuclearXS::G4GammaNuclearXS Initialise for Z < "
74 << MAXZGAMMAXS << G4endl;
75 }
77 if(ggXsection == nullptr) ggXsection = new G4PhotoNuclearCrossSection();
79}
80
82{
83 if(isMaster) {
84 delete data;
85 data = nullptr;
86 }
87}
88
89void G4GammaNuclearXS::CrossSectionDescription(std::ostream& outFile) const
90{
91 outFile << "G4GammaNuclearXS calculates the gamma nuclear\n"
92 << "cross-section for GDR energy region on nuclei using data from the high precision\n"
93 << "IAEA photonuclear database (2019). Then liniear connection\n"
94 <<"implemented with previous CHIPS photonuclear model\n";
95}
96
97G4bool
99 G4int, const G4Material*)
100{
101 return true;
102}
103
105 G4int, G4int,
106 const G4Element*, const G4Material*)
107{
108 return true;
109}
110
113 G4int ZZ, const G4Material* mat)
114{
115 const G4int Z = (ZZ >= MAXZGAMMAXS) ? MAXZGAMMAXS - 1 : ZZ;
116 auto pv = GetPhysicsVector(Z);
117
118 if(pv == nullptr) {
119 return ggXsection->GetElementCrossSection(aParticle, Z, mat);
120 }
121 const G4double emax = pv->GetMaxEnergy();
122 const G4double ekin = aParticle->GetKineticEnergy();
123 G4double xs = 0.0;
124 if(ekin <= emax) {
125 xs = pv->Value(ekin);
126 } else if(ekin >= rTransitionBound){
127 xs = ggXsection->GetElementCrossSection(aParticle, Z, mat);
128 } else {
129 const G4double rxs = xs150[Z];
130 const G4double lxs = pv->Value(emax);
131 xs = lxs + (ekin - emax)*(rxs - lxs)/(rTransitionBound-emax);
132 }
133
134#ifdef G4VERBOSE
135 if(verboseLevel > 1) {
136 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
137 << ", nElmXS(b)= " << xs/CLHEP::barn
138 << G4endl;
139 }
140#endif
141 return xs;
142}
143
146{
147 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
148 return GetElementCrossSection(&theGamma, ZZ);
149}
150
153{
154 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), ekin);
155 return GetIsoCrossSection(&theGamma, ZZ, A);
156}
157
159 const G4DynamicParticle* aParticle,
160 G4int ZZ, G4int A,
161 const G4Isotope*, const G4Element*,
162 const G4Material*)
163{
164 const G4int Z = (ZZ >= MAXZGAMMAXS) ? MAXZGAMMAXS - 1 : ZZ;
165 /*
166 G4cout << "IsoCrossSection Z= " << Z << " A= " << A
167 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
168 << " E(MeV)= " << ekin << G4endl;
169 */
170 auto pv = GetPhysicsVector(Z);
171 if(pv == nullptr) {
172 return ggXsection->GetIsoCrossSection(aParticle, Z, A);
173 }
174 const G4double ekin = aParticle->GetKineticEnergy();
175 const G4double emax = pv->GetMaxEnergy();
176 G4double xs = 0.0;
177
178 // compute isotope cross section if applicable
179 if(amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z] &&
180 ekin < rTransitionBound) {
181 auto pviso = data->GetComponentDataByIndex(Z, A - amin[Z]);
182 // isotope file exists
183 if(nullptr != pviso) {
184 const G4double emaxiso = pviso->GetMaxEnergy();
185 if(ekin <= emaxiso) {
186 xs = pviso->Value(ekin);
187 } else {
189 theGamma(gamma, G4ThreeVector(0,0,1.), rTransitionBound);
190 const G4double rxs = ggXsection->GetIsoCrossSection(&theGamma, Z, A);
191 const G4double lxs = pviso->Value(emaxiso);
192 xs = lxs + (ekin - emaxiso)*(rxs - lxs)/(rTransitionBound-emaxiso);
193 }
194#ifdef G4VERBOSE
195 if(verboseLevel > 1) {
196 G4cout << "G4GammaNuclearXS::IsoXS: Z= " << Z << " A= " << A
197 << " Ekin(MeV)= " << ekin/CLHEP::MeV
198 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
199 }
200#endif
201 return xs;
202 }
203 }
204
205 // use element x-section
206 // for the hydrogen target there is no element data
207 if(ekin <= emax && Z != 1) {
208 xs = pv->Value(ekin)*A/aeff[Z];
209
210 // CHIPS for high energy and for the hydrogen target
211 } else if(ekin >= rTransitionBound || Z == 1) {
212 if(Z <= 2 && ekin > 10.*GeV) {
213 xs = coeff[Z][A - amin[Z]]*
214 ggXsection->GetElementCrossSection(aParticle, Z, 0);
215 } else {
216 xs = ggXsection->GetIsoCrossSection(aParticle, Z, A);
217 }
218
219 // transition GDR to CHIPS
220 } else {
221 const G4double rxs = xs150[Z];
222 const G4double lxs = pv->Value(emax)*A/aeff[Z];
223 xs = lxs + (ekin - emax)*(rxs - lxs)/(rTransitionBound-emax);
224 }
225#ifdef G4VERBOSE
226 if(verboseLevel > 1) {
227 G4cout << "G4GammaNuclearXS::IsoXS: Z= " << Z << " A= " << A
228 << " Ekin(MeV)= " << ekin/CLHEP::MeV
229 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
230 }
231#endif
232 return xs;
233}
234
236 const G4Element* anElement, G4double kinEnergy, G4double)
237{
238 std::size_t nIso = anElement->GetNumberOfIsotopes();
239 const G4Isotope* iso = anElement->GetIsotope(0);
240
241 if(1 == nIso) { return iso; }
242
243 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
245 G4double sum = 0.0;
246 G4int j;
247 G4int Z = anElement->GetZasInt();
248
249 // condition to use only isotope abundance
250 if(amax[Z] == amin[Z] || kinEnergy > rTransitionBound || Z >= MAXZGAMMAXS ) {
251 for (j=0; j<(G4int)nIso; ++j) {
252 sum += abundVector[j];
253 if(q <= sum) {
254 iso = anElement->GetIsotope(j);
255 break;
256 }
257 }
258 return iso;
259 }
260 // use isotope cross sections
261 std::size_t nn = temp.size();
262 if(nn < nIso) { temp.resize(nIso, 0.); }
263
264 for (j=0; j<(G4int)nIso; ++j) {
265 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
266 // << " abund= " << abundVector[j] << G4endl;
267 sum += abundVector[j]*
268 IsoCrossSection(kinEnergy, Z, anElement->GetIsotope(j)->GetN());
269 temp[j] = sum;
270 }
271 sum *= q;
272 for (j = 0; j<(G4int)nIso; ++j) {
273 if(temp[j] >= sum) {
274 iso = anElement->GetIsotope(j);
275 break;
276 }
277 }
278 return iso;
279}
280
281void
283{
284 if(verboseLevel > 0){
285 G4cout << "G4GammaNuclearXS::BuildPhysicsTable for "
286 << p.GetParticleName() << G4endl;
287 }
288 if(p.GetParticleName() != "gamma") {
290 ed << p.GetParticleName() << " is a wrong particle type -"
291 << " only gamma is allowed";
292 G4Exception("G4GammaNuclearXS::BuildPhysicsTable(..)","had012",
293 FatalException, ed, "");
294 return;
295 }
296
297 if(nullptr == data) {
298#ifdef G4MULTITHREADED
299 G4MUTEXLOCK(&gNuclearXSMutex);
300 if(nullptr == data) {
301#endif
302 isMaster = true;
303 data = new G4ElementData();
304 data->SetName("PhotoNuclear");
305 FindDirectoryPath();
306#ifdef G4MULTITHREADED
307 }
308 G4MUTEXUNLOCK(&gNuclearXSMutex);
309#endif
310 }
311
312
313 // it is possible re-initialisation for the second run
314 // Upload data for elements used in geometry
316 if(isMaster) {
317 for ( auto & elm : *table ) {
318 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZGAMMAXS-1) );
319 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
320 }
321 }
322
323 // prepare isotope selection
324 std::size_t nIso = temp.size();
325 for ( auto & elm : *table ) {
326 std::size_t n = elm->GetNumberOfIsotopes();
327 if(n > nIso) { nIso = n; }
328 }
329 temp.resize(nIso, 0.0);
330}
331
332const G4String& G4GammaNuclearXS::FindDirectoryPath()
333{
334 // check environment variable
335 // build the complete string identifying the file with the data set
336 if(gDataDirectory.empty()) {
337 const char* path = G4FindDataDir("G4PARTICLEXSDATA");
338 if (nullptr != path) {
339 std::ostringstream ost;
340 ost << path << "/gamma/inel";
341 gDataDirectory = ost.str();
342 } else {
343 G4Exception("G4GammaNuclearXS::Initialise(..)","had013",
345 "Environment variable G4PARTICLEXSDATA is not defined");
346 }
347 }
348 return gDataDirectory;
349}
350
351void G4GammaNuclearXS::InitialiseOnFly(G4int Z)
352{
353#ifdef G4MULTITHREADED
354 G4MUTEXLOCK(&gNuclearXSMutex);
355 if(nullptr == data->GetElementData(Z)) {
356#endif
357 Initialise(Z);
358#ifdef G4MULTITHREADED
359 }
360 G4MUTEXUNLOCK(&gNuclearXSMutex);
361#endif
362}
363
364void G4GammaNuclearXS::Initialise(G4int Z)
365{
366 if(nullptr != data->GetElementData(Z)) { return; }
367
368 // upload data from file
369 std::ostringstream ost;
370 ost << FindDirectoryPath() << Z ;
371 G4PhysicsVector* v = RetrieveVector(ost, true, Z);
372
373 data->InitialiseForElement(Z, v);
374 /*
375 G4cout << "G4GammaNuclearXS::Initialise for Z= " << Z
376 << " A= " << Amean << " Amin= " << amin[Z]
377 << " Amax= " << amax[Z] << G4endl;
378 */
379 // upload isotope data
380 G4DynamicParticle theGamma(gamma, G4ThreeVector(1,0,0), rTransitionBound);
381 xs150[Z] = ggXsection->GetElementCrossSection(&theGamma, Z, 0);
382 if(amax[Z] > amin[Z]) {
383 G4int nmax = amax[Z]-amin[Z]+1;
384 data->InitialiseForComponent(Z, nmax);
385 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
386 std::ostringstream ost1;
387 ost1 << gDataDirectory << Z << "_" << A;
388 G4PhysicsVector* v1 = RetrieveVector(ost1, false, Z);
389 data->AddComponent(Z, A, v1);
390 if(Z<=2){
391 theGamma.SetKineticEnergy(10.*GeV);
392 G4double sig1 = ggXsection->GetIsoCrossSection(&theGamma, Z, A);
393 G4double sig2 = ggXsection->GetElementCrossSection(&theGamma, Z, 0);
394 if(sig2 > 0.) coeff[Z][A-amin[Z]]=(sig1/sig2);
395 else coeff[Z][A-amin[Z]]=1.;
396 }
397 }
398 }
399}
400
402G4GammaNuclearXS::RetrieveVector(std::ostringstream& ost, G4bool isElement, G4int Z)
403{
404 G4PhysicsVector* v = nullptr;
405
406 std::ifstream filein(ost.str().c_str());
407 if (!filein.is_open()) {
408 if(isElement) {
410 ed << "Data file <" << ost.str().c_str()
411 << "> is not opened!";
412 G4Exception("G4GammaNuclearXS::RetrieveVector(..)","had014",
413 FatalException, ed, "Check G4PARTICLEXSDATA");
414 }
415 } else {
416 if(verboseLevel > 1) {
417 G4cout << "File " << ost.str()
418 << " is opened by G4GammaNuclearXS" << G4endl;
419 }
420 // retrieve data from DB
421 if(std::find(std::begin(freeVectorException), std::end(freeVectorException), Z ) == std::end(freeVectorException) && isElement) {
422 v = new G4PhysicsLinearVector();
423 } else {
424 v = new G4PhysicsVector();
425 }
426 if(!v->Retrieve(filein, true)) {
428 ed << "Data file <" << ost.str().c_str()
429 << "> is not retrieved!";
430 G4Exception("G4GammaNuclearXS::RetrieveVector(..)","had015",
431 FatalException, ed, "Check G4PARTICLEXSDATA");
432 }
433 }
434 return v;
435}
436
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
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 * GetComponentDataByIndex(G4int Z, G4int idx)
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
G4PhysicsVector * GetElementData(G4int Z)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
void SetName(const G4String &nam)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZasInt() const
Definition: G4Element.hh:132
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 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
G4int GetN() const
Definition: G4Isotope.hh:93
const G4String & GetParticleName() const
G4double GetMaxEnergy() const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
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)