Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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"
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
52// factory
54//
56
57G4double G4NeutronInelasticXS::coeff[] = {1.0};
58G4ElementData* G4NeutronInelasticXS::data = nullptr;
59G4String G4NeutronInelasticXS::gDataDirectory = "";
60
61#ifdef G4MULTITHREADED
62 G4Mutex G4NeutronInelasticXS::neutronInelasticXSMutex = G4MUTEX_INITIALIZER;
63#endif
64
66 : G4VCrossSectionDataSet(Default_Name()),
67 neutron(G4Neutron::Neutron())
68{
69 verboseLevel = 0;
70 if(verboseLevel > 0){
71 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
72 << MAXZINEL << G4endl;
73 }
75 if(ggXsection == nullptr) ggXsection = new G4ComponentGGHadronNucleusXsc();
77 isMaster = false;
78 temp.resize(13,0.0);
79}
80
82{
83 if(isMaster) { delete data; data = nullptr; }
84}
85
86void G4NeutronInelasticXS::CrossSectionDescription(std::ostream& outFile) const
87{
88 outFile << "G4NeutronInelasticXS calculates the neutron inelastic 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 energy Glauber-Gribov cross section model is used\n";
93}
94
95G4bool
97 G4int, const G4Material*)
98{
99 return true;
100}
101
102G4bool
104 G4int, G4int,
105 const G4Element*, const G4Material*)
106{
107 return true;
108}
109
111 const G4DynamicParticle* aParticle,
112 G4int ZZ, const G4Material*)
113{
114 G4double xs = 0.0;
115 G4double ekin = aParticle->GetKineticEnergy();
116
117 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
118
119 auto pv = GetPhysicsVector(Z);
120 if(pv == nullptr) { return xs; }
121 // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin
122 // << " Z= " << Z << G4endl;
123
124 if(ekin <= pv->GetMaxEnergy()) {
125 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
126 } else {
127 xs = coeff[Z]*ggXsection->GetInelasticElementCrossSection(neutron,
128 ekin, Z, aeff[Z]);
129 }
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
142 const G4DynamicParticle* aParticle,
143 G4int Z, G4int A,
144 const G4Isotope*, const G4Element*,
145 const G4Material*)
146{
147 return IsoCrossSection(aParticle->GetKineticEnergy(),
148 aParticle->GetLogKineticEnergy(), Z, A);
149}
150
153 G4int ZZ, G4int A)
154{
155 G4double xs = 0.0;
156 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
157
158 /*
159 G4cout << "IsoCrossSection Z= " << Z << " A= " << A
160 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
161 << " E(MeV)= " << ekin << G4endl;
162 */
163 auto pv = GetPhysicsVector(Z);
164 if(pv == nullptr) { return xs; }
165
166 // compute isotope cross section if applicable
167 const G4double emax = pv->GetMaxEnergy();
168 if(ekin <= emax && amin[Z] > 0 && A >= amin[Z] && A <= amax[Z]) {
169 auto pviso = data->GetComponentDataByIndex(Z, A - amin[Z]);
170 if(pviso) {
171 xs = pviso->LogVectorValue(ekin, logekin);
172#ifdef G4VERBOSE
173 if(verboseLevel > 1) {
174 G4cout << "G4NeutronInelasticXS::IsoXS: Ekin(MeV)= "
175 << ekin/CLHEP::MeV
176 << " xs(b)= " << xs/CLHEP::barn
177 << " Z= " << Z << " A= " << A << G4endl;
178 }
179#endif
180 return xs;
181 }
182 }
183
184 // use element x-section
185 if(ekin <= emax) {
186 xs = pv->LogVectorValue(ekin, logekin);
187 } else {
188 xs = coeff[Z]*ggXsection->GetInelasticElementCrossSection(neutron,
189 ekin, Z, aeff[Z]);
190 }
191 xs *= A/aeff[Z];
192#ifdef G4VERBOSE
193 if(verboseLevel > 1) {
194 G4cout << "G4NeutronInelasticXS::IsoXS: Z= " << Z << " A= " << A
195 << " Ekin(MeV)= " << ekin/CLHEP::MeV
196 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
197 }
198#endif
199 return xs;
200}
201
203 const G4Element* anElement, G4double kinEnergy, G4double logE)
204{
205 size_t nIso = anElement->GetNumberOfIsotopes();
206 const G4Isotope* iso = anElement->GetIsotope(0);
207
208 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
209 if(1 == nIso) { return iso; }
210
211 // more than 1 isotope
212 G4int Z = anElement->GetZasInt();
213 //G4cout << "SelectIsotope Z= " << Z << G4endl;
214
215 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
217 G4double sum = 0.0;
218 size_t j;
219
220 // isotope wise cross section not available
221 if(0 == amin[Z] || Z >= MAXZINEL) {
222 for (j=0; j<nIso; ++j) {
223 sum += abundVector[j];
224 if(q <= sum) {
225 iso = anElement->GetIsotope(j);
226 break;
227 }
228 }
229 return iso;
230 }
231
232 // use isotope cross sections
233 size_t nn = temp.size();
234 if(nn < nIso) { temp.resize(nIso, 0.); }
235
236 for (j=0; j<nIso; ++j) {
237 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
238 // << " abund= " << abundVector[j] << G4endl;
239 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
240 anElement->GetIsotope(j)->GetN());
241 temp[j] = sum;
242 }
243 sum *= q;
244 for (j = 0; j<nIso; ++j) {
245 if(temp[j] >= sum) {
246 iso = anElement->GetIsotope(j);
247 break;
248 }
249 }
250 return iso;
251}
252
253void
255{
256 if(verboseLevel > 0) {
257 G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
258 << p.GetParticleName() << G4endl;
259 }
260 if(p.GetParticleName() != "neutron") {
262 ed << p.GetParticleName() << " is a wrong particle type -"
263 << " only neutron is allowed";
264 G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
265 FatalException, ed, "");
266 return;
267 }
268
269 if(nullptr == data) {
270#ifdef G4MULTITHREADED
271 G4MUTEXLOCK(&neutronInelasticXSMutex);
272 if(nullptr == data) {
273#endif
274 isMaster = true;
275 data = new G4ElementData();
276 data->SetName("NeutronInelastic");
277 FindDirectoryPath();
278#ifdef G4MULTITHREADED
279 }
280 G4MUTEXUNLOCK(&neutronInelasticXSMutex);
281#endif
282 }
283
284 // it is possible re-initialisation for the new run
285 if(isMaster) {
286
287 // Access to elements
288 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
289 size_t numOfCouples = theCoupleTable->GetTableSize();
290 for(size_t j=0; j<numOfCouples; ++j) {
291 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
292 auto elmVec = mat->GetElementVector();
293 size_t numOfElem = mat->GetNumberOfElements();
294 for (size_t ie = 0; ie < numOfElem; ++ie) {
295 G4int Z = std::max(1,std::min(((*elmVec)[ie])->GetZasInt(), MAXZINEL-1));
296 if(nullptr == data->GetElementData(Z)) { Initialise(Z); }
297 }
298 }
299 }
300}
301
302const G4String& G4NeutronInelasticXS::FindDirectoryPath()
303{
304 // check environment variable
305 // build the complete string identifying the file with the data set
306 if(gDataDirectory.empty()) {
307 char* path = std::getenv("G4PARTICLEXSDATA");
308 if (path) {
309 std::ostringstream ost;
310 ost << path << "/neutron/inel";
311 gDataDirectory = ost.str();
312 } else {
313 G4Exception("G4NeutronInelasticXS::Initialise(..)","had013",
315 "Environment variable G4PARTICLEXSDATA is not defined");
316 }
317 }
318 return gDataDirectory;
319}
320
321void G4NeutronInelasticXS::InitialiseOnFly(G4int Z)
322{
323#ifdef G4MULTITHREADED
324 G4MUTEXLOCK(&neutronInelasticXSMutex);
325 if(nullptr == data->GetElementData(Z)) {
326#endif
327 Initialise(Z);
328#ifdef G4MULTITHREADED
329 }
330 G4MUTEXUNLOCK(&neutronInelasticXSMutex);
331#endif
332}
333
334void G4NeutronInelasticXS::Initialise(G4int Z)
335{
336 if(nullptr != data->GetElementData(Z)) { return; }
337
338 // upload element data
339 std::ostringstream ost;
340 ost << FindDirectoryPath() << Z;
341 G4PhysicsVector* v = RetrieveVector(ost, true);
342 data->InitialiseForElement(Z, v);
343 /*
344 G4cout << "G4NeutronInelasticXS::Initialise for Z= " << Z
345 << " A= " << Amean << " Amin= " << amin[Z]
346 << " Amax= " << amax[Z] << G4endl;
347 */
348 // upload isotope data
349 if(amin[Z] > 0) {
350 size_t nmax = (size_t)(amax[Z]-amin[Z]+1);
351 data->InitialiseForComponent(Z, nmax);
352
353 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
354 std::ostringstream ost1;
355 ost1 << gDataDirectory << Z << "_" << A;
356 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
357 data->AddComponent(Z, A, v1);
358 }
359 }
360
361 // smooth transition
362 G4double sig1 = (*v)[v->GetVectorLength()-1];
363 G4double ehigh= v->GetMaxEnergy();
364 G4double sig2 = ggXsection->GetInelasticElementCrossSection(neutron,
365 ehigh, Z, aeff[Z]);
366 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
367}
368
370G4NeutronInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
371{
372 G4PhysicsLogVector* v = nullptr;
373 std::ifstream filein(ost.str().c_str());
374 if (!(filein)) {
375 if(warn) {
377 ed << "Data file <" << ost.str().c_str()
378 << "> is not opened!";
379 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had014",
380 FatalException, ed, "Check G4PARTICLEXSDATA");
381 }
382 } else {
383 if(verboseLevel > 1) {
384 G4cout << "File " << ost.str()
385 << " is opened by G4NeutronInelasticXS" << G4endl;
386 }
387 // retrieve data from DB
388 v = new G4PhysicsLogVector();
389 if(!v->Retrieve(filein, true)) {
391 ed << "Data file <" << ost.str().c_str()
392 << "> is not retrieved!";
393 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had015",
394 FatalException, ed, "Check G4PARTICLEXSDATA");
395 }
396 }
397 return v;
398}
#define G4_DECLARE_XS_FACTORY(cross_section)
double A(double temperature)
@ 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
const G4int MAXZINEL
#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
#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
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
G4PhysicsVector * GetElementData(G4int Z)
G4PhysicsVector * GetComponentDataByIndex(G4int Z, size_t idx)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
void SetName(const G4String &nam)
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4int GetZasInt() const
Definition: G4Element.hh:131
G4int GetN() const
Definition: G4Isotope.hh:93
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)
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
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii) final
G4double GetMaxEnergy() const
std::size_t GetVectorLength() const
G4double LogVectorValue(const G4double theEnergy, const G4double theLogEnergy) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
void SetForAllAtomsAndEnergies(G4bool val)