Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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"
47#include "Randomize.hh"
48#include "G4SystemOfUnits.hh"
49#include "G4IsotopeList.hh"
50#include "G4AutoLock.hh"
51
52#include <fstream>
53#include <sstream>
54
55G4PhysicsVector* G4NeutronElasticXS::data[] = {nullptr};
56G4double G4NeutronElasticXS::coeff[] = {0.0};
57G4String G4NeutronElasticXS::gDataDirectory = "";
58
59namespace
60{
61 G4Mutex nElasticXSMutex = G4MUTEX_INITIALIZER;
62}
63
65 : G4VCrossSectionDataSet(Default_Name()),
66 neutron(G4Neutron::Neutron())
67{
68 // verboseLevel = 0;
69 if(verboseLevel > 0){
70 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
71 << MAXZEL << G4endl;
72 }
74 if(ggXsection == nullptr) ggXsection = new G4ComponentGGHadronNucleusXsc();
76}
77
79{
80 if(isMaster) {
81 for(G4int i=0; i<MAXZEL; ++i) {
82 delete data[i];
83 data[i] = nullptr;
84 }
85 }
86}
87
88void G4NeutronElasticXS::CrossSectionDescription(std::ostream& outFile) const
89{
90 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
91 << "cross section on nuclei using data from the high precision\n"
92 << "neutron database. These data are simplified and smoothed over\n"
93 << "the resonance region in order to reduce CPU time.\n"
94 << "For high energies Glauber-Gribiv cross section is used.\n";
95}
96
97G4bool
99 G4int, const G4Material*)
100{
101 return true;
102}
103
105 G4int, G4int,
106 const G4Element*, const G4Material*)
107{
108 return false;
109}
110
113 G4int Z, const G4Material*)
114{
115 return ElementCrossSection(aParticle->GetKineticEnergy(), aParticle->GetLogKineticEnergy(), Z);
116}
117
121 const G4Element* elm,
122 const G4Material*)
123{
124 return ElementCrossSection(ekin, loge, elm->GetZasInt());
125}
126
128{
129 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
130 auto pv = GetPhysicsVector(Z);
131
132 G4double xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
133 : coeff[Z]*ggXsection->GetElasticElementCrossSection(neutron, ekin,
134 Z, aeff[Z]);
135
136#ifdef G4VERBOSE
137 if(verboseLevel > 1) {
138 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
139 << ", nElmXSel(b)= " << xs/CLHEP::barn
140 << G4endl;
141 }
142#endif
143 return xs;
144}
145
149 G4int Z, G4int A,
150 const G4Isotope*, const G4Element*,
151 const G4Material*)
152{
153 return ElementCrossSection(ekin, loge, Z)*A/aeff[Z];
154}
155
158 G4int Z, G4int A,
159 const G4Isotope*, const G4Element*,
160 const G4Material*)
161{
162 return ElementCrossSection(aParticle->GetKineticEnergy(),
163 aParticle->GetLogKineticEnergy(), Z)*A/aeff[Z];
164
165}
166
168 const G4Element* anElement, G4double, G4double)
169{
170 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
171 const G4Isotope* iso = anElement->GetIsotope(0);
172
173 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
174 if(1 == nIso) { return iso; }
175
176 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
178 G4double sum = 0.0;
179
180 // isotope wise cross section not used
181 for (G4int j=0; j<nIso; ++j) {
182 sum += abundVector[j];
183 if(q <= sum) {
184 iso = anElement->GetIsotope(j);
185 break;
186 }
187 }
188 return iso;
189}
190
191void
193{
194 if(verboseLevel > 0){
195 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
196 << p.GetParticleName() << G4endl;
197 }
198 if(p.GetParticleName() != "neutron") {
200 ed << p.GetParticleName() << " is a wrong particle type -"
201 << " only neutron is allowed";
202 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
203 FatalException, ed, "");
204 return;
205 }
206 if(0. == coeff[0]) {
207 G4AutoLock l(&nElasticXSMutex);
208 if(0. == coeff[0]) {
209 coeff[0] = 1.0;
210 isMaster = true;
211 FindDirectoryPath();
212 }
213 l.unlock();
214 }
215
216 // it is possible re-initialisation for the second run
217 if(isMaster) {
218
219 // Access to elements
221 for ( auto & elm : *table ) {
222 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZEL-1) );
223 if ( nullptr == data[Z] ) { Initialise(Z); }
224 }
225 }
226}
227
228const G4String& G4NeutronElasticXS::FindDirectoryPath()
229{
230 // check environment variable
231 // build the complete string identifying the file with the data set
232 if(gDataDirectory.empty()) {
233 const char* path = G4FindDataDir("G4PARTICLEXSDATA");
234 if (nullptr != path) {
235 std::ostringstream ost;
236 ost << path << "/neutron/el";
237 gDataDirectory = ost.str();
238 } else {
239 G4Exception("G4NeutronElasticXS::Initialise(..)","had013",
241 "Environment variable G4PARTICLEXSDATA is not defined");
242 }
243 }
244 return gDataDirectory;
245}
246
247void G4NeutronElasticXS::InitialiseOnFly(G4int Z)
248{
249 if(nullptr == data[Z]) {
250 G4AutoLock l(&nElasticXSMutex);
251 if(nullptr == data[Z]) { Initialise(Z); }
252 l.unlock();
253 }
254}
255
256void G4NeutronElasticXS::Initialise(G4int Z)
257{
258 if(data[Z] != nullptr) { return; }
259
260 // upload data from file
261 data[Z] = new G4PhysicsLogVector();
262
263 std::ostringstream ost;
264 ost << FindDirectoryPath() << Z ;
265 std::ifstream filein(ost.str().c_str());
266 if (!filein.is_open()) {
268 ed << "Data file <" << ost.str().c_str()
269 << "> is not opened!";
270 G4Exception("G4NeutronElasticXS::Initialise(..)","had014",
271 FatalException, ed, "Check G4PARTICLEXSDATA");
272 return;
273 }
274 if(verboseLevel > 1) {
275 G4cout << "file " << ost.str()
276 << " is opened by G4NeutronElasticXS" << G4endl;
277 }
278
279 // retrieve data from DB
280 if(!data[Z]->Retrieve(filein, true)) {
282 ed << "Data file <" << ost.str().c_str()
283 << "> is not retrieved!";
284 G4Exception("G4NeutronElasticXS::Initialise(..)","had015",
285 FatalException, ed, "Check G4PARTICLEXSDATA");
286 return;
287 }
288 // smooth transition
289 G4double sig1 = (*(data[Z]))[data[Z]->GetVectorLength()-1];
290 G4double ehigh = data[Z]->GetMaxEnergy();
291 G4double sig2 = ggXsection->GetElasticElementCrossSection(neutron,
292 ehigh, Z, aeff[Z]);
293 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
294}
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
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
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
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
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
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
std::size_t GetVectorLength() const
G4double GetElasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
void SetForAllAtomsAndEnergies(G4bool val)