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