Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
54
55#ifdef G4MULTITHREADED
56 G4Mutex G4NeutronCaptureXS::neutronCaptureXSMutex = G4MUTEX_INITIALIZER;
57#endif
58
60 : G4VCrossSectionDataSet(Default_Name()),
61 emax(20*CLHEP::MeV),elimit(1.0e-10*CLHEP::eV)
62{
63 // verboseLevel = 0;
64 if(verboseLevel > 0){
65 G4cout << "G4NeutronCaptureXS::G4NeutronCaptureXS: Initialise for Z < "
66 << MAXZCAPTURE << G4endl;
67 }
69}
70
72{
73 if(isMaster) { delete data; data = nullptr; }
74}
75
76void G4NeutronCaptureXS::CrossSectionDescription(std::ostream& outFile) const
77{
78 outFile << "G4NeutronCaptureXS calculates the neutron capture cross sections\n"
79 << "on nuclei using data from the high precision neutron database.\n"
80 << "These data are simplified and smoothed over the resonance region\n"
81 << "in order to reduce CPU time. G4NeutronCaptureXS is set to zero\n"
82 << "above 20 MeV for all targets. Cross section is zero also for Z>92.\n";
83}
84
85G4bool
87 G4int, const G4Material*)
88{
89 return true;
90}
91
92G4bool
94 G4int, G4int,
95 const G4Element*, const G4Material*)
96{
97 return true;
98}
99
102 G4int ZZ, const G4Material*)
103{
104 G4double xs = 0.0;
105 G4double ekin = aParticle->GetKineticEnergy();
106 if(ekin > emax) { return xs; }
107
108 G4int Z = std::min(ZZ, MAXZCAPTURE-1);
109 G4double logEkin = aParticle->GetLogKineticEnergy();
110 if(ekin < elimit) { ekin = elimit; logEkin = logElimit; }
111
112 auto pv = GetPhysicsVector(Z);
113 if(pv == nullptr) { return xs; }
114
115 const G4double e1 = pv->Energy(1);
116 xs = (ekin >= e1) ? pv->LogVectorValue(ekin, logEkin)
117 : (*pv)[1]*std::sqrt(e1/ekin);
118
119#ifdef G4VERBOSE
120 if(verboseLevel > 1){
121 G4cout << "Ekin= " << ekin/CLHEP::MeV
122 << " ElmXScap(b)= " << xs/CLHEP::barn << G4endl;
123 }
124#endif
125 return xs;
126}
127
130 G4int Z, G4int A,
131 const G4Isotope*, const G4Element*,
132 const G4Material*)
133{
134 return IsoCrossSection(aParticle->GetKineticEnergy(),
135 aParticle->GetLogKineticEnergy(),
136 Z, A);
137}
138
140 G4int ZZ, G4int A)
141{
142 G4double xs = 0.0;
143 if(eKin > emax) { return xs; }
144
145 G4int Z = std::min(ZZ, MAXZCAPTURE-1);
146 G4double ekin = eKin;
147 G4double logEkin = logE;
148 if(ekin < elimit) {
149 ekin = elimit;
150 logEkin = logElimit;
151 }
152
153 auto pv = GetPhysicsVector(Z);
154 if(pv == nullptr) { return xs; }
155
156 if(amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z]) {
158 if(pviso != nullptr) {
159 const G4double e1 = pviso->Energy(1);
160 xs = (ekin >= e1) ? pviso->LogVectorValue(ekin, logEkin)
161 : (*pviso)[1]*std::sqrt(e1/ekin);
162#ifdef G4VERBOSE
163 if(verboseLevel > 0) {
164 G4cout << "G4NeutronCaptureXS::IsoXS: Ekin(MeV)= " << ekin/MeV
165 << " xs(b)= " << xs/barn
166 << " Z= " << Z << " A= " << A << G4endl;
167 }
168#endif
169 return xs;
170 }
171 }
172 // isotope data are not available or applicable
173 const G4double e1 = pv->Energy(1);
174 xs = (ekin >= e1) ? pv->LogVectorValue(ekin, logEkin)
175 : (*pv)[1]*std::sqrt(e1/ekin);
176#ifdef G4VERBOSE
177 if(verboseLevel > 0) {
178 G4cout << "G4NeutronCaptureXS::IsoXS: Ekin(MeV)= " << ekin/MeV
179 << " xs(b)= " << xs/barn
180 << " Z= " << Z << " A= " << A << " no iso XS" << G4endl;
181 }
182#endif
183 return xs;
184}
185
186const G4Isotope*
188 G4double kinEnergy, G4double logE)
189{
190 size_t nIso = anElement->GetNumberOfIsotopes();
191 const G4Isotope* iso = anElement->GetIsotope(0);
192
193 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
194 if(1 == nIso) { return iso; }
195
196 // more than 1 isotope
197 G4int Z = anElement->GetZasInt();
198
199 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
201 G4double sum = 0.0;
202
203 // is there isotope wise cross section?
204 size_t j;
205 if(amax[Z] == amin[Z] || Z >= MAXZCAPTURE) {
206 for (j = 0; j<nIso; ++j) {
207 sum += abundVector[j];
208 if(q <= sum) {
209 iso = anElement->GetIsotope(j);
210 break;
211 }
212 }
213 return iso;
214 }
215 size_t nn = temp.size();
216 if(nn < nIso) { temp.resize(nIso, 0.); }
217
218 for (j=0; j<nIso; ++j) {
219 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
220 anElement->GetIsotope(j)->GetN());
221 temp[j] = sum;
222 }
223 sum *= q;
224 for (j = 0; j<nIso; ++j) {
225 if(temp[j] >= sum) {
226 iso = anElement->GetIsotope(j);
227 break;
228 }
229 }
230 return iso;
231}
232
233void
235{
236 if(verboseLevel > 0){
237 G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
238 << p.GetParticleName() << G4endl;
239 }
240 if(p.GetParticleName() != "neutron") {
242 ed << p.GetParticleName() << " is a wrong particle type -"
243 << " only neutron is allowed";
244 G4Exception("G4NeutronCaptureXS::BuildPhysicsTable(..)","had012",
245 FatalException, ed, "");
246 return;
247 }
248
249 if(nullptr == data) {
250#ifdef G4MULTITHREADED
251 G4MUTEXLOCK(&neutronCaptureXSMutex);
252 if(nullptr == data) {
253#endif
254 isMaster = true;
255 data = new G4ElementData();
256 data->SetName("NeutronCapture");
258#ifdef G4MULTITHREADED
259 }
260 G4MUTEXUNLOCK(&neutronCaptureXSMutex);
261#endif
262 }
263
264 // it is possible re-initialisation for the second run
266 if(isMaster) {
267
268 // Access to elements
269 for ( auto & elm : *table ) {
270 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZCAPTURE-1) );
271 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
272 }
273 }
274 // prepare isotope selection
275 size_t nIso = temp.size();
276 for ( auto & elm : *table ) {
277 size_t n = elm->GetNumberOfIsotopes();
278 if(n > nIso) { nIso = n; }
279 }
280 temp.resize(nIso, 0.0);
281}
282
284{
285 // check environment variable
286 // build the complete string identifying the file with the data set
287 if(gDataDirectory.empty()) {
288 char* path = std::getenv("G4PARTICLEXSDATA");
289 if (nullptr != path) {
290 std::ostringstream ost;
291 ost << path << "/neutron/cap";
292 gDataDirectory = ost.str();
293 } else {
294 G4Exception("G4NeutronCaptureXS::Initialise(..)","had013",
296 "Environment variable G4PARTICLEXSDATA is not defined");
297 }
298 }
299 return gDataDirectory;
300}
301
303{
304#ifdef G4MULTITHREADED
305 G4MUTEXLOCK(&neutronCaptureXSMutex);
306 if(nullptr == data->GetElementData(Z)) {
307#endif
308 Initialise(Z);
309#ifdef G4MULTITHREADED
310 }
311 G4MUTEXUNLOCK(&neutronCaptureXSMutex);
312#endif
313}
314
316{
317 if(nullptr != data->GetElementData(Z)) { return; }
318
319 // upload element data
320 std::ostringstream ost;
321 ost << FindDirectoryPath() << Z ;
322 G4PhysicsVector* v = RetrieveVector(ost, true);
324
325 // upload isotope data
326 if(amin[Z] < amax[Z]) {
327 G4int nmax = amax[Z] - amin[Z] + 1;
329
330 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
331 std::ostringstream ost1;
332 ost1 << gDataDirectory << Z << "_" << A;
333 v = RetrieveVector(ost1, false);
334 data->AddComponent(Z, A, v);
335 }
336 }
337}
338
340G4NeutronCaptureXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
341{
342 G4PhysicsLogVector* v = nullptr;
343 std::ifstream filein(ost.str().c_str());
344 if (!filein.is_open()) {
345 if(warn) {
347 ed << "Data file <" << ost.str().c_str()
348 << "> is not opened!";
349 G4Exception("G4NeutronCaptureXS::RetrieveVector(..)","had014",
350 FatalException, ed, "Check G4PARTICLEXSDATA");
351 }
352 } else {
353 if(verboseLevel > 1) {
354 G4cout << "File " << ost.str()
355 << " is opened by G4NeutronCaptureXS" << G4endl;
356 }
357 // retrieve data from DB
358 v = new G4PhysicsLogVector();
359 if(!v->Retrieve(filein, true)) {
361 ed << "Data file <" << ost.str().c_str()
362 << "> is not retrieved!";
363 G4Exception("G4NeutronCaptureXS::RetrieveVector(..)","had015",
364 FatalException, ed, "Check G4PARTICLEXSDATA");
365 }
366 }
367 return v;
368}
static const G4double e1[44]
static const G4double emax
std::vector< G4Element * > G4ElementTable
@ 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
static const G4int amax[]
static const G4int amin[]
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double barn
Definition: G4SIunits.hh:85
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double MeV
Definition: G4SIunits.hh:200
#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
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:397
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
void CrossSectionDescription(std::ostream &) const final
static G4String gDataDirectory
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
static const G4int MAXZCAPTURE
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)
void BuildPhysicsTable(const G4ParticleDefinition &) final
void InitialiseOnFly(G4int Z)
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) final
std::vector< G4double > temp
const G4String & FindDirectoryPath()
G4PhysicsVector * RetrieveVector(std::ostringstream &in, G4bool warn)
const G4PhysicsVector * GetPhysicsVector(G4int Z)
static G4ElementData * data
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
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double MeV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments