Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleInelasticXS.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: G4ParticleInelasticXS
32//
33// Author Ivantchenko, Geant4, 3-Aug-09
34//
35// Modifications:
36//
37
39#include "G4Neutron.hh"
40#include "G4DynamicParticle.hh"
41#include "G4ElementTable.hh"
42#include "G4Material.hh"
43#include "G4Element.hh"
44#include "G4PhysicsLogVector.hh"
49#include "G4Proton.hh"
50#include "Randomize.hh"
51#include "G4SystemOfUnits.hh"
52#include "G4IsotopeList.hh"
54#include "G4AutoLock.hh"
55
56#include <fstream>
57#include <sstream>
58#include <thread>
59
60G4ElementData* G4ParticleInelasticXS::data[] = {nullptr, nullptr, nullptr, nullptr, nullptr};
61G4double G4ParticleInelasticXS::coeff[MAXZINELP][5] = {{1.0}, {1.0}, {1.0}, {1.0}, {1.0}};
62G4String G4ParticleInelasticXS::gDataDirectory[] = {"", "", "", "", ""};
63
64static std::once_flag applyOnce;
65
66namespace
67{
68 G4Mutex pInelasticXSMutex = G4MUTEX_INITIALIZER;
69 G4String pname[5] = {"proton", "deuteron", "triton", "He3", "alpha"};
70}
71
73 : G4VCrossSectionDataSet("G4ParticleInelasticXS"),
74 particle(part),
75 elimit(20*CLHEP::MeV)
76{
77 if(nullptr == part) {
78 G4Exception("G4ParticleInelasticXS::G4ParticleInelasticXS(..)","had015",
79 FatalException, "NO particle definition in constructor");
80 } else {
81 verboseLevel = 0;
82 const G4String& particleName = particle->GetParticleName();
83 if(verboseLevel > 1) {
84 G4cout << "G4ParticleInelasticXS::G4ParticleInelasticXS for "
85 << particleName << " on atoms with Z < " << MAXZINELP << G4endl;
86 }
88 if(particleName == "proton") {
89 highEnergyXsection = xsr->GetComponentCrossSection("Glauber-Gribov");
90 if(highEnergyXsection == nullptr) {
91 highEnergyXsection = new G4ComponentGGHadronNucleusXsc();
92 }
93 } else {
94 highEnergyXsection =
95 xsr->GetComponentCrossSection("Glauber-Gribov Nucl-nucl");
96 if(highEnergyXsection == nullptr) {
97 highEnergyXsection = new G4ComponentGGNuclNuclXsc();
98 }
99 for (index=1; index<5; ++index) {
100 if (particleName == pname[index]) { break; }
101 }
102 if (index == 5) {
104 ed << particleName << " is a wrong particle type";
105 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
106 FatalException, ed, "");
107 }
108 if (1 < index) { SetMaxKinEnergy(25.6*CLHEP::PeV); }
109 }
110 }
112 if (data[0] == nullptr) {
113 for (G4int i=0; i<5; ++i) {
114 data[i] = new G4ElementData(MAXZINELP);
115 data[i]->SetName(pname[i] + "IonInel");
116 }
117 FindDirectoryPath();
118 }
119}
120
121void G4ParticleInelasticXS::CrossSectionDescription(std::ostream& outFile) const
122{
123 outFile << "G4ParticleInelasticXS calculates n, p, d, t, he3, he4 inelastic\n"
124 << "cross section on nuclei using data from the high precision\n"
125 << "neutron database. These data are simplified and smoothed over\n"
126 << "the resonance region in order to reduce CPU time.\n"
127 << "For high energy Glauber-Gribov cross section model is used.\n";
128}
129
130G4bool
132 G4int, const G4Material*)
133{
134 return true;
135}
136
137G4bool
139 G4int, G4int,
140 const G4Element*, const G4Material*)
141{
142 return true;
143}
144
147 G4int Z, const G4Material*)
148{
149 return ElementCrossSection(aParticle->GetKineticEnergy(),
150 aParticle->GetLogKineticEnergy(), Z);
151}
152
156 const G4Element* elm,
157 const G4Material*)
158{
159 return ElementCrossSection(ekin, loge, elm->GetZasInt());
160}
161
163{
164 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
165 auto pv = GetPhysicsVector(Z);
166
167 G4double xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
168 : coeff[Z][index]*highEnergyXsection->GetInelasticElementCrossSection(particle,
169 ekin, Z, aeff[Z]);
170
171#ifdef G4VERBOSE
172 if(verboseLevel > 1) {
173 G4cout << "ElmXS: Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
174 << " xs(bn)= " << xs/CLHEP::barn << " element data for "
175 << particle->GetParticleName()
176 << " idx= " << index << G4endl;
177 }
178#endif
179 return xs;
180}
181
185 G4int Z, G4int A, const G4Isotope*,
186 const G4Element*, const G4Material*)
187{
188 return IsoCrossSection(ekin, loge, Z, A);
189}
190
193 G4int Z, G4int A, const G4Isotope*,
194 const G4Element*, const G4Material*)
195{
196 return IsoCrossSection(aParticle->GetKineticEnergy(),
197 aParticle->GetLogKineticEnergy(), Z, A);
198}
199
202 G4int ZZ, G4int A)
203{
204 G4double xs = 0.0;
205 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
206 auto pv = GetPhysicsVector(Z);
207
208 // compute isotope cross section if applicable
209 if (ekin <= elimit && data[index]->GetNumberOfComponents(Z) > 0) {
210 auto pviso = data[index]->GetComponentDataByID(Z, A);
211 if(pviso != nullptr) {
212 xs = pviso->LogVectorValue(ekin, logE);
213#ifdef G4VERBOSE
214 if(verboseLevel > 1) {
215 G4cout << "G4ParticleInelasticXS::IsoXS: for "
216 << particle->GetParticleName() << " Ekin(MeV)= "
217 << ekin/CLHEP::MeV << " xs(b)= " << xs/CLHEP::barn
218 << " Z= " << Z << " A= " << A
219 << " idx= " << index << G4endl;
220 }
221#endif
222 return xs;
223 }
224 }
225 // use element x-section
226 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, logE)
227 : coeff[Z][index] *
228 highEnergyXsection->GetInelasticElementCrossSection(particle,
229 ekin, Z, aeff[Z]);
230 xs *= A/aeff[Z];
231#ifdef G4VERBOSE
232 if(verboseLevel > 1) {
233 G4cout << "IsoXS for " << particle->GetParticleName()
234 << " Target Z= " << Z << " A= " << A
235 << " Ekin(MeV)= " << ekin/CLHEP::MeV
236 << " xs(bn)= " << xs/CLHEP::barn
237 << " idx= " << index << G4endl;
238 }
239#endif
240 return xs;
241}
242
244 const G4Element* anElement, G4double kinEnergy, G4double logE)
245{
246 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
247 const G4Isotope* iso = anElement->GetIsotope(0);
248
249 if (1 == nIso) { return iso; }
250
251 // more than 1 isotope
252 G4int Z = anElement->GetZasInt();
253 if (nullptr == data[index]->GetElementData(Z)) { InitialiseOnFly(Z); }
254
255 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
257 G4double sum = 0.0;
258 G4int j;
259
260 // isotope wise cross section not available
261 if (Z >= MAXZINELP || 0 == data[index]->GetNumberOfComponents(Z)) {
262 for (j=0; j<nIso; ++j) {
263 sum += abundVector[j];
264 if(q <= sum) {
265 iso = anElement->GetIsotope(j);
266 break;
267 }
268 }
269 return iso;
270 }
271
272 G4int nn = (G4int)temp.size();
273 if (nn < nIso) { temp.resize(nIso, 0.); }
274
275 for (j=0; j<nIso; ++j) {
276 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
277 anElement->GetIsotope(j)->GetN());
278 temp[j] = sum;
279 }
280 sum *= q;
281 for (j=0; j<nIso; ++j) {
282 if (temp[j] >= sum) {
283 iso = anElement->GetIsotope(j);
284 break;
285 }
286 }
287 return iso;
288}
289
290void
292{
293 if (verboseLevel > 0){
294 G4cout << "G4ParticleInelasticXS::BuildPhysicsTable for "
295 << p.GetParticleName() << G4endl;
296 }
297 if (&p != particle) {
299 ed << p.GetParticleName() << " is a wrong particle type -"
300 << particle->GetParticleName() << " is expected";
301 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
302 FatalException, ed, "");
303 return;
304 }
305
306 // it is possible re-initialisation for the new run
308
309 // initialise static tables only once
310 std::call_once(applyOnce, [this]() { isInitializer = true; });
311
312 if (isInitializer) {
313 G4AutoLock l(&pInelasticXSMutex);
314
315 // Access to elements
316 for ( auto const & elm : *table ) {
317 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINELP-1) );
318 for (G4int i=0; i<5; ++i) {
319 if ( nullptr == (data[i])->GetElementData(Z) ) { Initialise(Z, i); }
320 }
321 }
322 l.unlock();
323 }
324 // prepare isotope selection
325 std::size_t nIso = temp.size();
326 for ( auto const & elm : *table ) {
327 std::size_t n = elm->GetNumberOfIsotopes();
328 if (n > nIso) { nIso = n; }
329 }
330 temp.resize(nIso, 0.0);
331}
332
333void G4ParticleInelasticXS::FindDirectoryPath()
334{
335 // build the complete string identifying the file with the data set
336 if (gDataDirectory[0].empty()) {
337 for (G4int i=0; i<5; ++i) {
338 std::ostringstream ost;
340 << pname[i] << "/inel";
341 gDataDirectory[i] = ost.str();
342 }
343 }
344}
345
346void G4ParticleInelasticXS::InitialiseOnFly(G4int Z)
347{
348 G4AutoLock l(&pInelasticXSMutex);
349 for (G4int i=0; i<5; ++i) {
350 if ( nullptr == (data[i])->GetElementData(Z) ) { Initialise(Z, i); }
351 }
352 l.unlock();
353}
354
355void G4ParticleInelasticXS::Initialise(G4int Z, G4int idx)
356{
357 if (nullptr != data[idx]->GetElementData(Z)) { return; }
358
359 // upload element data
360 std::ostringstream ost;
361 ost << gDataDirectory[idx] << Z ;
362 G4PhysicsVector* v = RetrieveVector(ost, true);
363 data[idx]->InitialiseForElement(Z, v);
364
365 // upload isotope data
366 G4bool noComp = true;
367 if (amin[Z] < amax[Z]) {
368
369 for (G4int A=amin[Z]; A<=amax[Z]; ++A) {
370 std::ostringstream ost1;
371 ost1 << gDataDirectory[idx] << Z << "_" << A;
372 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
373 if (nullptr != v1) {
374 if (noComp) {
375 G4int nmax = amax[Z] - A + 1;
376 data[idx]->InitialiseForComponent(Z, nmax);
377 noComp = false;
378 }
379 data[idx]->AddComponent(Z, A, v1);
380 }
381 }
382 }
383 // no components case
384 if (noComp) { data[idx]->InitialiseForComponent(Z, 0); }
385
386 // smooth transition
387 G4double sig1 = (*v)[v->GetVectorLength()-1];
388 G4double ehigh = v->GetMaxEnergy();
389 G4double sig2 = highEnergyXsection->GetInelasticElementCrossSection(
390 particle, ehigh, Z, aeff[Z]);
391 coeff[Z][idx] = (sig2 > 0.) ? sig1/sig2 : 1.0;
392}
393
395G4ParticleInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
396{
397 G4PhysicsLogVector* v = nullptr;
398 std::ifstream filein(ost.str().c_str());
399 if (!filein.is_open()) {
400 if(warn) {
402 ed << "Data file <" << ost.str().c_str()
403 << "> is not opened!";
404 G4Exception("G4ParticleInelasticXS::RetrieveVector(..)","had014",
405 FatalException, ed, "Check G4PARTICLEXSDATA");
406 }
407 } else {
408 if(verboseLevel > 1) {
409 G4cout << "File " << ost.str()
410 << " is opened by G4ParticleInelasticXS" << G4endl;
411 }
412 // retrieve data from DB
413 v = new G4PhysicsLogVector();
414 if(!v->Retrieve(filein, true)) {
416 ed << "Data file <" << ost.str().c_str()
417 << "> is not retrieved!";
418 G4Exception("G4ParticleInelasticXS::RetrieveVector(..)","had015",
419 FatalException, ed, "Check G4PARTICLEXSDATA");
420 }
421 }
422 return v;
423}
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
static G4CrossSectionDataSetRegistry * Instance()
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
G4PhysicsVector * GetComponentDataByID(G4int Z, G4int id) const
void SetName(const G4String &nam)
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
G4int GetZasInt() const
Definition G4Element.hh:120
static G4HadronicParameters * Instance()
const G4String & GetDirPARTICLEXS() const
G4int GetN() const
Definition G4Isotope.hh:83
const G4String & GetParticleName() const
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4double ComputeCrossSectionPerElement(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
G4ParticleInelasticXS(const G4ParticleDefinition *)
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4double IsoCrossSection(G4double ekin, G4double logE, G4int Z, G4int A)
G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)
void CrossSectionDescription(std::ostream &) const final
G4double GetMaxEnergy() const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
std::size_t GetVectorLength() const
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
void SetMaxKinEnergy(G4double value)
void SetForAllAtomsAndEnergies(G4bool val)