Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronInelasticXS.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: G4NeutronInelasticXS
32//
33// Author Ivantchenko, Geant4, 3-Aug-09
34//
35
37#include "G4Neutron.hh"
38#include "G4DynamicParticle.hh"
39#include "G4ElementTable.hh"
40#include "G4Material.hh"
41#include "G4Element.hh"
42#include "G4PhysicsLogVector.hh"
46#include "Randomize.hh"
47#include "G4Neutron.hh"
48#include "G4SystemOfUnits.hh"
49#include "G4IsotopeList.hh"
50#include "G4NuclearRadii.hh"
51#include "G4AutoLock.hh"
52
53#include <fstream>
54#include <sstream>
55#include <thread>
56
57G4double G4NeutronInelasticXS::coeff[] = {1.0};
58G4double G4NeutronInelasticXS::lowcoeff[] = {1.0};
59G4ElementData* G4NeutronInelasticXS::data = nullptr;
60G4ElementData* G4NeutronInelasticXS::dataR = nullptr;
61G4String G4NeutronInelasticXS::gDataDirectory = "";
62
63static std::once_flag applyOnce;
64
65namespace
66{
67 G4Mutex nInelasticXSMutex = G4MUTEX_INITIALIZER;
68}
69
72 neutron(G4Neutron::Neutron()),
73 lowElimit(1.0e-7*CLHEP::eV)
74{
75 verboseLevel = 0;
76 if (verboseLevel > 0) {
77 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
78 << MAXZINEL << G4endl;
79 }
80 loglowElimit = G4Log(lowElimit);
81 ggXsection =
83 if (ggXsection == nullptr)
84 ggXsection = new G4ComponentGGHadronNucleusXsc();
85
86 if (nullptr == data) {
87 data = new G4ElementData(MAXZINEL);
88 data->SetName("nInelastic");
89 dataR = new G4ElementData(MAXZINEL);
90 dataR->SetName("nRInelastic");
91 FindDirectoryPath();
92 for (G4int Z=1; Z<MAXZINEL; ++Z) { Initialise(Z); }
93 }
94
96}
97
98void G4NeutronInelasticXS::CrossSectionDescription(std::ostream& outFile) const
99{
100 outFile << "G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
101 << "cross section on nuclei using data from the high precision\n"
102 << "neutron database. These data are simplified and smoothed over\n"
103 << "the resonance region in order to reduce CPU time.\n"
104 << "For high energy Glauber-Gribov cross section model is used\n";
105}
106
107G4bool
109 G4int, const G4Material*)
110{
111 return true;
112}
113
114G4bool
116 G4int, G4int,
117 const G4Element*, const G4Material*)
118{
119 return true;
120}
121
124 G4int Z, const G4Material*)
125{
126 return ElementCrossSection(aParticle->GetKineticEnergy(),
127 aParticle->GetLogKineticEnergy(), Z);
128}
129
133 const G4Element* elm,
134 const G4Material*)
135{
136 return ElementCrossSection(ekin, loge, elm->GetZasInt());
137}
138
141{
142 G4int Z = std::min(ZZ, MAXZINEL-1);
143 G4double ekin = eKin;
144 G4double loge = logE;
145 G4double xs{0.0};
146 G4bool done{false};
147
148 // very low energy limit
149 if (ekin < lowElimit) {
150 ekin = lowElimit;
151 loge = loglowElimit;
152 }
153
154 // data from the resonance region
155 if (fRfilesEnabled) {
156 auto pv = GetPhysicsVectorR(Z);
157 if (nullptr != pv && ekin < inel_max_r_e[Z]) {
158 xs = pv->LogVectorValue(ekin, loge);
159 done = true;
160 }
161 }
162 // data above the resonance region pv should always defined
163 if (!done) {
164 auto pv = GetPhysicsVector(Z);
165 if (ekin <= pv->GetMaxEnergy()) {
166 xs = pv->LogVectorValue(ekin, loge);
167 }
168 else {
169 xs = coeff[Z]*ggXsection->GetInelasticElementCrossSection(neutron, ekin,
170 Z, aeff[Z]);
171 }
172 }
173
174#ifdef G4VERBOSE
175 if(verboseLevel > 1) {
176 G4cout << "G4NeutronInelasticXS::ElementCrossSection Z= " << Z
177 << " Ekin(MeV)= " << ekin/CLHEP::MeV
178 << ", ElmXSinel(b)= " << xs/CLHEP::barn
179 << G4endl;
180 }
181#endif
182 return xs;
183}
184
188 G4int Z, G4int A,
189 const G4Isotope*, const G4Element*,
190 const G4Material*)
191{
192 return IsoCrossSection(ekin, loge, Z, A);
193}
194
197 G4int Z, G4int A,
198 const G4Isotope*, const G4Element*,
199 const G4Material*)
200{
201 return IsoCrossSection(aParticle->GetKineticEnergy(),
202 aParticle->GetLogKineticEnergy(), Z, A);
203}
204
207 G4int ZZ, G4int A)
208{
209 G4double xs{0.0};
210 G4int Z = std::min(ZZ, MAXZINEL-1);
211 G4bool done{false};
212
213 // data from the resonance region
214 if (fRfilesEnabled) {
215 auto pv = GetPhysicsVectorR(Z);
216 if (nullptr != pv && ekin < inel_max_r_e[Z]) {
217 // use isotope x-section if possible
218 if (dataR->GetNumberOfComponents(Z) > 0) {
219 auto pviso = dataR->GetComponentDataByID(Z, A);
220 if (pviso != nullptr) {
221 xs = pviso->LogVectorValue(ekin, loge);
222 done = true;
223 }
224 }
225 // isotope data are not available or applicable
226 if (!done) {
227 xs = pv->LogVectorValue(ekin, loge);
228 done = true;
229 }
230 }
231 }
232 // data above the resonance region
233 if (!done) {
234 auto pv = GetPhysicsVector(Z);
235 // use isotope x-section if possible
236 if (data->GetNumberOfComponents(Z) > 0) {
237 auto pviso = data->GetComponentDataByID(Z, A);
238 if (pviso != nullptr && ekin <= pviso->GetMaxEnergy()) {
239 xs = pviso->LogVectorValue(ekin, loge);
240 done = true;
241 }
242 }
243 // isotope data are not available or applicable
244 // use element x-section
245 if (!done) {
246 if (ekin <= pv->GetMaxEnergy()) {
247 xs = pv->LogVectorValue(ekin, loge);
248 }
249 else {
250 xs = coeff[Z]*
251 ggXsection->GetInelasticElementCrossSection(neutron, ekin,
252 Z, aeff[Z]);
253 }
254 }
255 }
256
257#ifdef G4VERBOSE
258 if (verboseLevel > 1) {
259 G4cout << "G4NeutronInelasticXS::IsoXS: Z= " << Z << " A= " << A
260 << " Ekin(MeV)= " << ekin/CLHEP::MeV
261 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
262 }
263#endif
264 return xs;
265}
266
268 const G4Element* anElement, G4double kinEnergy, G4double logE)
269{
270 std::size_t nIso = anElement->GetNumberOfIsotopes();
271 const G4Isotope* iso = anElement->GetIsotope(0);
272 if(1 == nIso) { return iso; }
273
274 // more than 1 isotope
275 G4int Z = anElement->GetZasInt();
276 if (nullptr == data->GetElementData(Z)) { InitialiseOnFly(Z); }
277
278 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
280 G4double sum = 0.0;
281 std::size_t j;
282
283 // isotope wise cross section not available
284 if (Z >= MAXZINEL || 0 == data->GetNumberOfComponents(Z)) {
285 for (j=0; j<nIso; ++j) {
286 sum += abundVector[j];
287 if(q <= sum) {
288 iso = anElement->GetIsotope((G4int)j);
289 break;
290 }
291 }
292 return iso;
293 }
294
295 // use isotope cross sections
296 auto nn = temp.size();
297 if(nn < nIso) { temp.resize(nIso, 0.); }
298
299 for (j=0; j<nIso; ++j) {
300 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
301 anElement->GetIsotope((G4int)j)->GetN());
302 temp[j] = sum;
303 }
304 sum *= q;
305 for (j = 0; j<nIso; ++j) {
306 if (temp[j] >= sum) {
307 iso = anElement->GetIsotope((G4int)j);
308 break;
309 }
310 }
311 return iso;
312}
313
314void
316{
317 if (verboseLevel > 0) {
318 G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
319 << p.GetParticleName() << G4endl;
320 }
321 if (p.GetParticleName() != "neutron") {
323 ed << p.GetParticleName() << " is a wrong particle type -"
324 << " only neutron is allowed";
325 G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
326 FatalException, ed, "");
327 return;
328 }
329
331
332 // it is possible re-initialisation for the new run
334
335 // initialise static tables only once
336 std::call_once(applyOnce, [this]() { isInitializer = true; });
337
338 if (isInitializer) {
339 G4AutoLock l(&nInelasticXSMutex);
340
341 // Upload data for elements used in geometry
342 for ( auto const & elm : *table ) {
343 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINEL-1) );
344 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
345 }
346 l.unlock();
347 }
348 // prepare isotope selection
349 std::size_t nIso = temp.size();
350 for ( auto const & elm : *table ) {
351 std::size_t n = elm->GetNumberOfIsotopes();
352 if (n > nIso) { nIso = n; }
353 }
354 temp.resize(nIso, 0.0);
355}
356
357const G4String& G4NeutronInelasticXS::FindDirectoryPath()
358{
359 // build the complete string identifying the file with the data set
360 if (gDataDirectory.empty()) {
361 std::ostringstream ost;
362 ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/neutron/";
363 gDataDirectory = ost.str();
364 }
365 return gDataDirectory;
366}
367
368void G4NeutronInelasticXS::InitialiseOnFly(G4int Z)
369{
370 G4AutoLock l(&nInelasticXSMutex);
371 Initialise(Z);
372 l.unlock();
373}
374
375void G4NeutronInelasticXS::Initialise(G4int Z)
376{
377 if (nullptr != data->GetElementData(Z)) { return; }
378
379 // upload element data
380 std::ostringstream ost;
381 ost << FindDirectoryPath() << "inel" << Z;
382 G4PhysicsVector* v = RetrieveVector(ost, true);
383 data->InitialiseForElement(Z, v);
384
385 G4PhysicsVector* vr = nullptr;
386 if (fRfilesEnabled) {
387 std::ostringstream ostr;
388 ostr << FindDirectoryPath() << "Rinel" << Z;
389 vr = RetrieveVector(ostr, false);
390 dataR->InitialiseForElement(Z, vr);
391 }
392
393 // upload isotope data
394 G4bool noComp = true;
395 G4bool noCompR = true;
396 if (amin[Z] < amax[Z]) {
397
398 for (G4int A=amin[Z]; A<=amax[Z]; ++A) {
399 std::ostringstream ost1;
400 ost1 << gDataDirectory << "inel" << Z << "_" << A;
401 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
402 if (nullptr != v1) {
403 if (noComp) {
404 G4int nmax = amax[Z] - A + 1;
405 data->InitialiseForComponent(Z, nmax);
406 noComp = false;
407 }
408 data->AddComponent(Z, A, v1);
409 }
410 if (nullptr != vr) {
411 std::ostringstream ost2;
412 ost2 << gDataDirectory << "Rinel" << Z << "_" << A;
413 G4PhysicsVector* v2 = RetrieveVector(ost2, false);
414 if (nullptr != v2) {
415 if (noCompR) {
416 G4int nmax = amax[Z] - A + 1;
417 dataR->InitialiseForComponent(Z, nmax);
418 noCompR = false;
419 }
420 dataR->AddComponent(Z, A, v2);
421 }
422 }
423 }
424 }
425 // no components case
426 if (noComp) { data->InitialiseForComponent(Z, 0); }
427 if (noCompR) { dataR->InitialiseForComponent(Z, 0); }
428
429 // smooth transition
430 G4double sig1 = (*v)[v->GetVectorLength()-1];
431 G4double ehigh= v->GetMaxEnergy();
432 G4double sig2 = ggXsection->GetInelasticElementCrossSection(neutron,
433 ehigh, Z, aeff[Z]);
434 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
435}
436
438G4NeutronInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
439{
440 G4PhysicsLogVector* v = nullptr;
441 std::ifstream filein(ost.str().c_str());
442 if (!filein.is_open()) {
443 if (warn) {
445 ed << "Data file <" << ost.str().c_str()
446 << "> is not opened!";
447 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had014",
448 FatalException, ed, "Check G4PARTICLEXSDATA");
449 }
450 } else {
451 if (verboseLevel > 1) {
452 G4cout << "File " << ost.str()
453 << " is opened by G4NeutronInelasticXS" << G4endl;
454 }
455 // retrieve data from DB
456 v = new G4PhysicsLogVector();
457 if (!v->Retrieve(filein, true)) {
459 ed << "Data file <" << ost.str().c_str()
460 << "> is not retrieved!";
461 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had015",
462 FatalException, ed, "Check G4PARTICLEXSDATA");
463 }
464 }
465 return v;
466}
G4TemplateAutoLock< G4Mutex > G4AutoLock
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4double G4Log(G4double x)
Definition G4Log.hh:169
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
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
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
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
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
G4int GetZasInt() const
Definition G4Element.hh:120
static G4HadronicParameters * Instance()
const G4String & GetDirPARTICLEXS() const
G4int GetN() const
Definition G4Isotope.hh:83
G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
static const char * Default_Name()
G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) 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 *) final
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
const G4String & GetParticleName() const
G4double GetMaxEnergy() const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
std::size_t GetVectorLength() const
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)