Geant4 11.2.2
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
56G4PhysicsVector* G4NeutronElasticXS::data[] = {nullptr};
57G4double G4NeutronElasticXS::coeff[] = {0.0};
58G4String G4NeutronElasticXS::gDataDirectory = "";
59G4bool G4NeutronElasticXS::fLock = true;
60
61namespace
62{
63 G4Mutex nElasticXSMutex = G4MUTEX_INITIALIZER;
64}
65
67 : G4VCrossSectionDataSet(Default_Name()),
68 neutron(G4Neutron::Neutron())
69{
70 // verboseLevel = 0;
71 if (verboseLevel > 0){
72 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
73 << MAXZEL << G4endl;
74 }
75 ggXsection =
77 if (ggXsection == nullptr)
78 ggXsection = new G4ComponentGGHadronNucleusXsc();
80}
81
83{
84 if (isFirst) {
85 for(G4int i=0; i<MAXZEL; ++i) {
86 delete data[i];
87 data[i] = nullptr;
88 }
89 }
90}
91
92void G4NeutronElasticXS::CrossSectionDescription(std::ostream& outFile) const
93{
94 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
95 << "cross section on nuclei using data from the high precision\n"
96 << "neutron database. These data are simplified and smoothed over\n"
97 << "the resonance region in order to reduce CPU time.\n"
98 << "For high energies Glauber-Gribiv cross section is used.\n";
99}
100
101G4bool
103 G4int, const G4Material*)
104{
105 return true;
106}
107
109 G4int, G4int,
110 const G4Element*, const G4Material*)
111{
112 return false;
113}
114
117 G4int Z, const G4Material*)
118{
119 return ElementCrossSection(aParticle->GetKineticEnergy(),
120 aParticle->GetLogKineticEnergy(), Z);
121}
122
126 const G4Element* elm,
127 const G4Material*)
128{
129 return ElementCrossSection(ekin, loge, elm->GetZasInt());
130}
131
133{
134 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
135 auto pv = GetPhysicsVector(Z);
136
137 G4double xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
138 : coeff[Z]*ggXsection->GetElasticElementCrossSection(neutron, ekin,
139 Z, aeff[Z]);
140
141#ifdef G4VERBOSE
142 if(verboseLevel > 1) {
143 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
144 << ", nElmXSel(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 ElementCrossSection(ekin, loge, Z)*A/aeff[Z];
159}
160
163 G4int Z, G4int A,
164 const G4Isotope*, const G4Element*,
165 const G4Material*)
166{
167 return ElementCrossSection(aParticle->GetKineticEnergy(),
168 aParticle->GetLogKineticEnergy(), Z)*A/aeff[Z];
169
170}
171
173 const G4Element* anElement, G4double, G4double)
174{
175 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
176 const G4Isotope* iso = anElement->GetIsotope(0);
177
178 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
179 if(1 == nIso) { return iso; }
180
181 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
183 G4double sum = 0.0;
184
185 // isotope wise cross section not used
186 for (G4int j=0; j<nIso; ++j) {
187 sum += abundVector[j];
188 if(q <= sum) {
189 iso = anElement->GetIsotope(j);
190 break;
191 }
192 }
193 return iso;
194}
195
196void
198{
199 if(verboseLevel > 0){
200 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
201 << p.GetParticleName() << G4endl;
202 }
203 if(p.GetParticleName() != "neutron") {
205 ed << p.GetParticleName() << " is a wrong particle type -"
206 << " only neutron is allowed";
207 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
208 FatalException, ed, "");
209 return;
210 }
211 if (fLock || isFirst) {
212 G4AutoLock l(&nElasticXSMutex);
213 if (fLock) {
214 isFirst = true;
215 fLock = false;
216 FindDirectoryPath();
217 }
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 l.unlock();
226 }
227}
228
229const G4String& G4NeutronElasticXS::FindDirectoryPath()
230{
231 // build the complete string identifying the file with the data set
232 if (gDataDirectory.empty()) {
233 std::ostringstream ost;
234 ost << G4HadronicParameters::Instance()->GetDirPARTICLEXS() << "/neutron/el";
235 gDataDirectory = ost.str();
236 }
237 return gDataDirectory;
238}
239
240void G4NeutronElasticXS::InitialiseOnFly(G4int Z)
241{
242 G4AutoLock l(&nElasticXSMutex);
243 Initialise(Z);
244 l.unlock();
245}
246
247void G4NeutronElasticXS::Initialise(G4int Z)
248{
249 if(data[Z] != nullptr) { return; }
250
251 // upload data from file
252 data[Z] = new G4PhysicsLogVector();
253
254 std::ostringstream ost;
255 ost << FindDirectoryPath() << Z ;
256 std::ifstream filein(ost.str().c_str());
257 if (!filein.is_open()) {
259 ed << "Data file <" << ost.str().c_str()
260 << "> is not opened!";
261 G4Exception("G4NeutronElasticXS::Initialise(..)","had014",
262 FatalException, ed, "Check G4PARTICLEXSDATA");
263 return;
264 }
265 if(verboseLevel > 1) {
266 G4cout << "file " << ost.str()
267 << " is opened by G4NeutronElasticXS" << G4endl;
268 }
269
270 // retrieve data from DB
271 if(!data[Z]->Retrieve(filein, true)) {
273 ed << "Data file <" << ost.str().c_str()
274 << "> is not retrieved!";
275 G4Exception("G4NeutronElasticXS::Initialise(..)","had015",
276 FatalException, ed, "Check G4PARTICLEXSDATA");
277 return;
278 }
279 // smooth transition
280 G4double sig1 = (*(data[Z]))[data[Z]->GetVectorLength()-1];
281 G4double ehigh = data[Z]->GetMaxEnergy();
282 G4double sig2 = ggXsection->GetElasticElementCrossSection(neutron,
283 ehigh, Z, aeff[Z]);
284 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
285}
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
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
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
G4double GetElasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
void SetForAllAtomsAndEnergies(G4bool val)