Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronInelasticXS.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: G4NeutronInelasticXS
32//
33// Author Ivantchenko, Geant4, 3-Aug-09
34//
35
37#include "G4Neutron.hh"
38#include "G4DynamicParticle.hh"
39#include "G4ElementTable.hh"
40#include "G4Material.hh"
41#include "G4Element.hh"
42#include "G4PhysicsLogVector.hh"
46#include "Randomize.hh"
47#include "G4Neutron.hh"
48#include "G4SystemOfUnits.hh"
49#include "G4IsotopeList.hh"
50#include "G4AutoLock.hh"
51
52#include <fstream>
53#include <sstream>
54#include <thread>
55
56G4double G4NeutronInelasticXS::coeff[] = {1.0};
57G4ElementData* G4NeutronInelasticXS::data = nullptr;
58G4String G4NeutronInelasticXS::gDataDirectory = "";
59
60static std::once_flag applyOnce;
61
62namespace
63{
64 G4Mutex nInelasticXSMutex = G4MUTEX_INITIALIZER;
65}
66
68 : G4VCrossSectionDataSet(Default_Name()),
69 neutron(G4Neutron::Neutron()),
70 elimit(20*CLHEP::MeV)
71{
72 verboseLevel = 0;
73 if (verboseLevel > 0){
74 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
75 << MAXZINEL << G4endl;
76 }
77 if (nullptr == data) {
78 data = new G4ElementData(MAXZINEL);
79 data->SetName("nInelastic");
80 FindDirectoryPath();
81 }
82 ggXsection =
84 if(ggXsection == nullptr)
85 ggXsection = new G4ComponentGGHadronNucleusXsc();
87}
88
89void G4NeutronInelasticXS::CrossSectionDescription(std::ostream& outFile) const
90{
91 outFile << "G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
92 << "cross section on nuclei using data from the high precision\n"
93 << "neutron database. These data are simplified and smoothed over\n"
94 << "the resonance region in order to reduce CPU time.\n"
95 << "For high energy Glauber-Gribov cross section model is used\n";
96}
97
98G4bool
100 G4int, const G4Material*)
101{
102 return true;
103}
104
105G4bool
107 G4int, G4int,
108 const G4Element*, const G4Material*)
109{
110 return true;
111}
112
115 G4int Z, const G4Material*)
116{
117 return ElementCrossSection(aParticle->GetKineticEnergy(),
118 aParticle->GetLogKineticEnergy(), Z);
119}
120
124 const G4Element* elm,
125 const G4Material*)
126{
127 return ElementCrossSection(ekin, loge, elm->GetZasInt());
128}
129
132{
133 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
134 auto pv = GetPhysicsVector(Z);
135
136 G4double xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
137 : coeff[Z]*ggXsection->GetInelasticElementCrossSection(neutron, ekin,
138 Z, aeff[Z]);
139
140#ifdef G4VERBOSE
141 if(verboseLevel > 1) {
142 G4cout << "G4NeutronInelasticXS::ElementCrossSection Z= " << Z
143 << " Ekin(MeV)= " << ekin/CLHEP::MeV
144 << ", ElmXSinel(b)= " << xs/CLHEP::barn
145 << G4endl;
146 }
147#endif
148 return xs;
149}
150
154 G4int Z, G4int A,
155 const G4Isotope*, const G4Element*,
156 const G4Material*)
157{
158 return IsoCrossSection(ekin, loge, Z, A);
159}
160
163 G4int Z, G4int A,
164 const G4Isotope*, const G4Element*,
165 const G4Material*)
166{
167 return IsoCrossSection(aParticle->GetKineticEnergy(),
168 aParticle->GetLogKineticEnergy(), Z, A);
169}
170
173 G4int ZZ, G4int A)
174{
175 G4double xs = 0.0;
176 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
177
178 auto pv = GetPhysicsVector(Z);
179 /*
180 G4cout << "G4NeutronInelasticXS::IsoCrossSection Z= "
181 << Z << " A= " << A << G4endl;
182 G4cout << " Amin= " << amin[Z] << " Amax= " << amax[Z]
183 << " E(MeV)= " << ekin << " Ncomp="
184 << data->GetNumberOfComponents(Z) << G4endl;
185 */
186
187 // compute isotope cross section if applicable
188 if (ekin <= elimit && data->GetNumberOfComponents(Z) > 0) {
189 auto pviso = data->GetComponentDataByID(Z, A);
190 if (nullptr != pviso) {
191 xs = pviso->LogVectorValue(ekin, logekin);
192#ifdef G4VERBOSE
193 if(verboseLevel > 1) {
194 G4cout << "G4NeutronInelasticXS::IsoXS: Ekin(MeV)= "
195 << ekin/CLHEP::MeV
196 << " xs(b)= " << xs/CLHEP::barn
197 << " Z= " << Z << " A= " << A << G4endl;
198 }
199#endif
200 return xs;
201 }
202 }
203
204 // use element x-section
205 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, logekin)
206 : coeff[Z]*ggXsection->GetInelasticElementCrossSection(neutron, ekin,
207 Z, aeff[Z]);
208 xs *= A/aeff[Z];
209#ifdef G4VERBOSE
210 if(verboseLevel > 1) {
211 G4cout << "G4NeutronInelasticXS::IsoXS: Z= " << Z << " A= " << A
212 << " Ekin(MeV)= " << ekin/CLHEP::MeV
213 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
214 }
215#endif
216 return xs;
217}
218
220 const G4Element* anElement, G4double kinEnergy, G4double logE)
221{
222 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
223 const G4Isotope* iso = anElement->GetIsotope(0);
224 if(1 == nIso) { return iso; }
225
226 // more than 1 isotope
227 G4int Z = anElement->GetZasInt();
228 if (nullptr == data->GetElementData(Z)) { InitialiseOnFly(Z); }
229
230 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
232 G4double sum = 0.0;
233 G4int j;
234
235 // isotope wise cross section not available
236 if (Z >= MAXZINEL || 0 == data->GetNumberOfComponents(Z)) {
237 for (j=0; j<nIso; ++j) {
238 sum += abundVector[j];
239 if(q <= sum) {
240 iso = anElement->GetIsotope(j);
241 break;
242 }
243 }
244 return iso;
245 }
246
247 // use isotope cross sections
248 G4int nn = (G4int)temp.size();
249 if(nn < nIso) { temp.resize(nIso, 0.); }
250
251 for (j=0; j<nIso; ++j) {
252 // G4cout << j << "-th isotope " << anElement->GetIsotope(j)->GetN()
253 // << " abund= " << abundVector[j] << G4endl;
254 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
255 anElement->GetIsotope(j)->GetN());
256 temp[j] = sum;
257 }
258 sum *= q;
259 for (j = 0; j<nIso; ++j) {
260 if (temp[j] >= sum) {
261 iso = anElement->GetIsotope(j);
262 break;
263 }
264 }
265 return iso;
266}
267
268void
270{
271 if (verboseLevel > 0) {
272 G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
273 << p.GetParticleName() << G4endl;
274 }
275 if (p.GetParticleName() != "neutron") {
277 ed << p.GetParticleName() << " is a wrong particle type -"
278 << " only neutron is allowed";
279 G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
280 FatalException, ed, "");
281 return;
282 }
283 // it is possible re-initialisation for the new run
285
286 // initialise static tables only once
287 std::call_once(applyOnce, [this]() { isInitializer = true; });
288
289 if (isInitializer) {
290 G4AutoLock l(&nInelasticXSMutex);
291
292 // Upload data for elements used in geometry
293 for ( auto const & elm : *table ) {
294 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINEL-1) );
295 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
296 }
297 l.unlock();
298 }
299 // prepare isotope selection
300 std::size_t nIso = temp.size();
301 for ( auto const & elm : *table ) {
302 std::size_t n = elm->GetNumberOfIsotopes();
303 if (n > nIso) { nIso = n; }
304 }
305 temp.resize(nIso, 0.0);
306}
307
308const G4String& G4NeutronInelasticXS::FindDirectoryPath()
309{
310 // build the complete string identifying the file with the data set
311 if (gDataDirectory.empty()) {
312 std::ostringstream ost;
313 ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/neutron/inel";
314 gDataDirectory = ost.str();
315 }
316 return gDataDirectory;
317}
318
319void G4NeutronInelasticXS::InitialiseOnFly(G4int Z)
320{
321 G4AutoLock l(&nInelasticXSMutex);
322 Initialise(Z);
323 l.unlock();
324}
325
326void G4NeutronInelasticXS::Initialise(G4int Z)
327{
328 if (nullptr != data->GetElementData(Z)) { return; }
329
330 // upload element data
331 std::ostringstream ost;
332 ost << FindDirectoryPath() << Z;
333 G4PhysicsVector* v = RetrieveVector(ost, true);
334 data->InitialiseForElement(Z, v);
335 if (verboseLevel > 1) {
336 G4cout << "G4NeutronInelasticXS::Initialise for Z= " << Z
337 << " A= " << aeff[Z] << " Amin= " << amin[Z]
338 << " Amax= " << amax[Z] << G4endl;
339 }
340 // upload isotope data
341 G4bool noComp = true;
342 if (amin[Z] < amax[Z]) {
343
344 for (G4int A=amin[Z]; A<=amax[Z]; ++A) {
345 std::ostringstream ost1;
346 ost1 << gDataDirectory << Z << "_" << A;
347 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
348 if (nullptr != v1) {
349 if (noComp) {
350 G4int nmax = amax[Z] - A + 1;
351 data->InitialiseForComponent(Z, nmax);
352 noComp = false;
353 }
354 data->AddComponent(Z, A, v1);
355 }
356 }
357 }
358 // no components case
359 if (noComp) { data->InitialiseForComponent(Z, 0); }
360
361 // smooth transition
362 G4double sig1 = (*v)[v->GetVectorLength()-1];
363 G4double ehigh= v->GetMaxEnergy();
364 G4double sig2 = ggXsection->GetInelasticElementCrossSection(neutron,
365 ehigh, Z, aeff[Z]);
366 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
367}
368
370G4NeutronInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
371{
372 G4PhysicsLogVector* v = nullptr;
373 std::ifstream filein(ost.str().c_str());
374 if (!filein.is_open()) {
375 if(warn) {
377 ed << "Data file <" << ost.str().c_str()
378 << "> is not opened!";
379 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had014",
380 FatalException, ed, "Check G4PARTICLEXSDATA");
381 }
382 } else {
383 if(verboseLevel > 1) {
384 G4cout << "File " << ost.str()
385 << " is opened by G4NeutronInelasticXS" << G4endl;
386 }
387 // retrieve data from DB
388 v = new G4PhysicsLogVector();
389 if(!v->Retrieve(filein, true)) {
391 ed << "Data file <" << ost.str().c_str()
392 << "> is not retrieved!";
393 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had015",
394 FatalException, ed, "Check G4PARTICLEXSDATA");
395 }
396 }
397 return v;
398}
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#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
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
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 ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
void CrossSectionDescription(std::ostream &) const final
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) final
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
const G4String & GetParticleName() const
G4double GetMaxEnergy() const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
std::size_t GetVectorLength() const
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
void SetForAllAtomsAndEnergies(G4bool val)