Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionHP.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// V. Ivanchenko September 2023
27//
28// G4CrossSectionHP is a generic class implementing
29// cross sections for neutrons, protons and light ions
30// It is an alternative to code developed by J.P. Wellisch & T.Koi
31//
32
33#include <fstream>
34#include <sstream>
35#include <thread>
36
37#include "G4CrossSectionHP.hh"
38#include "G4Material.hh"
39#include "G4Element.hh"
41#include "G4ElementTable.hh"
42#include "G4IsotopeList.hh"
47#include "G4Pow.hh"
48#include "G4Log.hh"
49#include "Randomize.hh"
50#include "G4RandomDirection.hh"
51#include "G4SystemOfUnits.hh"
52#include "G4Neutron.hh"
53#include "G4NucleiProperties.hh"
54#include "G4AutoLock.hh"
55
56namespace
57{
59}
60
62 const G4String& nameData,
63 const G4String& nameDir, G4double emaxHP,
64 G4int zmin, G4int zmax)
65 : G4VCrossSectionDataSet(nameData),
66 fParticle(p),
67 fNeutron(G4Neutron::Neutron()),
68 fManagerHP(G4ParticleHPManager::GetInstance()),
69 emax(emaxHP),
70 emaxT(fManagerHP->GetMaxEnergyDoppler()),
71 elimit(1.0e-11*CLHEP::eV),
72 logElimit(G4Log(elimit)),
73 minZ(zmin),
74 maxZ(zmax),
75 fDataName(nameData),
76 fDataDirectory(nameDir)
77{
78 // verboseLevel = 1;
79 if (verboseLevel > 1) {
80 G4cout << "G4CrossSectionHP::G4CrossSectionHP: Initialise for "
81 << fDataName << " " << minZ << " < Z < " << maxZ
82 << " EmaxT(MeV)=" << emaxT << G4endl;
83 G4cout << "Data directory: " << fDataDirectory << G4endl;
84 }
86 auto data = ptr->GetElementDataByName(fDataName);
87 if (nullptr == data) {
88 data = new G4ElementData(maxZ - minZ + 1);
89 data->SetName(fDataName);
90 }
91 fData = data;
92}
93
95 G4int, G4int,
96 const G4Element*,
97 const G4Material*)
98{
99 return (dp->GetKineticEnergy() <= emax);
100}
101
103 G4int Z, G4int A,
104 const G4Isotope*,
105 const G4Element*,
106 const G4Material* mat)
107{
108 G4double ekin = dp->GetKineticEnergy();
109 G4double loge = dp->GetLogKineticEnergy();
110 if (ekin < elimit) {
111 ekin = elimit;
112 loge = logElimit;
113 }
114 if (mat != fCurrentMat) { PrepareCache(mat); }
115
116 return IsoCrossSection(ekin, loge, Z, A, mat->GetTemperature());
117}
118
122 G4int Z, G4int A,
123 const G4Isotope*,
124 const G4Element*,
125 const G4Material* mat)
126{
127 G4double ekin = e;
128 G4double loge = le;
129 if (ekin < elimit) {
130 ekin = elimit;
131 loge = logElimit;
132 }
133 if (mat != fCurrentMat) { PrepareCache(mat); }
134
135 return IsoCrossSection(ekin, loge, Z, A, mat->GetTemperature());
136}
137
138G4double G4CrossSectionHP::IsoCrossSection(const G4double ekin,
139 const G4double logek,
140 const G4int Z, const G4int A,
141 const G4double T)
142{
143 //G4cout << "G4CrossSectionHP::IsoCrossSection Z=" << Z << " A=" << A
144 // << " E(MeV)=" << ekin/MeV << " T=" << T << " " << GetName() << G4endl;
145 G4double xs = 0.0;
146 if (ekin > emax || Z > maxZ || Z < minZ || ekin < elimit) { return xs; }
147
148 const G4PhysicsVector* pv0 = fData->GetElementData(Z - minZ);
149 if (nullptr == pv0) {
150 InitialiseOnFly(Z);
151 pv0 = fData->GetElementData(Z - minZ);
152 }
153 if (nullptr == pv0) { return xs; }
154 const G4PhysicsVector* pv = fData->GetComponentDataByID(Z - minZ, A);
155 if (nullptr == pv) { return xs; }
156
157 // no Doppler broading
158 G4double factT = T/CLHEP::STP_Temperature;
159 if (ekin >= emaxT*factT || fManagerHP->GetNeglectDoppler()) {
160 xs = pv->LogFreeVectorValue(ekin, logek);
161
162 } else {
163
164 // Doppler broading
165 G4double e0 = CLHEP::k_Boltzmann*T;
166 G4double mass = fParticle->GetPDGMass();
168
169 // projectile
170 G4LorentzVector lv(0., 0., std::sqrt(ekin*(ekin + 2*mass)), mass + ekin);
171
172 // limits of integration
173 const G4double lim = 1.01;
174 const G4int nmin = 3;
175 G4int ii = 0;
176 const G4int nn = 20;
177 G4double xs2 = 0.0;
178
179 for (G4int i=0; i<nn; ++i) {
180 G4double erand = G4RandGamma::shoot(2.0, e0);
181 auto mom = G4RandomDirection()*std::sqrt(2*massTarget*erand);
182 fLV.set(mom.x(), mom.y(), mom.z(), massTarget + erand);
183 fBoost = fLV.boostVector();
184 fLV = lv.boost(fBoost);
185 if (fLV.pz() <= 0.0) { continue; }
186 ++ii;
187 G4double e = fLV.e() - mass;
188 G4double y = pv->Value(e, index);
189 xs += y;
190 xs2 += y*y;
191 if (ii >= nmin && ii*xs2 <= lim*xs*xs) { break; }
192 }
193 if (ii > 0) { xs /= (G4double)(ii); }
194 }
195#ifdef G4VERBOSE
196 if (verboseLevel > 1) {
197 G4cout << "G4CrossSectionHP::IsoXS " << fDataName
198 << " Z= " << Z << " A= " << A
199 << " Ekin(MeV)= " << ekin/MeV << " xs(b)= " << xs/barn
200 << " size=" << fZA.size() << G4endl;
201 }
202#endif
203
204 // save cross section into struct
205 for (std::size_t i=0; i<fZA.size(); ++i) {
206 if (Z == fZA[i].first && A == fZA[i].second) {
207 fIsoXS[i] = xs;
208 break;
209 }
210 }
211 return xs;
212}
213
216{
217 std::size_t nIso = elm->GetNumberOfIsotopes();
218 const G4Isotope* iso = elm->GetIsotope(0);
219
220 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
221 if(1 == nIso) { return iso; }
222
223 // more than 1 isotope
224 G4int Z = elm->GetZasInt();
225 if (Z >= minZ && Z <= maxZ && nullptr == fData->GetElementData(Z - minZ)) {
226 InitialiseOnFly(Z);
227 }
228
229 const G4double* abundVector = elm->GetRelativeAbundanceVector();
231 G4double sum = 0.0;
232
233 // is there cached isotope wise cross section?
234 std::size_t j;
235 if (Z < minZ || Z > maxZ || !CheckCache(Z) ||
236 0 == fData->GetNumberOfComponents(Z - minZ)) {
237 for (j = 0; j<nIso; ++j) {
238 sum += abundVector[j];
239 if(q <= sum) {
240 iso = elm->GetIsotope((G4int)j);
241 break;
242 }
243 }
244 return iso;
245 }
246 std::size_t nn = fTemp.size();
247 if (nn < nIso) { fTemp.resize(nIso, 0.); }
248
249 // reuse cache
250 for (j=0; j<nIso; ++j) {
251 sum += abundVector[j]*
252 GetCrossSection(Z - minZ, elm->GetIsotope((G4int)j)->GetN());
253 fTemp[j] = sum;
254 }
255 sum *= q;
256 for (j = 0; j<nIso; ++j) {
257 if (fTemp[j] >= sum) {
258 iso = elm->GetIsotope((G4int)j);
259 break;
260 }
261 }
262 return iso;
263}
264
266{
267 if (verboseLevel > 1){
268 G4cout << "G4CrossSectionHP::BuildPhysicsTable for "
269 << p.GetParticleName() << " and " << fDataName << G4endl;
270 }
271
272 // it is possible re-initialisation for the second run
274
275 // Access to elements
276 for ( auto const & elm : *table ) {
277 G4int Z = elm->GetZasInt();
278 if (Z >= minZ && Z <= maxZ &&
279 nullptr == fData->GetElementData(Z - minZ)) {
280 InitialiseOnFly(Z);
281 }
282 }
283
284 // prepare isotope selection
285 std::size_t nmax = 0;
286 std::size_t imax = 0;
287 for ( auto const & mat : *G4Material::GetMaterialTable() ) {
288 std::size_t n = 0;
289 for ( auto const & elm : *mat->GetElementVector() ) {
290 std::size_t niso = elm->GetNumberOfIsotopes();
291 n += niso;
292 if(niso > imax) { imax = niso; }
293 }
294 if (n > nmax) { nmax = n; }
295 }
296 fTemp.resize(imax, 0.0);
297 fZA.clear();
298 fZA.reserve(nmax);
299 fIsoXS.resize(nmax, 0.0);
300}
301
303{
304 if (fManagerHP->GetVerboseLevel() == 0 || fPrinted)
305 return;
306
307 //
308 // Dump element based cross section
309 // range 10e-5 eV to 20 MeV
310 // 10 point per decade
311 // in barn
312 //
313 fPrinted = true;
314 G4cout << G4endl;
315 G4cout << "HP Cross Section " << fDataName << " for "
316 << fParticle->GetParticleName() << G4endl;
317 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
318 G4cout << G4endl;
319 G4cout << "Name of Element" << G4endl;
320 G4cout << "Energy[eV] XS[barn]" << G4endl;
321 G4cout << G4endl;
322
324 for ( auto const & elm : *table ) {
325 G4int Z = elm->GetZasInt();
326 if (Z >= minZ && Z <= maxZ && nullptr != fData->GetElementData(Z - minZ)) {
327 G4cout << "---------------------------------------------------" << G4endl;
328 G4cout << elm->GetName() << G4endl;
329 std::size_t n = fData->GetNumberOfComponents(Z);
330 for (size_t i=0; i<n; ++i) {
331 G4cout << "## Z=" << Z << " A=" << fData->GetComponentID(Z - minZ, i);
332 G4cout << *(fData->GetComponentDataByIndex(Z - minZ, i)) << G4endl;
333 }
334 }
335 }
336}
337
338void G4CrossSectionHP::PrepareCache(const G4Material* mat)
339{
340 fCurrentMat = mat;
341 fZA.clear();
342 for ( auto const & elm : *(mat->GetElementVector()) ) {
343 G4int Z = elm->GetZasInt();
344 for ( auto const & iso : *(elm->GetIsotopeVector()) ) {
345 G4int A = iso->GetN();
346 fZA.emplace_back(Z, A);
347 }
348 }
349 fIsoXS.resize(fZA.size(), 0.0);
350}
351
352void G4CrossSectionHP::InitialiseOnFly(const G4int Z)
353{
354 G4AutoLock l(&theHPXS);
355 Initialise(Z);
356 l.unlock();
357}
358
359void G4CrossSectionHP::Initialise(const G4int Z)
360{
361 if (fManagerHP->GetVerboseLevel() > 1) {
362 G4cout << " G4CrossSectionHP::Initialise: Z=" << Z
363 << " for " << fDataName
364 << " minZ=" << minZ << " maxZ=" << maxZ << G4endl;
365 }
366 if (Z < minZ || Z > maxZ || nullptr != fData->GetElementData(Z - minZ)) {
367 return;
368 }
369
370 // add empty vector to avoid double initialisation
371 fData->InitialiseForElement(Z - minZ, new G4PhysicsVector());
372
373 G4String tnam = "temp";
374 G4bool noComp = true;
375 for (G4int A=amin[Z]; A<=amax[Z]; ++A) {
376 std::ostringstream ost;
377 ost << fDataDirectory;
378 // first check special cases
379 if (6 == Z && 12 == A && fParticle == fNeutron) {
380 ost << Z << "_nat_" << elementName[Z];
381 } else if (18 == Z && 40 != A) {
382 continue;
383 } else if (27 == Z && 62 == A) {
384 ost << Z << "_62m1_" << elementName[Z];
385 } else if (47 == Z && 106 == A) {
386 ost << Z << "_106m1_" << elementName[Z];
387 } else if (48 == Z && 115 == A) {
388 ost << Z << "_115m1_" << elementName[Z];
389 } else if (52 == Z && 127 == A) {
390 ost << Z << "_127m1_" << elementName[Z];
391 } else if (52 == Z && 129 == A) {
392 ost << Z << "_129m1_" << elementName[Z];
393 } else if (52 == Z && 131 == A) {
394 ost << Z << "_131m1_" << elementName[Z];
395 } else if (61 == Z && 145 == A) {
396 ost << Z << "_147_" << elementName[Z];
397 } else if (67 == Z && 166 == A) {
398 ost << Z << "_166m1_" << elementName[Z];
399 } else if (73 == Z && 180 == A) {
400 ost << Z << "_180m1_" << elementName[Z];
401 } else if ((Z == 85 && A == 210) || (Z == 86 && A == 222) || (Z == 87 && A == 223)) {
402 ost << "84_209_" << elementName[84];
403 } else {
404 // the main file name
405 ost << Z << "_" << A << "_" << elementName[Z];
406 }
407 std::istringstream theXSData(tnam, std::ios::in);
408 fManagerHP->GetDataStream(ost.str().c_str(), theXSData);
409 if (theXSData) {
410 G4int i1, i2, n;
411 theXSData >> i1 >> i2 >> n;
412 if (fManagerHP->GetVerboseLevel() > 1) {
413 G4cout << "## G4CrossSectionHP::Initialise for Z=" << Z
414 << " A=" << A << " Npoints=" << n << G4endl;
415 }
416 G4double x, y;
418 for (G4int i=0; i<n; ++i) {
419 theXSData >> x >> y;
420 x *= CLHEP::eV;
421 y *= CLHEP::barn;
422 //G4cout << " e=" << x << " xs=" << y << G4endl;
423 v->PutValues((std::size_t)i, x, y);
424 }
425 v->EnableLogBinSearch(binSearch);
426 if (noComp) {
427 G4int nmax = amax[Z] - A + 1;
428 fData->InitialiseForComponent(Z - minZ, nmax);
429 noComp = false;
430 }
431 fData->AddComponent(Z - minZ, A, v);
432 }
433 }
434 if (noComp) { fData->InitialiseForComponent(Z - minZ, 0); }
435}
436
std::vector< G4Element * > G4ElementTable
G4double G4Log(G4double x)
Definition G4Log.hh:227
G4ThreeVector G4RandomDirection()
#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
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void set(double x, double y, double z, double t)
G4double ComputeIsoCrossSection(G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) override
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionHP(const G4ParticleDefinition *, const G4String &nameData, const G4String &nameDir, G4double emaxHP, G4int zmin, G4int zmax)
void DumpPhysicsTable(const G4ParticleDefinition &) override
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) override
const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy, G4double logE) override
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
static G4ElementDataRegistry * Instance()
G4PhysicsVector * GetElementData(G4int Z) const
void InitialiseForComponent(G4int Z, G4int nComponents=0)
G4PhysicsVector * GetComponentDataByIndex(G4int Z, std::size_t idx) const
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
G4int GetComponentID(G4int Z, std::size_t idx) const
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
G4PhysicsVector * GetComponentDataByID(G4int Z, G4int id) const
std::size_t GetNumberOfComponents(G4int Z) const
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
G4int GetN() const
Definition G4Isotope.hh:83
const G4ElementVector * GetElementVector() const
G4double GetTemperature() const
static G4MaterialTable * GetMaterialTable()
static G4double GetNuclearMass(const G4double A, const G4double Z)
const G4String & GetParticleName() const
void GetDataStream(const G4String &, std::istringstream &iss)
G4bool GetNeglectDoppler() const
void PutValues(const std::size_t index, const G4double energy, const G4double value)
void EnableLogBinSearch(const G4int n=1)
G4double LogFreeVectorValue(const G4double energy, const G4double theLogEnergy) const
G4double Value(const G4double energy, std::size_t &lastidx) const