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