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