Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronCaptureXS.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: G4NeutronCaptureXS
32//
33// Author Ivantchenko, Geant4, 3-Aug-09
34//
35// Modifications:
36//
37
38#include <fstream>
39#include <sstream>
40#include <thread>
41
42#include "G4SystemOfUnits.hh"
43#include "G4NeutronCaptureXS.hh"
44#include "G4Material.hh"
45#include "G4Element.hh"
46#include "G4PhysicsLogVector.hh"
47#include "G4DynamicParticle.hh"
48#include "G4ElementTable.hh"
49#include "G4IsotopeList.hh"
51#include "Randomize.hh"
52#include "G4Log.hh"
53#include "G4AutoLock.hh"
54
55G4ElementData* G4NeutronCaptureXS::data = nullptr;
56G4ElementData* G4NeutronCaptureXS::dataR = nullptr;
57G4String G4NeutronCaptureXS::gDataDirectory = "";
58
59static std::once_flag applyOnce;
60
61namespace
62{
63 G4Mutex neutronCaptureXSMutex = G4MUTEX_INITIALIZER;
64 const G4int MAXZCAPTURE = 92;
65}
66
69 emax(20*CLHEP::MeV), elimit(1.0e-5*CLHEP::eV)
70{
71 verboseLevel = 0;
72 if (verboseLevel > 0) {
73 G4cout << "G4NeutronCaptureXS::G4NeutronCaptureXS: Initialise for Z < "
74 << MAXZCAPTURE << G4endl;
75 }
76 logElimit = G4Log(elimit);
77 if (nullptr == data) {
78 data = new G4ElementData(MAXZCAPTURE+1);
79 data->SetName("nCapture");
80 dataR = new G4ElementData(MAXZCAPTURE+1);
81 dataR->SetName("nRCapture");
82 FindDirectoryPath();
83 }
84}
85
86void G4NeutronCaptureXS::CrossSectionDescription(std::ostream& outFile) const
87{
88 outFile << "G4NeutronCaptureXS calculates the neutron capture cross sections\n"
89 << "on nuclei using data from the high precision neutron database.\n"
90 << "These data are simplified and smoothed over the resonance region\n"
91 << "in order to reduce CPU time. G4NeutronCaptureXS is set to zero\n"
92 << "above 20 MeV for all targets. For Z > 92 the cross section of\n"
93 << "Uranium is used.\n";
94}
95
96G4bool
98 G4int, const G4Material*)
99{
100 return true;
101}
102
103G4bool
105 G4int, G4int,
106 const G4Element*, const G4Material*)
107{
108 return true;
109}
110
113 G4int Z, const G4Material*)
114{
115 G4double xs = 0.0;
116 G4double ekin = aParticle->GetKineticEnergy();
117 if (ekin < emax) {
118 xs = ElementCrossSection(ekin, aParticle->GetLogKineticEnergy(), Z);
119 }
120 return xs;
121}
122
126 const G4Element* elm,
127 const G4Material*)
128{
129 G4double xs = 0.0;
130 if (ekin < emax) {
131 xs = ElementCrossSection(ekin, loge, elm->GetZasInt());
132 }
133 return xs;
134}
135
138{
139 G4int Z = std::min(ZZ, MAXZCAPTURE);
140 G4double ekin = eKin;
141 G4double logEkin = logE;
142 if (ekin < elimit) {
143 ekin = elimit;
144 logEkin = logElimit;
145 }
146
147 G4double xs;
148 G4bool done{false};
149
150 // data from the resonance region
151 if (fRfilesEnabled) {
152 auto pv = GetPhysicsVectorR(Z);
153 if (nullptr != pv && ekin < cap_max_r_e[Z]) {
154 const G4double e0 = pv->Energy(0);
155 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
156 : (*pv)[0]*std::sqrt(e0/ekin);
157 done = true;
158 }
159 }
160 // data above the resonance region
161 if (!done) {
162 auto pv = GetPhysicsVector(Z);
163 const G4double e0 = pv->Energy(0);
164 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
165 : (*pv)[0]*std::sqrt(e0/ekin);
166 }
167
168#ifdef G4VERBOSE
169 if (verboseLevel > 1){
170 G4cout << "Ekin= " << ekin/CLHEP::MeV
171 << " ElmXScap(b)= " << xs/CLHEP::barn << G4endl;
172 }
173#endif
174 return xs;
175}
176
180 G4int Z, G4int A,
181 const G4Isotope*, const G4Element*,
182 const G4Material*)
183{
184 return IsoCrossSection(ekin, loge, Z, A);
185}
186
189 G4int Z, G4int A,
190 const G4Isotope*, const G4Element*,
191 const G4Material*)
192{
193 return IsoCrossSection(aParticle->GetKineticEnergy(),
194 aParticle->GetLogKineticEnergy(),
195 Z, A);
196}
197
199 G4int ZZ, G4int A)
200{
201 G4double xs = 0.0;
202 if (eKin > emax) { return xs; }
203
204 G4int Z = std::min(ZZ, MAXZCAPTURE);
205 G4double ekin = eKin;
206 G4double logEkin = logE;
207 if (ekin < elimit) {
208 ekin = elimit;
209 logEkin = logElimit;
210 }
211
212 G4bool done{false};
213
214 // data from the resonance region
215 if (fRfilesEnabled) {
216 auto pv = GetPhysicsVectorR(Z);
217 if (nullptr != pv && ekin < cap_max_r_e[Z]) {
218 // use isotope x-section if possible
219 if (dataR->GetNumberOfComponents(Z) > 0) {
220 auto pviso = dataR->GetComponentDataByID(Z, A);
221 if (pviso != nullptr) {
222 const G4double e0 = pviso->Energy(0);
223 xs = (ekin >= e0) ? pviso->LogVectorValue(ekin, logEkin)
224 : (*pviso)[0]*std::sqrt(e0/ekin);
225 done = true;
226 }
227 }
228 // isotope data are not available or applicable
229 if (!done) {
230 const G4double e0 = pv->Energy(0);
231 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
232 : (*pv)[0]*std::sqrt(e0/ekin);
233 done = true;
234 }
235 }
236 }
237 // data above the resonance region
238 if (!done) {
239 auto pv = GetPhysicsVector(Z);
240 // use isotope x-section if possible
241 if (data->GetNumberOfComponents(Z) > 0) {
242 auto pviso = data->GetComponentDataByID(Z, A);
243 if (pviso != nullptr) {
244 const G4double e0 = pviso->Energy(0);
245 xs = (ekin >= e0) ? pviso->LogVectorValue(ekin, logEkin)
246 : (*pviso)[0]*std::sqrt(e0/ekin);
247 done = true;
248 }
249 }
250 // isotope data are not available or applicable
251 if (!done) {
252 const G4double e0 = pv->Energy(0);
253 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
254 : (*pv)[0]*std::sqrt(e0/ekin);
255 }
256 }
257#ifdef G4VERBOSE
258 if (verboseLevel > 0) {
259 G4cout << "G4NeutronCaptureXS::IsoXS: Ekin(MeV)= " << ekin/MeV
260 << " xs(b)= " << xs/CLHEP::barn
261 << " Z= " << Z << " A= " << A << " no iso XS" << G4endl;
262 }
263#endif
264 return xs;
265}
266
267const G4Isotope*
269 G4double kinEnergy, G4double logE)
270{
271 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
272 const G4Isotope* iso = anElement->GetIsotope(0);
273
274 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
275 if(1 == nIso) { return iso; }
276
277 // more than 1 isotope
278 G4int Z = anElement->GetZasInt();
279 if (nullptr == data->GetElementData(Z)) { InitialiseOnFly(Z); }
280
281 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
283 G4double sum = 0.0;
284
285 // is there isotope wise cross section?
286 G4int j;
287 if (Z > MAXZCAPTURE || 0 == data->GetNumberOfComponents(Z)) {
288 for (j = 0; j<nIso; ++j) {
289 sum += abundVector[j];
290 if (q <= sum) {
291 iso = anElement->GetIsotope(j);
292 break;
293 }
294 }
295 return iso;
296 }
297 G4int nn = (G4int)temp.size();
298 if (nn < nIso) { temp.resize(nIso, 0.); }
299
300 for (j=0; j<nIso; ++j) {
301 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
302 anElement->GetIsotope(j)->GetN());
303 temp[j] = sum;
304 }
305 sum *= q;
306 for (j = 0; j<nIso; ++j) {
307 if (temp[j] >= sum) {
308 iso = anElement->GetIsotope(j);
309 break;
310 }
311 }
312 return iso;
313}
314
315void
317{
318 if (verboseLevel > 0){
319 G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
320 << p.GetParticleName() << G4endl;
321 }
322 if (p.GetParticleName() != "neutron") {
324 ed << p.GetParticleName() << " is a wrong particle type -"
325 << " only neutron is allowed";
326 G4Exception("G4NeutronCaptureXS::BuildPhysicsTable(..)","had012",
327 FatalException, ed, "");
328 return;
329 }
330
332
333 // it is possible re-initialisation for the second run
335
336 // initialise static tables only once
337 std::call_once(applyOnce, [this]() { isInitializer = true; });
338
339 if (isInitializer) {
340 G4AutoLock l(&neutronCaptureXSMutex);
341 // Access to elements
342 for ( auto const & elm : *table ) {
343 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZCAPTURE) );
344 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
345 }
346 l.unlock();
347 }
348
349 // prepare isotope selection
350 std::size_t nIso = temp.size();
351 for ( auto const & elm : *table ) {
352 std::size_t n = elm->GetNumberOfIsotopes();
353 if (n > nIso) { nIso = n; }
354 }
355 temp.resize(nIso, 0.0);
356}
357
358const G4String& G4NeutronCaptureXS::FindDirectoryPath()
359{
360 // build the complete string identifying the file with the data set
361 if(gDataDirectory.empty()) {
362 std::ostringstream ost;
363 ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/neutron/";
364 gDataDirectory = ost.str();
365 }
366 return gDataDirectory;
367}
368
369void G4NeutronCaptureXS::InitialiseOnFly(G4int Z)
370{
371 G4AutoLock l(&neutronCaptureXSMutex);
372 Initialise(Z);
373 l.unlock();
374}
375
376void G4NeutronCaptureXS::Initialise(G4int Z)
377{
378 if (nullptr != data->GetElementData(Z)) { return; }
379
380 // upload element data
381 std::ostringstream ost;
382 ost << FindDirectoryPath() << "cap" << Z;
383 G4PhysicsVector* v = RetrieveVector(ost, true);
384 data->InitialiseForElement(Z, v);
385 G4PhysicsVector* vr = nullptr;
386 if (fRfilesEnabled) {
387 std::ostringstream ostr;
388 ostr << FindDirectoryPath() << "Rcap" << 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 for (G4int A=amin[Z]; A<=amax[Z]; ++A) {
398 std::ostringstream ost1;
399 ost1 << gDataDirectory << "cap" << Z << "_" << A;
400 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
401 if (nullptr != v1) {
402 if (noComp) {
403 G4int nmax = amax[Z] - A + 1;
404 data->InitialiseForComponent(Z, nmax);
405 noComp = false;
406 }
407 data->AddComponent(Z, A, v1);
408 }
409 if (nullptr != vr) {
410 std::ostringstream ost2;
411 ost2 << gDataDirectory << "Rcap" << Z << "_" << A;
412 G4PhysicsVector* v2 = RetrieveVector(ost2, false);
413 if (nullptr != v2) {
414 if (noCompR) {
415 G4int nmax = amax[Z] - A + 1;
416 dataR->InitialiseForComponent(Z, nmax);
417 noCompR = false;
418 }
419 dataR->AddComponent(Z, A, v2);
420 }
421 }
422 }
423 }
424 // no components case
425 if (noComp) { data->InitialiseForComponent(Z, 0); }
426 if (noCompR && nullptr != vr) { dataR->InitialiseForComponent(Z, 0); }
427}
428
430G4NeutronCaptureXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
431{
432 G4PhysicsLogVector* v = nullptr;
433 std::ifstream filein(ost.str().c_str());
434 if (!filein.is_open()) {
435 if (warn) {
437 ed << "Data file <" << ost.str().c_str()
438 << "> is not opened!";
439 G4Exception("G4NeutronCaptureXS::RetrieveVector(..)","had014",
440 FatalException, ed, "Check G4PARTICLEXSDATA");
441 }
442 } else {
443 if (verboseLevel > 1) {
444 G4cout << "File " << ost.str()
445 << " is opened by G4NeutronCaptureXS" << G4endl;
446 }
447 // retrieve data from DB
448 v = new G4PhysicsLogVector();
449 if (!v->Retrieve(filein, true)) {
451 ed << "Data file <" << ost.str().c_str()
452 << "> is not retrieved!";
453 G4Exception("G4NeutronCaptureXS::RetrieveVector(..)","had015",
454 FatalException, ed, "Check G4PARTICLEXSDATA");
455 }
456 }
457 return v;
458}
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
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 ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
void CrossSectionDescription(std::ostream &) const final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
static const char * Default_Name()
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) final
const G4String & GetParticleName() const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4VCrossSectionDataSet(const G4String &nam="")