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