Geant4 10.7.0
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"
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"
52
53#include <fstream>
54#include <sstream>
55
56G4ElementData* G4ParticleInelasticXS::data[] = {nullptr};
57G4double G4ParticleInelasticXS::coeff[MAXZINELP][5] = {{1.0}, {1.0}, {1.0}, {1.0}, {1.0}};
58G4String G4ParticleInelasticXS::gDataDirectory[] = {""};
59
60#ifdef G4MULTITHREADED
61 G4Mutex G4ParticleInelasticXS::particleInelasticXSMutex = G4MUTEX_INITIALIZER;
62#endif
63
65 : G4VCrossSectionDataSet("G4ParticleInelasticXS"),
66 highEnergyXsection(nullptr),
67 particle(part),
68 index(0),
69 isMaster(false)
70{
71 if(!part) {
72 G4Exception("G4ParticleInelasticXS::G4ParticleInelasticXS(..)","had015",
73 FatalException, "NO particle definition in constructor");
74 } else {
75 verboseLevel = 0;
76 const G4String& particleName = particle->GetParticleName();
77 if(verboseLevel > 1) {
78 G4cout << "G4ParticleInelasticXS::G4ParticleInelasticXS for "
79 << particleName << " on atoms with Z < " << MAXZINELP << G4endl;
80 }
81 if(particleName == "proton") {
82 highEnergyXsection = G4CrossSectionDataSetRegistry::Instance()->GetComponentCrossSection("Glauber-Gribov");
83 if(highEnergyXsection == nullptr) {
84 highEnergyXsection = new G4ComponentGGHadronNucleusXsc();
85 }
86 } else {
87 highEnergyXsection = G4CrossSectionDataSetRegistry::Instance()->GetComponentCrossSection("Glauber-Gribov Nucl-nucl");
88 if(highEnergyXsection == nullptr) {
89 highEnergyXsection = new G4ComponentGGNuclNuclXsc();
90 }
91 if(particleName == "deuteron") index = 1;
92 else if(particleName == "triton") index = 2;
93 else if(particleName == "He3") index = 3;
94 else if(particleName == "alpha") index = 4;
95 else {
97 ed << particleName << " is a wrong particle type";
98 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
99 FatalException, ed, "");
100 }
101 }
102 }
104 temp.resize(13,0.0);
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
140 const G4DynamicParticle* aParticle,
141 G4int ZZ, const G4Material*)
142{
143 G4double xs = 0.0;
144 G4double ekin = aParticle->GetKineticEnergy();
145
146 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
147
148 auto pv = GetPhysicsVector(Z);
149 if(nullptr == pv) { return xs; }
150 // G4cout << "G4ParticleInelasticXS::GetCrossSection e= " << ekin
151 // << " Z= " << Z << G4endl;
152
153 if(ekin <= pv->GetMaxEnergy()) {
154 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
155 } else {
156 xs = coeff[Z][index]*highEnergyXsection->GetInelasticElementCrossSection(particle,
157 ekin, Z, aeff[Z]);
158 }
159
160#ifdef G4VERBOSE
161 if(verboseLevel > 1) {
162 G4cout << "ElmXS: Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
163 << " xs(bn)= " << xs/CLHEP::barn << " element data for "
164 << particle->GetParticleName()
165 << " idx= " << index << G4endl;
166 }
167#endif
168 return xs;
169}
170
172 const G4DynamicParticle* aParticle,
173 G4int Z, G4int A,
174 const G4Isotope*, const G4Element*, const G4Material*)
175{
176 return IsoCrossSection(aParticle->GetKineticEnergy(),
177 aParticle->GetLogKineticEnergy(),Z, A);
178}
179
182 G4int ZZ, G4int A)
183{
184 G4double xs = 0.0;
185 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
186 /*
187 G4cout << "G4ParticleInelasticXS: IsoCrossSection Z= "
188 << Z << " A= " << A
189 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
190 << " E(MeV)= " << ekin << G4endl;
191 */
192 auto pv = GetPhysicsVector(Z);
193 if(pv == nullptr) { return xs; }
194
195 // compute isotope cross section if applicable
196 const G4double emax = pv->GetMaxEnergy();
197 if(ekin <= emax && amin[Z]>0 && A >= amin[Z] && A <= amax[Z]) {
198 auto pviso = data[index]->GetComponentDataByIndex(Z, A - amin[Z]);
199 if(pviso != nullptr) {
200 xs = pviso->LogVectorValue(ekin, logE);
201#ifdef G4VERBOSE
202 if(verboseLevel > 1) {
203 G4cout << "G4ParticleInelasticXS::IsoXS: for "
204 << particle->GetParticleName() << " Ekin(MeV)= "
205 << ekin/CLHEP::MeV << " xs(b)= " << xs/CLHEP::barn
206 << " Z= " << Z << " A= " << A
207 << " idx= " << index << G4endl;
208 }
209#endif
210 return xs;
211 }
212 }
213 // use element x-section
214 if(ekin <= emax) {
215 xs = pv->LogVectorValue(ekin, logE);
216 } else {
217 xs = coeff[Z][index] *
218 highEnergyXsection->GetInelasticElementCrossSection(particle,
219 ekin, Z, aeff[Z]);
220 }
221 xs *= A/aeff[Z];
222#ifdef G4VERBOSE
223 if(verboseLevel > 1) {
224 G4cout << "IsoXS for " << particle->GetParticleName()
225 << " Target Z= " << Z << " A= " << A
226 << " Ekin(MeV)= " << ekin/CLHEP::MeV
227 << " xs(bn)= " << xs/CLHEP::barn
228 << " idx= " << index << G4endl;
229 }
230#endif
231 return xs;
232}
233
235 const G4Element* anElement, G4double kinEnergy, G4double logE)
236{
237 size_t nIso = anElement->GetNumberOfIsotopes();
238 const G4Isotope* iso = anElement->GetIsotope(0);
239
240 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
241 if(1 == nIso) { return iso; }
242
243 // more than 1 isotope
244 G4int Z = anElement->GetZasInt();
245 //G4cout << "SelectIsotope Z= " << Z << G4endl;
246
247 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
249 G4double sum = 0.0;
250 size_t j;
251
252 // isotope wise cross section not available
253 if(0 == 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 size_t nn = temp.size();
265 if(nn < nIso) { temp.resize(nIso, 0.); }
266
267 for (j=0; j<nIso; ++j) {
268 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
269 // << " abund= " << abundVector[j] << G4endl;
270 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
271 anElement->GetIsotope(j)->GetN());
272 temp[j] = sum;
273 }
274 sum *= q;
275 for (j=0; j<nIso; ++j) {
276 if(temp[j] >= sum) {
277 iso = anElement->GetIsotope(j);
278 break;
279 }
280 }
281 return iso;
282}
283
284void
286{
287 if(verboseLevel > 0){
288 G4cout << "G4ParticleInelasticXS::BuildPhysicsTable for "
289 << p.GetParticleName() << G4endl;
290 }
291 if(&p != particle) {
293 ed << p.GetParticleName() << " is a wrong particle type -"
294 << particle->GetParticleName() << " is expected";
295 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
296 FatalException, ed, "");
297 return;
298 }
299
300 if(data[index] == nullptr) {
301#ifdef G4MULTITHREADED
302 G4MUTEXLOCK(&particleInelasticXSMutex);
303 if(data[index] == nullptr) {
304#endif
305 isMaster = true;
306 data[index] = new G4ElementData();
307 data[index]->SetName(particle->GetParticleName() + "Inelastic");
308 FindDirectoryPath();
309#ifdef G4MULTITHREADED
310 }
311 G4MUTEXUNLOCK(&particleInelasticXSMutex);
312#endif
313 }
314
315 // it is possible re-initialisation for the new run
316 if(isMaster) {
317
318 // Access to elements
319 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
320 size_t numOfCouples = theCoupleTable->GetTableSize();
321 for(size_t j=0; j<numOfCouples; ++j) {
322 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
323 auto elmVec = mat->GetElementVector();
324 size_t numOfElem = mat->GetNumberOfElements();
325 for (size_t ie = 0; ie < numOfElem; ++ie) {
326 G4int Z = std::max(1,std::min(((*elmVec)[ie])->GetZasInt(), MAXZINELP-1));
327 if(nullptr == data[index]->GetElementData(Z)) { Initialise(Z); }
328 }
329 }
330 }
331}
332
333const G4String& G4ParticleInelasticXS::FindDirectoryPath()
334{
335 // check environment variable
336 // build the complete string identifying the file with the data set
337 if(gDataDirectory[index].empty()) {
338 char* path = std::getenv("G4PARTICLEXSDATA");
339 if (path) {
340 std::ostringstream ost;
341 ost << path << "/" << particle->GetParticleName() << "/inel";
342 gDataDirectory[index] = ost.str();
343 } else {
344 G4Exception("G4NeutronInelasticXS::Initialise(..)","had013",
346 "Environment variable G4PARTICLEXSDATA is not defined");
347 }
348 }
349 return gDataDirectory[index];
350}
351
352void G4ParticleInelasticXS::InitialiseOnFly(G4int Z)
353{
354#ifdef G4MULTITHREADED
355 G4MUTEXLOCK(&particleInelasticXSMutex);
356 if(nullptr == data[index]->GetElementData(Z)) {
357#endif
358 Initialise(Z);
359#ifdef G4MULTITHREADED
360 }
361 G4MUTEXUNLOCK(&particleInelasticXSMutex);
362#endif
363}
364
365void G4ParticleInelasticXS::Initialise(G4int Z)
366{
367 if(nullptr != data[index]->GetElementData(Z)) { return; }
368
369 // upload element data
370 std::ostringstream ost;
371 ost << FindDirectoryPath() << Z ;
372 G4PhysicsVector* v = RetrieveVector(ost, true);
373 data[index]->InitialiseForElement(Z, v);
374 /*
375 G4cout << "G4ParticleInelasticXS::Initialise for Z= " << Z
376 << " idx= " << index
377 << " Amin= " << amin[Z]
378 << " Amax= " << amax[Z]
379 << " " << FindDirectoryPath() << G4endl;
380 */
381 // upload isotope data
382 if(amin[Z] > 0) {
383 size_t nmax = (size_t)(amax[Z]-amin[Z]+1);
384 data[index]->InitialiseForComponent(Z, nmax);
385
386 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
387 std::ostringstream ost1;
388 ost1 << FindDirectoryPath() << Z << "_" << A;
389 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
390 data[index]->AddComponent(Z, A, v1);
391 /*
392 G4cout << " Isotope x-section Z= " << Z << " A= " << A
393 << " v1= " << v1 << G4endl;
394 */
395 }
396 }
397 // smooth transition
398 G4double sig1 = (*v)[v->GetVectorLength()-1];
399 G4double ehigh = v->GetMaxEnergy();
400 G4double sig2 = highEnergyXsection->GetInelasticElementCrossSection(
401 particle, ehigh, Z, aeff[Z]);
402 coeff[Z][index] = (sig2 > 0.) ? sig1/sig2 : 1.0;
403 /*
404 G4cout << "G4ParticleInelasticXS: index= " << index
405 << " Z= " << Z << " coeff= " << coeff[Z][index] << G4endl;
406 */
407}
408
410G4ParticleInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
411{
412 G4PhysicsLogVector* v = nullptr;
413 std::ifstream filein(ost.str().c_str());
414 if (!(filein)) {
415 if(warn) {
417 ed << "Data file <" << ost.str().c_str()
418 << "> is not opened!";
419 G4Exception("G4ParticleInelasticXS::RetrieveVector(..)","had014",
420 FatalException, ed, "Check G4PARTICLEXSDATA");
421 }
422 } else {
423 if(verboseLevel > 1) {
424 G4cout << "File " << ost.str()
425 << " is opened by G4ParticleInelasticXS" << G4endl;
426 }
427 // retrieve data from DB
428 v = new G4PhysicsLogVector();
429 if(!v->Retrieve(filein, true)) {
431 ed << "Data file <" << ost.str().c_str()
432 << "> is not retrieved!";
433 G4Exception("G4ParticleInelasticXS::RetrieveVector(..)","had015",
434 FatalException, ed, "Check G4PARTICLEXSDATA");
435 }
436 }
437 return v;
438}
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 MAXZINELP
#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 * 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
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
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)
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
void CrossSectionDescription(std::ostream &) const final
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)