Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronElasticXS.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: G4NeutronElasticXS
32//
33// Author Ivantchenko, Geant4, 3-Aug-09
34//
35// Modifications:
36//
37
38#include "G4NeutronElasticXS.hh"
39#include "G4Neutron.hh"
40#include "G4DynamicParticle.hh"
41#include "G4ElementTable.hh"
42#include "G4Material.hh"
43#include "G4Element.hh"
44#include "G4PhysicsLogVector.hh"
48#include "Randomize.hh"
49#include "G4SystemOfUnits.hh"
50#include "G4IsotopeList.hh"
51#include "G4AutoLock.hh"
52
53#include <fstream>
54#include <sstream>
55
56G4ElementData* G4NeutronElasticXS::data = nullptr;
57G4ElementData* G4NeutronElasticXS::dataR = nullptr;
58G4double G4NeutronElasticXS::coeff[] = {1.0};
59G4String G4NeutronElasticXS::gDataDirectory = "";
60
61static std::once_flag applyOnce;
62
63namespace
64{
65 G4Mutex nElasticXSMutex = G4MUTEX_INITIALIZER;
66}
67
70 neutron(G4Neutron::Neutron())
71{
72 if (verboseLevel > 0) {
73 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
74 << MAXZEL << G4endl;
75 }
76 ggXsection =
78 if (ggXsection == nullptr)
79 ggXsection = new G4ComponentGGHadronNucleusXsc();
81 if (nullptr == data) {
82 data = new G4ElementData(MAXZEL);
83 data->SetName("nElastic");
84 dataR = new G4ElementData(MAXZEL);
85 dataR->SetName("nRElastic");
86 FindDirectoryPath();
87 }
88}
89
90void G4NeutronElasticXS::CrossSectionDescription(std::ostream& outFile) const
91{
92 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
93 << "cross section on nuclei using data from the high precision\n"
94 << "neutron database. These data are simplified and smoothed over\n"
95 << "the resonance region in order to reduce CPU time.\n"
96 << "For high energies Glauber-Gribiv cross section is used.\n";
97}
98
99G4bool
101 G4int, const G4Material*)
102{
103 return true;
104}
105
107 G4int, G4int,
108 const G4Element*, const G4Material*)
109{
110 return false;
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 = std::min(ZZ, MAXZEL-1);
134 G4double xs;
135 G4bool done{false};
136
137 // data from the resonance region
138 if (fRfilesEnabled) {
139 auto pv = GetPhysicsVectorR(Z);
140 if (nullptr != pv && ekin < el_max_r_e[Z]) {
141 xs = pv->LogVectorValue(ekin, loge);
142 done = true;
143 }
144 }
145 // data above the resonance region
146 if (!done) {
147 auto pv = GetPhysicsVector(Z);
148 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
149 : coeff[Z]*ggXsection->GetElasticElementCrossSection(neutron, ekin,
150 Z, aeff[Z]);
151 }
152
153#ifdef G4VERBOSE
154 if (verboseLevel > 1) {
155 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
156 << ", nElmXSel(b)= " << xs/CLHEP::barn
157 << G4endl;
158 }
159#endif
160 return xs;
161}
162
166 G4int ZZ, G4int A,
167 const G4Isotope*, const G4Element*,
168 const G4Material*)
169{
170 G4int Z = std::min(ZZ, MAXZEL-1);
171 return ElementCrossSection(ekin, loge, Z)*A/aeff[Z];
172}
173
176 G4int ZZ, G4int A,
177 const G4Isotope*, const G4Element*,
178 const G4Material*)
179{
180 G4int Z = std::min(ZZ, MAXZEL-1);
181 return ElementCrossSection(aParticle->GetKineticEnergy(),
182 aParticle->GetLogKineticEnergy(), Z)*A/aeff[Z];
183
184}
185
187 const G4Element* anElement, G4double, G4double)
188{
189 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
190 const G4Isotope* iso = anElement->GetIsotope(0);
191
192 if(1 == nIso) { return iso; }
193
194 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
196 G4double sum = 0.0;
197
198 // isotope wise cross section not used
199 for (G4int j=0; j<nIso; ++j) {
200 sum += abundVector[j];
201 if(q <= sum) {
202 iso = anElement->GetIsotope(j);
203 break;
204 }
205 }
206 return iso;
207}
208
209void
211{
212 if(verboseLevel > 0){
213 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
214 << p.GetParticleName() << G4endl;
215 }
216 if(p.GetParticleName() != "neutron") {
218 ed << p.GetParticleName() << " is a wrong particle type -"
219 << " only neutron is allowed";
220 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
221 FatalException, ed, "");
222 return;
223 }
224
226
227 // initialise static tables only once
228 std::call_once(applyOnce, [this]() { isInitializer = true; });
229
230 if (isInitializer) {
231 G4AutoLock l(&nElasticXSMutex);
232
233 // Access to elements
235 for ( auto const & elm : *table ) {
236 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZEL-1) );
237 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
238 }
239 l.unlock();
240 }
241}
242
243const G4String& G4NeutronElasticXS::FindDirectoryPath()
244{
245 // build the complete string identifying the file with the data set
246 if (gDataDirectory.empty()) {
247 std::ostringstream ost;
248 ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/neutron/";
249 gDataDirectory = ost.str();
250 }
251 return gDataDirectory;
252}
253
254void G4NeutronElasticXS::InitialiseOnFly(G4int Z)
255{
256 G4AutoLock l(&nElasticXSMutex);
257 Initialise(Z);
258 l.unlock();
259}
260
261void G4NeutronElasticXS::Initialise(G4int Z)
262{
263 if (nullptr != data->GetElementData(Z)) { return; }
264
265 // upload element data
266 std::ostringstream ost;
267 ost << FindDirectoryPath() << "el" << Z;
268 G4PhysicsVector* v = RetrieveVector(ost, true);
269 data->InitialiseForElement(Z, v);
270
271 G4PhysicsVector* vr = nullptr;
272 if (fRfilesEnabled) {
273 std::ostringstream ostr;
274 ostr << FindDirectoryPath() << "Rel" << Z;
275 vr = RetrieveVector(ostr, false);
276 dataR->InitialiseForElement(Z, vr);
277 }
278
279 // smooth transition
280 G4double sig1 = (*v)[v->GetVectorLength()-1];
281 G4double ehigh = v->GetMaxEnergy();
282 G4double sig2 =
283 ggXsection->GetElasticElementCrossSection(neutron, ehigh, Z, aeff[Z]);
284 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
285}
286
288G4NeutronElasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
289{
290 G4PhysicsLogVector* v = nullptr;
291 std::ifstream filein(ost.str().c_str());
292 if (!filein.is_open()) {
293 if (warn) {
295 ed << "Data file <" << ost.str().c_str()
296 << "> is not opened!";
297 G4Exception("G4NeutronElasticXS::RetrieveVector(..)","had014",
298 FatalException, ed, "Check G4PARTICLEXSDATA");
299 }
300 } else {
301 if (verboseLevel > 1) {
302 G4cout << "File " << ost.str()
303 << " is opened by G4NeutronElasticXS" << G4endl;
304 }
305 // retrieve data from DB
306 v = new G4PhysicsLogVector();
307 if (!v->Retrieve(filein, true)) {
309 ed << "Data file <" << ost.str().c_str()
310 << "> is not retrieved!";
311 G4Exception("G4NeutronElasticXS::RetrieveVector(..)","had015",
312 FatalException, ed, "Check G4PARTICLEXSDATA");
313 }
314 }
315 return v;
316}
G4TemplateAutoLock< G4Mutex > G4AutoLock
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
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
G4int GetZasInt() const
Definition G4Element.hh:120
static G4HadronicParameters * Instance()
const G4String & GetDirPARTICLEXS() const
G4bool IsElementApplicable(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
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
static const char * Default_Name()
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
void CrossSectionDescription(std::ostream &) const final
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) final
const G4String & GetParticleName() const
G4double GetMaxEnergy() const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
std::size_t GetVectorLength() const
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)