Geant4-11
G4NeutronInelasticXS.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: G4NeutronInelasticXS
32//
33// Author Ivantchenko, Geant4, 3-Aug-09
34//
35
37#include "G4Neutron.hh"
38#include "G4DynamicParticle.hh"
39#include "G4ElementTable.hh"
40#include "G4Material.hh"
41#include "G4Element.hh"
42#include "G4PhysicsLogVector.hh"
45#include "Randomize.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4IsotopeList.hh"
48
49#include <fstream>
50#include <sstream>
51
55
56#ifdef G4MULTITHREADED
57 G4Mutex G4NeutronInelasticXS::neutronInelasticXSMutex = G4MUTEX_INITIALIZER;
58#endif
59
61 : G4VCrossSectionDataSet(Default_Name()),
63{
64 verboseLevel = 0;
65 if(verboseLevel > 0){
66 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
67 << MAXZINEL << G4endl;
68 }
72}
73
75{
76 if(isMaster) { delete data; data = nullptr; }
77}
78
79void G4NeutronInelasticXS::CrossSectionDescription(std::ostream& outFile) const
80{
81 outFile << "G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
82 << "cross section on nuclei using data from the high precision\n"
83 << "neutron database. These data are simplified and smoothed over\n"
84 << "the resonance region in order to reduce CPU time.\n"
85 << "For high energy Glauber-Gribov cross section model is used\n";
86}
87
88G4bool
90 G4int, const G4Material*)
91{
92 return true;
93}
94
95G4bool
97 G4int, G4int,
98 const G4Element*, const G4Material*)
99{
100 return true;
101}
102
104 const G4DynamicParticle* aParticle,
105 G4int ZZ, const G4Material*)
106{
107 G4double xs = 0.0;
108 G4double ekin = aParticle->GetKineticEnergy();
109
110 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
111
112 auto pv = GetPhysicsVector(Z);
113 if(pv == nullptr) { return xs; }
114 // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin
115 // << " Z= " << Z << G4endl;
116
117 if(ekin <= pv->GetMaxEnergy()) {
118 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
119 } else {
121 ekin, Z, aeff[Z]);
122 }
123
124#ifdef G4VERBOSE
125 if(verboseLevel > 1) {
126 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
127 << ", ElmXSinel(b)= " << xs/CLHEP::barn
128 << G4endl;
129 }
130#endif
131 return xs;
132}
133
135 const G4DynamicParticle* aParticle,
136 G4int Z, G4int A,
137 const G4Isotope*, const G4Element*,
138 const G4Material*)
139{
140 return IsoCrossSection(aParticle->GetKineticEnergy(),
141 aParticle->GetLogKineticEnergy(), Z, A);
142}
143
146 G4int ZZ, G4int A)
147{
148 G4double xs = 0.0;
149 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
150
151 /*
152 G4cout << "IsoCrossSection Z= " << Z << " A= " << A
153 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
154 << " E(MeV)= " << ekin << G4endl;
155 */
156 auto pv = GetPhysicsVector(Z);
157 if(pv == nullptr) { return xs; }
158
159 // compute isotope cross section if applicable
160 const G4double emax = pv->GetMaxEnergy();
161 if(ekin <= emax && amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z]) {
162 auto pviso = data->GetComponentDataByIndex(Z, A - amin[Z]);
163 if(pviso) {
164 xs = pviso->LogVectorValue(ekin, logekin);
165#ifdef G4VERBOSE
166 if(verboseLevel > 1) {
167 G4cout << "G4NeutronInelasticXS::IsoXS: Ekin(MeV)= "
168 << ekin/CLHEP::MeV
169 << " xs(b)= " << xs/CLHEP::barn
170 << " Z= " << Z << " A= " << A << G4endl;
171 }
172#endif
173 return xs;
174 }
175 }
176
177 // use element x-section
178 if(ekin <= emax) {
179 xs = pv->LogVectorValue(ekin, logekin);
180 } else {
182 ekin, Z, aeff[Z]);
183 }
184 xs *= A/aeff[Z];
185#ifdef G4VERBOSE
186 if(verboseLevel > 1) {
187 G4cout << "G4NeutronInelasticXS::IsoXS: Z= " << Z << " A= " << A
188 << " Ekin(MeV)= " << ekin/CLHEP::MeV
189 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
190 }
191#endif
192 return xs;
193}
194
196 const G4Element* anElement, G4double kinEnergy, G4double logE)
197{
198 size_t nIso = anElement->GetNumberOfIsotopes();
199 const G4Isotope* iso = anElement->GetIsotope(0);
200
201 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
202 if(1 == nIso) { return iso; }
203
204 // more than 1 isotope
205 G4int Z = anElement->GetZasInt();
206 //G4cout << "SelectIsotope Z= " << Z << G4endl;
207
208 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
210 G4double sum = 0.0;
211 size_t j;
212
213 // isotope wise cross section not available
214 if(amax[Z] == amin[Z] || Z >= MAXZINEL) {
215 for (j=0; j<nIso; ++j) {
216 sum += abundVector[j];
217 if(q <= sum) {
218 iso = anElement->GetIsotope(j);
219 break;
220 }
221 }
222 return iso;
223 }
224
225 // use isotope cross sections
226 size_t nn = temp.size();
227 if(nn < nIso) { temp.resize(nIso, 0.); }
228
229 for (j=0; j<nIso; ++j) {
230 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
231 // << " abund= " << abundVector[j] << G4endl;
232 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
233 anElement->GetIsotope(j)->GetN());
234 temp[j] = sum;
235 }
236 sum *= q;
237 for (j = 0; j<nIso; ++j) {
238 if(temp[j] >= sum) {
239 iso = anElement->GetIsotope(j);
240 break;
241 }
242 }
243 return iso;
244}
245
246void
248{
249 if(verboseLevel > 0) {
250 G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
251 << p.GetParticleName() << G4endl;
252 }
253 if(p.GetParticleName() != "neutron") {
255 ed << p.GetParticleName() << " is a wrong particle type -"
256 << " only neutron is allowed";
257 G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
258 FatalException, ed, "");
259 return;
260 }
261
262 if(nullptr == data) {
263#ifdef G4MULTITHREADED
264 G4MUTEXLOCK(&neutronInelasticXSMutex);
265 if(nullptr == data) {
266#endif
267 isMaster = true;
268 data = new G4ElementData();
269 data->SetName("NeutronInelastic");
271#ifdef G4MULTITHREADED
272 }
273 G4MUTEXUNLOCK(&neutronInelasticXSMutex);
274#endif
275 }
276
277 // it is possible re-initialisation for the new run
279 if(isMaster) {
280 // Upload data for elements used in geometry
281 for ( auto & elm : *table ) {
282 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINEL-1) );
283 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
284 }
285 }
286 // prepare isotope selection
287 size_t nIso = temp.size();
288 for ( auto & elm : *table ) {
289 size_t n = elm->GetNumberOfIsotopes();
290 if(n > nIso) { nIso = n; }
291 }
292 temp.resize(nIso, 0.0);
293}
294
296{
297 // check environment variable
298 // build the complete string identifying the file with the data set
299 if(gDataDirectory.empty()) {
300 char* path = std::getenv("G4PARTICLEXSDATA");
301 if (nullptr != path) {
302 std::ostringstream ost;
303 ost << path << "/neutron/inel";
304 gDataDirectory = ost.str();
305 } else {
306 G4Exception("G4NeutronInelasticXS::Initialise(..)","had013",
308 "Environment variable G4PARTICLEXSDATA is not defined");
309 }
310 }
311 return gDataDirectory;
312}
313
315{
316#ifdef G4MULTITHREADED
317 G4MUTEXLOCK(&neutronInelasticXSMutex);
318 if(nullptr == data->GetElementData(Z)) {
319#endif
320 Initialise(Z);
321#ifdef G4MULTITHREADED
322 }
323 G4MUTEXUNLOCK(&neutronInelasticXSMutex);
324#endif
325}
326
328{
329 if(nullptr != data->GetElementData(Z)) { return; }
330
331 // upload element data
332 std::ostringstream ost;
333 ost << FindDirectoryPath() << Z;
334 G4PhysicsVector* v = RetrieveVector(ost, true);
336 /*
337 G4cout << "G4NeutronInelasticXS::Initialise for Z= " << Z
338 << " A= " << Amean << " Amin= " << amin[Z]
339 << " Amax= " << amax[Z] << G4endl;
340 */
341 // upload isotope data
342 if(amin[Z] < amax[Z]) {
343 G4int nmax = amax[Z] - amin[Z] + 1;
345
346 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
347 std::ostringstream ost1;
348 ost1 << gDataDirectory << Z << "_" << A;
349 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
350 data->AddComponent(Z, A, v1);
351 }
352 }
353
354 // smooth transition
355 G4double sig1 = (*v)[v->GetVectorLength()-1];
356 G4double ehigh= v->GetMaxEnergy();
358 ehigh, Z, aeff[Z]);
359 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
360}
361
363G4NeutronInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
364{
365 G4PhysicsLogVector* v = nullptr;
366 std::ifstream filein(ost.str().c_str());
367 if (!filein.is_open()) {
368 if(warn) {
370 ed << "Data file <" << ost.str().c_str()
371 << "> is not opened!";
372 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had014",
373 FatalException, ed, "Check G4PARTICLEXSDATA");
374 }
375 } else {
376 if(verboseLevel > 1) {
377 G4cout << "File " << ost.str()
378 << " is opened by G4NeutronInelasticXS" << G4endl;
379 }
380 // retrieve data from DB
381 v = new G4PhysicsLogVector();
382 if(!v->Retrieve(filein, true)) {
384 ed << "Data file <" << ost.str().c_str()
385 << "> is not retrieved!";
386 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had015",
387 FatalException, ed, "Check G4PARTICLEXSDATA");
388 }
389 }
390 return v;
391}
static const G4double emax
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 G4int amax[]
static const G4double aeff[]
static const G4int amin[]
#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
G4PhysicsVector * GetComponentDataByIndex(G4int Z, G4int idx)
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
G4PhysicsVector * GetElementData(G4int Z)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
void SetName(const G4String &nam)
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
G4int GetZasInt() const
Definition: G4Element.hh:132
G4int GetN() const
Definition: G4Isotope.hh:93
static G4ElementData * data
const G4String & FindDirectoryPath()
static const G4int MAXZINEL
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
std::vector< G4double > temp
const G4ParticleDefinition * neutron
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
G4VComponentCrossSection * ggXsection
void CrossSectionDescription(std::ostream &) const final
const G4PhysicsVector * GetPhysicsVector(G4int Z)
G4PhysicsVector * RetrieveVector(std::ostringstream &in, G4bool warn)
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
static G4String gDataDirectory
static G4double coeff[MAXZINEL]
G4double GetElementCrossSection(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
const G4String & GetParticleName() const
G4double GetMaxEnergy() const
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
std::size_t GetVectorLength() const
G4double GetInelasticElementCrossSection(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