Geant4 10.7.0
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// Author V.Ivantchenko, Geant4, 22 October 2020
34//
35// Modifications:
36//
37
38#include "G4GammaNuclearXS.hh"
39#include "G4Gamma.hh"
40#include "G4DynamicParticle.hh"
42#include "G4Material.hh"
43#include "G4Element.hh"
44#include "G4PhysicsLogVector.hh"
47#include "Randomize.hh"
48#include "G4SystemOfUnits.hh"
49#include "G4IsotopeList.hh"
50
51#include <fstream>
52#include <sstream>
53
54// factory
56//
58
59G4PhysicsVector* G4GammaNuclearXS::data[] = {nullptr};
60G4double G4GammaNuclearXS::coeff[] = {0.0};
61G4String G4GammaNuclearXS::gDataDirectory = "";
62
63#ifdef G4MULTITHREADED
64 G4Mutex G4GammaNuclearXS::gNuclearXSMutex = G4MUTEX_INITIALIZER;
65#endif
66
68 : G4VCrossSectionDataSet(Default_Name()),
69 ggXsection(nullptr),
70 gamma(G4Gamma::Gamma()),
71 isMaster(false)
72{
73 // verboseLevel = 0;
74 if(verboseLevel > 0){
75 G4cout << "G4GammaNuclearXS::G4GammaNuclearXS Initialise for Z < "
76 << MAXZEL << G4endl;
77 }
79 if(ggXsection == nullptr) ggXsection = new G4PhotoNuclearCrossSection();
81}
82
84{
85 if(isMaster) {
86 for(G4int i=0; i<MAXZEL; ++i) {
87 delete data[i];
88 data[i] = nullptr;
89 }
90 }
91}
92
93void G4GammaNuclearXS::CrossSectionDescription(std::ostream& outFile) const
94{
95 outFile << "G4GammaNuclearXS calculates the gamma nuclear\n"
96 << "cross section on nuclei using data from the high precision\n"
97 << "LEND gamma database. The data are simplified and smoothed over\n"
98 << "the resonance region in order to reduce CPU time.\n"
99 << "For high energies Glauber-Gribiv cross section is used.\n";
100}
101
102G4bool
104 G4int, const G4Material*)
105{
106 return true;
107}
108
110 G4int, G4int,
111 const G4Element*, const G4Material*)
112{
113 return true;
114}
115
118 G4int ZZ, const G4Material* mat)
119{
120 G4double xs = 0.0;
121 G4double ekin = aParticle->GetKineticEnergy();
122
123 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
124
125 auto pv = GetPhysicsVector(Z);
126 if(pv == nullptr) { return xs; }
127 // G4cout << "G4GammaNuclearXS::GetCrossSection e= " << ekin
128 // << " Z= " << Z << G4endl;
129
130 if(ekin <= pv->GetMaxEnergy()) {
131 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
132 } else {
133 xs = coeff[Z]*ggXsection->GetElementCrossSection(aParticle, Z, mat);
134 }
135
136#ifdef G4VERBOSE
137 if(verboseLevel > 1) {
138 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
139 << ", nElmXS(b)= " << xs/CLHEP::barn
140 << G4endl;
141 }
142#endif
143 return xs;
144}
145
147 const G4DynamicParticle* aParticle,
148 G4int Z, G4int A,
149 const G4Isotope*, const G4Element*,
150 const G4Material* mat)
151{
152 return GetElementCrossSection(aParticle, Z, mat) * A/aeff[Z];
153}
154
156 const G4Element* anElement, G4double, G4double)
157{
158 size_t nIso = anElement->GetNumberOfIsotopes();
159 const G4Isotope* iso = anElement->GetIsotope(0);
160
161 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
162 if(1 == nIso) { return iso; }
163
164 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
166 G4double sum = 0.0;
167 size_t j;
168
169 // isotope wise cross section not used
170 for (j=0; j<nIso; ++j) {
171 sum += abundVector[j];
172 if(q <= sum) {
173 iso = anElement->GetIsotope(j);
174 break;
175 }
176 }
177 return iso;
178}
179
180void
182{
183 if(verboseLevel > 0){
184 G4cout << "G4GammaNuclearXS::BuildPhysicsTable for "
185 << p.GetParticleName() << G4endl;
186 }
187 if(p.GetParticleName() != "gamma") {
189 ed << p.GetParticleName() << " is a wrong particle type -"
190 << " only gamma is allowed";
191 G4Exception("G4GammaNuclearXS::BuildPhysicsTable(..)","had012",
192 FatalException, ed, "");
193 return;
194 }
195 if(0. == coeff[0]) {
196#ifdef G4MULTITHREADED
197 G4MUTEXLOCK(&gNuclearXSMutex);
198 if(0. == coeff[0]) {
199#endif
200 coeff[0] = 1.0;
201 isMaster = true;
202 FindDirectoryPath();
203#ifdef G4MULTITHREADED
204 }
205 G4MUTEXUNLOCK(&gNuclearXSMutex);
206#endif
207 }
208
209 // it is possible re-initialisation for the second run
210 if(isMaster) {
211
212 // Access to elements
213 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
214 size_t numOfCouples = theCoupleTable->GetTableSize();
215 for(size_t j=0; j<numOfCouples; ++j) {
216 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
217 auto elmVec = mat->GetElementVector();
218 size_t numOfElem = mat->GetNumberOfElements();
219 for (size_t ie = 0; ie < numOfElem; ++ie) {
220 G4int Z = std::max(1,std::min(((*elmVec)[ie])->GetZasInt(), MAXZEL-1));
221 if(data[Z] == nullptr) { Initialise(Z); }
222 }
223 }
224 }
225}
226
227const G4String& G4GammaNuclearXS::FindDirectoryPath()
228{
229 // check environment variable
230 // build the complete string identifying the file with the data set
231 if(gDataDirectory.empty()) {
232 char* path = std::getenv("G4PARTICLEXSDATA");
233 if (path) {
234 std::ostringstream ost;
235 ost << path << "/gamma/inel";
236 gDataDirectory = ost.str();
237 } else {
238 G4Exception("G4GammaNuclearXS::Initialise(..)","had013",
240 "Environment variable G4PARTICLEXSDATA is not defined");
241 }
242 }
243 return gDataDirectory;
244}
245
246void G4GammaNuclearXS::InitialiseOnFly(G4int Z)
247{
248#ifdef G4MULTITHREADED
249 G4MUTEXLOCK(&gNuclearXSMutex);
250 if(data[Z] == nullptr) {
251#endif
252 Initialise(Z);
253#ifdef G4MULTITHREADED
254 }
255 G4MUTEXUNLOCK(&gNuclearXSMutex);
256#endif
257}
258
259void G4GammaNuclearXS::Initialise(G4int Z)
260{
261 if(data[Z] != nullptr) { return; }
262
263 // upload data from file
264 data[Z] = new G4PhysicsLogVector();
265
266 std::ostringstream ost;
267 ost << FindDirectoryPath() << Z ;
268 std::ifstream filein(ost.str().c_str());
269 if (!(filein)) {
271 ed << "Data file <" << ost.str().c_str()
272 << "> is not opened!";
273 G4Exception("G4GammaNuclearXS::Initialise(..)","had014",
274 FatalException, ed, "Check G4PARTICLEXSDATA");
275 return;
276 }
277 if(verboseLevel > 1) {
278 G4cout << "file " << ost.str()
279 << " is opened by G4GammaNuclearXS" << G4endl;
280 }
281
282 // retrieve data from DB
283 if(!data[Z]->Retrieve(filein, true)) {
285 ed << "Data file <" << ost.str().c_str()
286 << "> is not retrieved!";
287 G4Exception("G4GammaNuclearXS::Initialise(..)","had015",
288 FatalException, ed, "Check G4PARTICLEXSDATA");
289 return;
290 }
291 // smooth transition
292 G4Material* mat(nullptr);
293 G4ThreeVector mom(0.0,0.0,1.0);
294 G4DynamicParticle dp(gamma, mom, data[Z]->GetMaxEnergy());
295 G4double sig1 = (*(data[Z]))[data[Z]->GetVectorLength()-1];
296 G4double sig2 = ggXsection->GetElementCrossSection(&dp, Z, mat);
297 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
298}
#define G4_DECLARE_XS_FACTORY(cross_section)
double A(double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
const G4int MAXZEL
#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
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
static G4CrossSectionDataSetRegistry * Instance()
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *mat) final
void CrossSectionDescription(std::ostream &) const final
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
const G4String & GetParticleName() const
std::size_t GetVectorLength() const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
void SetForAllAtomsAndEnergies(G4bool val)