Geant4 11.2.2
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;
56G4String G4NeutronCaptureXS::gDataDirectory = "";
57
58static std::once_flag applyOnce;
59
60namespace
61{
62 G4Mutex neutronCaptureXSMutex = G4MUTEX_INITIALIZER;
63}
64
66 : G4VCrossSectionDataSet(Default_Name()),
67 emax(20*CLHEP::MeV), elimit(1.0e-10*CLHEP::eV)
68{
69 verboseLevel = 0;
70 if (verboseLevel > 0) {
71 G4cout << "G4NeutronCaptureXS::G4NeutronCaptureXS: Initialise for Z < "
72 << MAXZCAPTURE << G4endl;
73 }
74 logElimit = G4Log(elimit);
75 if (nullptr == data) {
76 data = new G4ElementData(MAXZCAPTURE);
77 data->SetName("nCapture");
78 FindDirectoryPath();
79 }
80}
81
82void G4NeutronCaptureXS::CrossSectionDescription(std::ostream& outFile) const
83{
84 outFile << "G4NeutronCaptureXS calculates the neutron capture cross sections\n"
85 << "on nuclei using data from the high precision neutron database.\n"
86 << "These data are simplified and smoothed over the resonance region\n"
87 << "in order to reduce CPU time. G4NeutronCaptureXS is set to zero\n"
88 << "above 20 MeV for all targets. For Z > 92 the cross section of\n"
89 << "Uranium is used.\n";
90}
91
92G4bool
98
99G4bool
101 G4int, G4int,
102 const G4Element*, const G4Material*)
103{
104 return true;
105}
106
109 G4int Z, const G4Material*)
110{
111 G4double xs = 0.0;
112 G4double ekin = aParticle->GetKineticEnergy();
113 if (ekin < emax) {
114 xs = ElementCrossSection(ekin, aParticle->GetLogKineticEnergy(), Z);
115 }
116 return xs;
117}
118
122 const G4Element* elm,
123 const G4Material*)
124{
125 G4double xs = 0.0;
126 if (ekin < emax) {
127 xs = ElementCrossSection(ekin, loge, elm->GetZasInt());
128 }
129 return xs;
130}
131
134{
135 G4int Z = std::min(ZZ, MAXZCAPTURE-1);
136 G4double logEkin = loge;
137 if (ekin < elimit) { ekin = elimit; logEkin = logElimit; }
138
139 auto pv = GetPhysicsVector(Z);
140 const G4double e1 = pv->Energy(1);
141 G4double xs = (ekin >= e1) ? pv->LogVectorValue(ekin, logEkin)
142 : (*pv)[1]*std::sqrt(e1/ekin);
143
144#ifdef G4VERBOSE
145 if (verboseLevel > 1){
146 G4cout << "Ekin= " << ekin/CLHEP::MeV
147 << " ElmXScap(b)= " << xs/CLHEP::barn << G4endl;
148 }
149#endif
150 return xs;
151}
152
156 G4int Z, G4int A,
157 const G4Isotope*, const G4Element*,
158 const G4Material*)
159{
160 return IsoCrossSection(ekin, loge, Z, A);
161}
162
165 G4int Z, G4int A,
166 const G4Isotope*, const G4Element*,
167 const G4Material*)
168{
169 return IsoCrossSection(aParticle->GetKineticEnergy(),
170 aParticle->GetLogKineticEnergy(),
171 Z, A);
172}
173
175 G4int ZZ, G4int A)
176{
177 G4double xs = 0.0;
178 if (eKin > emax) { return xs; }
179
180 G4int Z = std::min(ZZ, MAXZCAPTURE-1);
181 G4double ekin = eKin;
182 G4double logEkin = logE;
183 if (ekin < elimit) {
184 ekin = elimit;
185 logEkin = logElimit;
186 }
187
188 auto pv = GetPhysicsVector(Z);
189 if (pv == nullptr) { return xs; }
190
191 if (data->GetNumberOfComponents(Z) > 0) {
192 G4PhysicsVector* pviso = data->GetComponentDataByID(Z, A);
193 if(pviso != nullptr) {
194 const G4double e1 = pviso->Energy(1);
195 xs = (ekin >= e1) ? pviso->LogVectorValue(ekin, logEkin)
196 : (*pviso)[1]*std::sqrt(e1/ekin);
197#ifdef G4VERBOSE
198 if(verboseLevel > 0) {
199 G4cout << "G4NeutronCaptureXS::IsoXS: Ekin(MeV)= " << ekin/MeV
200 << " xs(b)= " << xs/barn
201 << " Z= " << Z << " A= " << A << G4endl;
202 }
203#endif
204 return xs;
205 }
206 }
207 // isotope data are not available or applicable
208 const G4double e1 = pv->Energy(1);
209 xs = (ekin >= e1) ? pv->LogVectorValue(ekin, logEkin)
210 : (*pv)[1]*std::sqrt(e1/ekin);
211#ifdef G4VERBOSE
212 if (verboseLevel > 0) {
213 G4cout << "G4NeutronCaptureXS::IsoXS: Ekin(MeV)= " << ekin/MeV
214 << " xs(b)= " << xs/barn
215 << " Z= " << Z << " A= " << A << " no iso XS" << G4endl;
216 }
217#endif
218 return xs;
219}
220
221const G4Isotope*
223 G4double kinEnergy, G4double logE)
224{
225 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
226 const G4Isotope* iso = anElement->GetIsotope(0);
227
228 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
229 if(1 == nIso) { return iso; }
230
231 // more than 1 isotope
232 G4int Z = anElement->GetZasInt();
233 if (nullptr == data->GetElementData(Z)) { InitialiseOnFly(Z); }
234
235 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
237 G4double sum = 0.0;
238
239 // is there isotope wise cross section?
240 G4int j;
241 if (Z >= MAXZCAPTURE || 0 == data->GetNumberOfComponents(Z)) {
242 for (j = 0; j<nIso; ++j) {
243 sum += abundVector[j];
244 if(q <= sum) {
245 iso = anElement->GetIsotope(j);
246 break;
247 }
248 }
249 return iso;
250 }
251 G4int nn = (G4int)temp.size();
252 if (nn < nIso) { temp.resize(nIso, 0.); }
253
254 for (j=0; j<nIso; ++j) {
255 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
256 anElement->GetIsotope(j)->GetN());
257 temp[j] = sum;
258 }
259 sum *= q;
260 for (j = 0; j<nIso; ++j) {
261 if (temp[j] >= sum) {
262 iso = anElement->GetIsotope(j);
263 break;
264 }
265 }
266 return iso;
267}
268
269void
271{
272 if (verboseLevel > 0){
273 G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
274 << p.GetParticleName() << G4endl;
275 }
276 if (p.GetParticleName() != "neutron") {
278 ed << p.GetParticleName() << " is a wrong particle type -"
279 << " only neutron is allowed";
280 G4Exception("G4NeutronCaptureXS::BuildPhysicsTable(..)","had012",
281 FatalException, ed, "");
282 return;
283 }
284
285 // it is possible re-initialisation for the second run
287
288 // initialise static tables only once
289 std::call_once(applyOnce, [this]() { isInitializer = true; });
290
291 if (isInitializer) {
292 G4AutoLock l(&neutronCaptureXSMutex);
293 // Access to elements
294 for ( auto const & elm : *table ) {
295 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZCAPTURE-1) );
296 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
297 }
298 l.unlock();
299 }
300
301 // prepare isotope selection
302 std::size_t nIso = temp.size();
303 for ( auto const & elm : *table ) {
304 std::size_t n = elm->GetNumberOfIsotopes();
305 if (n > nIso) { nIso = n; }
306 }
307 temp.resize(nIso, 0.0);
308}
309
310const G4String& G4NeutronCaptureXS::FindDirectoryPath()
311{
312 // build the complete string identifying the file with the data set
313 if(gDataDirectory.empty()) {
314 std::ostringstream ost;
315 ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/neutron/cap";
316 gDataDirectory = ost.str();
317 }
318 return gDataDirectory;
319}
320
321void G4NeutronCaptureXS::InitialiseOnFly(G4int Z)
322{
323 G4AutoLock l(&neutronCaptureXSMutex);
324 Initialise(Z);
325 l.unlock();
326}
327
328void G4NeutronCaptureXS::Initialise(G4int Z)
329{
330 if (nullptr != data->GetElementData(Z)) { return; }
331
332 // upload element data
333 std::ostringstream ost;
334 ost << FindDirectoryPath() << Z ;
335 G4PhysicsVector* v = RetrieveVector(ost, true);
336 data->InitialiseForElement(Z, v);
337
338 // upload isotope data
339 G4bool noComp = true;
340 if (amin[Z] < amax[Z]) {
341 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
342 std::ostringstream ost1;
343 ost1 << gDataDirectory << Z << "_" << A;
344 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
345 if (nullptr != v1) {
346 if (noComp) {
347 G4int nmax = amax[Z] - A + 1;
348 data->InitialiseForComponent(Z, nmax);
349 noComp = false;
350 }
351 data->AddComponent(Z, A, v1);
352 }
353 }
354 }
355 // no components case
356 if (noComp) { data->InitialiseForComponent(Z, 0); }
357}
358
360G4NeutronCaptureXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
361{
362 G4PhysicsLogVector* v = nullptr;
363 std::ifstream filein(ost.str().c_str());
364 if (!filein.is_open()) {
365 if (warn) {
367 ed << "Data file <" << ost.str().c_str()
368 << "> is not opened!";
369 G4Exception("G4NeutronCaptureXS::RetrieveVector(..)","had014",
370 FatalException, ed, "Check G4PARTICLEXSDATA");
371 }
372 } else {
373 if (verboseLevel > 1) {
374 G4cout << "File " << ost.str()
375 << " is opened by G4NeutronCaptureXS" << G4endl;
376 }
377 // retrieve data from DB
378 v = new G4PhysicsLogVector();
379 if (!v->Retrieve(filein, true)) {
381 ed << "Data file <" << ost.str().c_str()
382 << "> is not retrieved!";
383 G4Exception("G4NeutronCaptureXS::RetrieveVector(..)","had015",
384 FatalException, ed, "Check G4PARTICLEXSDATA");
385 }
386 }
387 return v;
388}
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:227
#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
G4PhysicsVector * GetElementData(G4int Z) const
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
G4PhysicsVector * GetComponentDataByID(G4int Z, G4int id) const
std::size_t GetNumberOfComponents(G4int Z) const
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
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
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) final
const G4String & GetParticleName() const
G4double Energy(const std::size_t index) const
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)