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