57G4String G4NeutronCaptureXS::gDataDirectory =
"";
59static std::once_flag applyOnce;
64 const G4int MAXZCAPTURE = 92;
73 G4cout <<
"G4NeutronCaptureXS::G4NeutronCaptureXS: Initialise for Z < "
76 logElimit =
G4Log(elimit);
77 if (
nullptr == data) {
79 data->SetName(
"nCapture");
81 dataR->SetName(
"nRCapture");
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";
139 G4int Z = std::min(ZZ, MAXZCAPTURE);
151 if (fRfilesEnabled) {
152 auto pv = GetPhysicsVectorR(Z);
153 if (
nullptr != pv && ekin < cap_max_r_e[Z]) {
155 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
156 : (*pv)[0]*std::sqrt(e0/ekin);
162 auto pv = GetPhysicsVector(Z);
164 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
165 : (*pv)[0]*std::sqrt(e0/ekin);
170 G4cout <<
"Ekin= " << ekin/CLHEP::MeV
171 <<
" ElmXScap(b)= " << xs/CLHEP::barn <<
G4endl;
202 if (eKin > emax) {
return xs; }
204 G4int Z = std::min(ZZ, MAXZCAPTURE);
215 if (fRfilesEnabled) {
216 auto pv = GetPhysicsVectorR(Z);
217 if (
nullptr != pv && ekin < cap_max_r_e[Z]) {
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);
231 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
232 : (*pv)[0]*std::sqrt(e0/ekin);
239 auto pv = GetPhysicsVector(Z);
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);
253 xs = (ekin >= e0) ? pv->LogVectorValue(ekin, logEkin)
254 : (*pv)[0]*std::sqrt(e0/ekin);
259 G4cout <<
"G4NeutronCaptureXS::IsoXS: Ekin(MeV)= " << ekin/MeV
260 <<
" xs(b)= " << xs/CLHEP::barn
261 <<
" Z= " << Z <<
" A= " <<
A <<
" no iso XS" <<
G4endl;
275 if(1 == nIso) {
return iso; }
279 if (
nullptr == data->GetElementData(Z)) { InitialiseOnFly(Z); }
287 if (Z > MAXZCAPTURE || 0 == data->GetNumberOfComponents(Z)) {
288 for (j = 0; j<nIso; ++j) {
289 sum += abundVector[j];
298 if (nn < nIso) { temp.resize(nIso, 0.); }
300 for (j=0; j<nIso; ++j) {
306 for (j = 0; j<nIso; ++j) {
307 if (temp[j] >= sum) {
319 G4cout <<
"G4NeutronCaptureXS::BuildPhysicsTable for "
325 <<
" only neutron is allowed";
326 G4Exception(
"G4NeutronCaptureXS::BuildPhysicsTable(..)",
"had012",
337 std::call_once(applyOnce, [
this]() { isInitializer =
true; });
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); }
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; }
355 temp.resize(nIso, 0.0);
358const G4String& G4NeutronCaptureXS::FindDirectoryPath()
361 if(gDataDirectory.empty()) {
362 std::ostringstream ost;
364 gDataDirectory = ost.str();
366 return gDataDirectory;
369void G4NeutronCaptureXS::InitialiseOnFly(
G4int Z)
376void G4NeutronCaptureXS::Initialise(
G4int Z)
378 if (
nullptr != data->GetElementData(Z)) {
return; }
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);
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);
403 G4int nmax = amax[Z] -
A + 1;
404 data->InitialiseForComponent(Z, nmax);
407 data->AddComponent(Z,
A, v1);
410 std::ostringstream ost2;
411 ost2 << gDataDirectory <<
"Rcap" << Z <<
"_" <<
A;
412 G4PhysicsVector* v2 = RetrieveVector(ost2,
false);
415 G4int nmax = amax[Z] -
A + 1;
416 dataR->InitialiseForComponent(Z, nmax);
419 dataR->AddComponent(Z,
A, v2);
425 if (noComp) { data->InitialiseForComponent(Z, 0); }
426 if (noCompR &&
nullptr != vr) { dataR->InitialiseForComponent(Z, 0); }
430G4NeutronCaptureXS::RetrieveVector(std::ostringstream& ost,
G4bool warn)
432 G4PhysicsLogVector* v =
nullptr;
433 std::ifstream filein(ost.str().c_str());
434 if (!filein.is_open()) {
437 ed <<
"Data file <" << ost.str().c_str()
438 <<
"> is not opened!";
439 G4Exception(
"G4NeutronCaptureXS::RetrieveVector(..)",
"had014",
444 G4cout <<
"File " << ost.str()
445 <<
" is opened by G4NeutronCaptureXS" <<
G4endl;
448 v =
new G4PhysicsLogVector();
451 ed <<
"Data file <" << ost.str().c_str()
452 <<
"> is not retrieved!";
453 G4Exception(
"G4NeutronCaptureXS::RetrieveVector(..)",
"had015",
G4TemplateAutoLock< G4Mutex > G4AutoLock
std::vector< G4Element * > G4ElementTable
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4double G4Log(G4double x)
#define G4MUTEX_INITIALIZER
G4GLOB_DLL std::ostream G4cout
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4double * GetRelativeAbundanceVector() const
std::size_t GetNumberOfIsotopes() const
const G4Isotope * GetIsotope(G4int iso) const
static const G4ElementTable * GetElementTable()
static G4HadronicParameters * Instance()
G4bool UseRFilesForXS() const
const G4String & GetDirPARTICLEXS() const
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="")