Geant4-11
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
51#include <fstream>
52#include <sstream>
53
57
58#ifdef G4MULTITHREADED
59 G4Mutex G4NeutronElasticXS::neutronElasticXSMutex = G4MUTEX_INITIALIZER;
60#endif
61
63 : G4VCrossSectionDataSet(Default_Name()),
65{
66 // verboseLevel = 0;
67 if(verboseLevel > 0){
68 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
69 << MAXZEL << G4endl;
70 }
74}
75
77{
78 if(isMaster) {
79 for(G4int i=0; i<MAXZEL; ++i) {
80 delete data[i];
81 data[i] = nullptr;
82 }
83 }
84}
85
86void G4NeutronElasticXS::CrossSectionDescription(std::ostream& outFile) const
87{
88 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
89 << "cross section on nuclei using data from the high precision\n"
90 << "neutron database. These data are simplified and smoothed over\n"
91 << "the resonance region in order to reduce CPU time.\n"
92 << "For high energies Glauber-Gribiv cross section is used.\n";
93}
94
95G4bool
97 G4int, const G4Material*)
98{
99 return true;
100}
101
103 G4int, G4int,
104 const G4Element*, const G4Material*)
105{
106 return false;
107}
108
111 G4int ZZ, const G4Material*)
112{
113 G4double xs = 0.0;
114 G4double ekin = aParticle->GetKineticEnergy();
115
116 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
117
118 auto pv = GetPhysicsVector(Z);
119 if(pv == nullptr) { return xs; }
120 // G4cout << "G4NeutronElasticXS::GetCrossSection e= " << ekin
121 // << " Z= " << Z << G4endl;
122
123 if(ekin <= pv->Energy(1)) {
124 xs = (*pv)[1];
125 } else if(ekin <= pv->GetMaxEnergy()) {
126 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
127 } else {
129 ekin, Z, aeff[Z]);
130 }
131
132#ifdef G4VERBOSE
133 if(verboseLevel > 1) {
134 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
135 << ", nElmXSel(b)= " << xs/CLHEP::barn
136 << G4endl;
137 }
138#endif
139 return xs;
140}
141
143 const G4DynamicParticle* aParticle,
144 G4int Z, G4int A,
145 const G4Isotope*, const G4Element*,
146 const G4Material* mat)
147{
148 return GetElementCrossSection(aParticle, Z, mat) * A/aeff[Z];
149}
150
152 const G4Element* anElement, G4double, G4double)
153{
154 size_t nIso = anElement->GetNumberOfIsotopes();
155 const G4Isotope* iso = anElement->GetIsotope(0);
156
157 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
158 if(1 == nIso) { return iso; }
159
160 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
162 G4double sum = 0.0;
163 size_t j;
164
165 // isotope wise cross section not used
166 for (j=0; j<nIso; ++j) {
167 sum += abundVector[j];
168 if(q <= sum) {
169 iso = anElement->GetIsotope(j);
170 break;
171 }
172 }
173 return iso;
174}
175
176void
178{
179 if(verboseLevel > 0){
180 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
181 << p.GetParticleName() << G4endl;
182 }
183 if(p.GetParticleName() != "neutron") {
185 ed << p.GetParticleName() << " is a wrong particle type -"
186 << " only neutron is allowed";
187 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
188 FatalException, ed, "");
189 return;
190 }
191 if(0. == coeff[0]) {
192#ifdef G4MULTITHREADED
193 G4MUTEXLOCK(&neutronElasticXSMutex);
194 if(0. == coeff[0]) {
195#endif
196 coeff[0] = 1.0;
197 isMaster = true;
199#ifdef G4MULTITHREADED
200 }
201 G4MUTEXUNLOCK(&neutronElasticXSMutex);
202#endif
203 }
204
205 // it is possible re-initialisation for the second run
206 if(isMaster) {
207
208 // Access to elements
210 for ( auto & elm : *table ) {
211 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZEL-1) );
212 if ( nullptr == data[Z] ) { Initialise(Z); }
213 }
214 }
215}
216
218{
219 // check environment variable
220 // build the complete string identifying the file with the data set
221 if(gDataDirectory.empty()) {
222 char* path = std::getenv("G4PARTICLEXSDATA");
223 if (nullptr != path) {
224 std::ostringstream ost;
225 ost << path << "/neutron/el";
226 gDataDirectory = ost.str();
227 } else {
228 G4Exception("G4NeutronElasticXS::Initialise(..)","had013",
230 "Environment variable G4PARTICLEXSDATA is not defined");
231 }
232 }
233 return gDataDirectory;
234}
235
237{
238#ifdef G4MULTITHREADED
239 G4MUTEXLOCK(&neutronElasticXSMutex);
240 if(data[Z] == nullptr) {
241#endif
242 Initialise(Z);
243#ifdef G4MULTITHREADED
244 }
245 G4MUTEXUNLOCK(&neutronElasticXSMutex);
246#endif
247}
248
250{
251 if(data[Z] != nullptr) { return; }
252
253 // upload data from file
254 data[Z] = new G4PhysicsLogVector();
255
256 std::ostringstream ost;
257 ost << FindDirectoryPath() << Z ;
258 std::ifstream filein(ost.str().c_str());
259 if (!filein.is_open()) {
261 ed << "Data file <" << ost.str().c_str()
262 << "> is not opened!";
263 G4Exception("G4NeutronElasticXS::Initialise(..)","had014",
264 FatalException, ed, "Check G4PARTICLEXSDATA");
265 return;
266 }
267 if(verboseLevel > 1) {
268 G4cout << "file " << ost.str()
269 << " is opened by G4NeutronElasticXS" << G4endl;
270 }
271
272 // retrieve data from DB
273 if(!data[Z]->Retrieve(filein, true)) {
275 ed << "Data file <" << ost.str().c_str()
276 << "> is not retrieved!";
277 G4Exception("G4NeutronElasticXS::Initialise(..)","had015",
278 FatalException, ed, "Check G4PARTICLEXSDATA");
279 return;
280 }
281 // smooth transition
282 G4double sig1 = (*(data[Z]))[data[Z]->GetVectorLength()-1];
283 G4double ehigh = data[Z]->GetMaxEnergy();
285 ehigh, Z, aeff[Z]);
286 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
287}
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static const G4double aeff[]
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
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:397
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
void InitialiseOnFly(G4int Z)
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
G4PhysicsVector * GetPhysicsVector(G4int Z)
void BuildPhysicsTable(const G4ParticleDefinition &) final
const G4String & FindDirectoryPath()
static G4double coeff[MAXZEL]
static G4PhysicsVector * data[MAXZEL]
const G4ParticleDefinition * neutron
static const G4int MAXZEL
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
G4VComponentCrossSection * ggXsection
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
static G4String gDataDirectory
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)
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double MeV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments