Geant4-11
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
54#include <fstream>
55#include <sstream>
56
58G4double G4ParticleInelasticXS::coeff[MAXZINELP][5] = {{1.0}, {1.0}, {1.0}, {1.0}, {1.0}};
60
61#ifdef G4MULTITHREADED
62 G4Mutex G4ParticleInelasticXS::particleInelasticXSMutex = G4MUTEX_INITIALIZER;
63#endif
64
66 : G4VCrossSectionDataSet("G4ParticleInelasticXS"),
67 particle(part)
68{
69 if(nullptr == part) {
70 G4Exception("G4ParticleInelasticXS::G4ParticleInelasticXS(..)","had015",
71 FatalException, "NO particle definition in constructor");
72 } else {
73 verboseLevel = 0;
74 const G4String& particleName = particle->GetParticleName();
75 if(verboseLevel > 1) {
76 G4cout << "G4ParticleInelasticXS::G4ParticleInelasticXS for "
77 << particleName << " on atoms with Z < " << MAXZINELP << G4endl;
78 }
79 if(particleName == "proton") {
81 if(highEnergyXsection == nullptr) {
83 }
84 } else {
86 if(highEnergyXsection == nullptr) {
88 }
89 if(particleName == "deuteron") index = 1;
90 else if(particleName == "triton") index = 2;
91 else if(particleName == "He3") index = 3;
92 else if(particleName == "alpha") index = 4;
93 else {
95 ed << particleName << " is a wrong particle type";
96 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
97 FatalException, ed, "");
98 }
99 }
100 }
102}
103
105{
106 if(isMaster) {
107 delete data[index];
108 data[index] = nullptr;
109 }
110}
111
112void G4ParticleInelasticXS::CrossSectionDescription(std::ostream& outFile) const
113{
114 outFile << "G4ParticleInelasticXS calculates n, p, d, t, he3, he4 inelastic\n"
115 << "cross section on nuclei using data from the high precision\n"
116 << "neutron database. These data are simplified and smoothed over\n"
117 << "the resonance region in order to reduce CPU time.\n"
118 << "For high energy Glauber-Gribov cross section model is used\n";
119}
120
121G4bool
123 G4int, const G4Material*)
124{
125 return true;
126}
127
128G4bool
130 G4int, G4int,
131 const G4Element*, const G4Material*)
132{
133 return true;
134}
135
137 const G4DynamicParticle* aParticle,
138 G4int ZZ, const G4Material*)
139{
140 G4double xs = 0.0;
141 G4double ekin = aParticle->GetKineticEnergy();
142
143 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
144
145 auto pv = GetPhysicsVector(Z);
146 if(nullptr == pv) { return xs; }
147 // G4cout << "G4ParticleInelasticXS::GetCrossSection e= " << ekin
148 // << " Z= " << Z << G4endl;
149
150 if(ekin <= pv->GetMaxEnergy()) {
151 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
152 } else {
154 ekin, Z, aeff[Z]);
155 }
156
157#ifdef G4VERBOSE
158 if(verboseLevel > 1) {
159 G4cout << "ElmXS: Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
160 << " xs(bn)= " << xs/CLHEP::barn << " element data for "
162 << " idx= " << index << G4endl;
163 }
164#endif
165 return xs;
166}
167
169 const G4DynamicParticle* aParticle,
170 G4int Z, G4int A,
171 const G4Isotope*, const G4Element*, const G4Material*)
172{
173 return IsoCrossSection(aParticle->GetKineticEnergy(),
174 aParticle->GetLogKineticEnergy(),Z, A);
175}
176
179 G4int ZZ, G4int A)
180{
181 G4double xs = 0.0;
182 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
183 /*
184 G4cout << "G4ParticleInelasticXS: IsoCrossSection Z= "
185 << Z << " A= " << A
186 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
187 << " E(MeV)= " << ekin << G4endl;
188 */
189 auto pv = GetPhysicsVector(Z);
190 if(pv == nullptr) { return xs; }
191
192 // compute isotope cross section if applicable
193 const G4double emax = pv->GetMaxEnergy();
194 if(ekin <= emax && amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z]) {
195 auto pviso = data[index]->GetComponentDataByIndex(Z, A - amin[Z]);
196 if(pviso != nullptr) {
197 xs = pviso->LogVectorValue(ekin, logE);
198#ifdef G4VERBOSE
199 if(verboseLevel > 1) {
200 G4cout << "G4ParticleInelasticXS::IsoXS: for "
201 << particle->GetParticleName() << " Ekin(MeV)= "
202 << ekin/CLHEP::MeV << " xs(b)= " << xs/CLHEP::barn
203 << " Z= " << Z << " A= " << A
204 << " idx= " << index << G4endl;
205 }
206#endif
207 return xs;
208 }
209 }
210 // use element x-section
211 if(ekin <= emax) {
212 xs = pv->LogVectorValue(ekin, logE);
213 } else {
214 xs = coeff[Z][index] *
216 ekin, Z, aeff[Z]);
217 }
218 xs *= A/aeff[Z];
219#ifdef G4VERBOSE
220 if(verboseLevel > 1) {
221 G4cout << "IsoXS for " << particle->GetParticleName()
222 << " Target Z= " << Z << " A= " << A
223 << " Ekin(MeV)= " << ekin/CLHEP::MeV
224 << " xs(bn)= " << xs/CLHEP::barn
225 << " idx= " << index << G4endl;
226 }
227#endif
228 return xs;
229}
230
232 const G4Element* anElement, G4double kinEnergy, G4double logE)
233{
234 size_t nIso = anElement->GetNumberOfIsotopes();
235 const G4Isotope* iso = anElement->GetIsotope(0);
236
237 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
238 if(1 == nIso) { return iso; }
239
240 // more than 1 isotope
241 G4int Z = anElement->GetZasInt();
242 //G4cout << "SelectIsotope Z= " << Z << G4endl;
243
244 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
246 G4double sum = 0.0;
247 size_t j;
248
249 // isotope wise cross section not available
250 if(amax[Z] == amin[Z] || Z >= MAXZINELP) {
251 for (j=0; j<nIso; ++j) {
252 sum += abundVector[j];
253 if(q <= sum) {
254 iso = anElement->GetIsotope(j);
255 break;
256 }
257 }
258 return iso;
259 }
260
261 size_t nn = temp.size();
262 if(nn < nIso) { temp.resize(nIso, 0.); }
263
264 for (j=0; j<nIso; ++j) {
265 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
266 // << " abund= " << abundVector[j] << G4endl;
267 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
268 anElement->GetIsotope(j)->GetN());
269 temp[j] = sum;
270 }
271 sum *= q;
272 for (j=0; j<nIso; ++j) {
273 if(temp[j] >= sum) {
274 iso = anElement->GetIsotope(j);
275 break;
276 }
277 }
278 return iso;
279}
280
281void
283{
284 if(verboseLevel > 0){
285 G4cout << "G4ParticleInelasticXS::BuildPhysicsTable for "
286 << p.GetParticleName() << G4endl;
287 }
288 if(&p != particle) {
290 ed << p.GetParticleName() << " is a wrong particle type -"
291 << particle->GetParticleName() << " is expected";
292 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
293 FatalException, ed, "");
294 return;
295 }
296
297 G4int fact = (p.GetParticleName() == "proton") ? 1 : 256;
298 SetMaxKinEnergy(G4HadronicParameters::Instance()->GetMaxEnergy() * fact);
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");
309#ifdef G4MULTITHREADED
310 }
311 G4MUTEXUNLOCK(&particleInelasticXSMutex);
312#endif
313 }
314
315 // it is possible re-initialisation for the new run
317 if(isMaster) {
318
319 // Access to elements
320 for ( auto & elm : *table ) {
321 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINELP-1) );
322 if ( nullptr == (data[index])->GetElementData(Z) ) { Initialise(Z); }
323 }
324 }
325 // prepare isotope selection
326 size_t nIso = temp.size();
327 for ( auto & elm : *table ) {
328 size_t n = elm->GetNumberOfIsotopes();
329 if(n > nIso) { nIso = n; }
330 }
331 temp.resize(nIso, 0.0);
332}
333
335{
336 // check environment variable
337 // build the complete string identifying the file with the data set
338 if(gDataDirectory[index].empty()) {
339 char* path = std::getenv("G4PARTICLEXSDATA");
340 if (nullptr != path) {
341 std::ostringstream ost;
342 ost << path << "/" << particle->GetParticleName() << "/inel";
343 gDataDirectory[index] = ost.str();
344 } else {
345 G4Exception("G4NeutronInelasticXS::Initialise(..)","had013",
347 "Environment variable G4PARTICLEXSDATA is not defined");
348 }
349 }
350 return gDataDirectory[index];
351}
352
354{
355#ifdef G4MULTITHREADED
356 G4MUTEXLOCK(&particleInelasticXSMutex);
357 if(nullptr == data[index]->GetElementData(Z)) {
358#endif
359 Initialise(Z);
360#ifdef G4MULTITHREADED
361 }
362 G4MUTEXUNLOCK(&particleInelasticXSMutex);
363#endif
364}
365
367{
368 if(nullptr != data[index]->GetElementData(Z)) { return; }
369
370 // upload element data
371 std::ostringstream ost;
372 ost << FindDirectoryPath() << Z ;
373 G4PhysicsVector* v = RetrieveVector(ost, true);
375 /*
376 G4cout << "G4ParticleInelasticXS::Initialise for Z= " << Z
377 << " idx= " << index
378 << " Amin= " << amin[Z]
379 << " Amax= " << amax[Z]
380 << " " << FindDirectoryPath() << G4endl;
381 */
382 // upload isotope data
383 if(amin[Z] < amax[Z]) {
384 G4int nmax = amax[Z] - amin[Z] + 1;
386
387 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
388 std::ostringstream ost1;
389 ost1 << FindDirectoryPath() << Z << "_" << A;
390 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
391 data[index]->AddComponent(Z, A, v1);
392 /*
393 G4cout << " Isotope x-section Z= " << Z << " A= " << A
394 << " v1= " << v1 << G4endl;
395 */
396 }
397 }
398 // smooth transition
399 G4double sig1 = (*v)[v->GetVectorLength()-1];
400 G4double ehigh = v->GetMaxEnergy();
402 particle, ehigh, Z, aeff[Z]);
403 coeff[Z][index] = (sig2 > 0.) ? sig1/sig2 : 1.0;
404 /*
405 G4cout << "G4ParticleInelasticXS: index= " << index
406 << " Z= " << Z << " coeff= " << coeff[Z][index] << G4endl;
407 */
408}
409
411G4ParticleInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
412{
413 G4PhysicsLogVector* v = nullptr;
414 std::ifstream filein(ost.str().c_str());
415 if (!filein.is_open()) {
416 if(warn) {
418 ed << "Data file <" << ost.str().c_str()
419 << "> is not opened!";
420 G4Exception("G4ParticleInelasticXS::RetrieveVector(..)","had014",
421 FatalException, ed, "Check G4PARTICLEXSDATA");
422 }
423 } else {
424 if(verboseLevel > 1) {
425 G4cout << "File " << ost.str()
426 << " is opened by G4ParticleInelasticXS" << G4endl;
427 }
428 // retrieve data from DB
429 v = new G4PhysicsLogVector();
430 if(!v->Retrieve(filein, true)) {
432 ed << "Data file <" << ost.str().c_str()
433 << "> is not retrieved!";
434 G4Exception("G4ParticleInelasticXS::RetrieveVector(..)","had015",
435 FatalException, ed, "Check G4PARTICLEXSDATA");
436 }
437 }
438 return v;
439}
static const G4double emax
std::vector< G4Element * > G4ElementTable
@ 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
static const G4int amax[]
static const G4double aeff[]
static const G4int amin[]
#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
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:397
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
const G4String & FindDirectoryPath()
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
G4VComponentCrossSection * highEnergyXsection
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
static const G4int MAXZINELP
G4ParticleInelasticXS(const G4ParticleDefinition *)
const G4PhysicsVector * GetPhysicsVector(G4int Z)
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
static G4ElementData * data[5]
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
static G4double coeff[MAXZINELP][5]
G4double IsoCrossSection(G4double ekin, G4double logE, G4int Z, G4int A)
std::vector< G4double > temp
G4PhysicsVector * RetrieveVector(std::ostringstream &in, G4bool warn)
const G4ParticleDefinition * particle
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) final
static G4String gDataDirectory[5]
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)
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double MeV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments