Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPInelasticData.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// particle_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 070523 add neglecting doppler broadening on the fly. T. Koi
31// 070613 fix memory leaking by T. Koi
32// 071002 enable cross section dump by T. Koi
33// 080428 change checking point of "neglecting doppler broadening" flag
34// from GetCrossSection to BuildPhysicsTable by T. Koi
35// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
36//
37// P. Arce, June-2014 Conversion neutron_hp to particle_hp
38//
40
41#include "G4ElementTable.hh"
43#include "G4Neutron.hh"
44#include "G4NucleiProperties.hh"
45#include "G4ParticleHPData.hh"
47#include "G4Pow.hh"
48
51{
52 const char* dataDirVariable;
53 G4String particleName;
54 if (projectile == G4Neutron::Neutron()) {
55 dataDirVariable = "G4NEUTRONHPDATA";
56 }
57 else if (projectile == G4Proton::Proton()) {
58 dataDirVariable = "G4PROTONHPDATA";
59 particleName = "Proton";
60 }
61 else if (projectile == G4Deuteron::Deuteron()) {
62 dataDirVariable = "G4DEUTERONHPDATA";
63 particleName = "Deuteron";
64 }
65 else if (projectile == G4Triton::Triton()) {
66 dataDirVariable = "G4TRITONHPDATA";
67 particleName = "Triton";
68 }
69 else if (projectile == G4He3::He3()) {
70 dataDirVariable = "G4HE3HPDATA";
71 particleName = "He3";
72 }
73 else if (projectile == G4Alpha::Alpha()) {
74 dataDirVariable = "G4ALPHAHPDATA";
75 particleName = "Alpha";
76 }
77 else {
78 G4String message(
79 "G4ParticleHPInelasticData may only be called for neutron, proton, deuteron, triton, He3 or "
80 "alpha, while it is called for "
81 + projectile->GetParticleName());
82 throw G4HadronicException(__FILE__, __LINE__, message.c_str());
83 }
84 G4String dataName = projectile->GetParticleName() + "HPInelasticXS";
85 dataName.at(0) = (char)std::toupper(dataName.at(0));
86 SetName(dataName);
87
88 if ((G4FindDataDir(dataDirVariable) == nullptr) && (G4FindDataDir("G4PARTICLEHPDATA") == nullptr))
89 {
90 G4String message("Please setenv G4PARTICLEHPDATA (recommended) or, at least setenv "
91 + G4String(dataDirVariable) + " to point to the "
92 + projectile->GetParticleName() + " cross-section files.");
93 throw G4HadronicException(__FILE__, __LINE__, message.c_str());
94 }
95
96 G4String dirName;
97 if (G4FindDataDir(dataDirVariable) != nullptr) {
98 dirName = G4FindDataDir(dataDirVariable);
99 }
100 else {
101 G4String baseName = G4FindDataDir("G4PARTICLEHPDATA");
102 dirName = baseName + "/" + particleName;
103 }
104#ifdef G4VERBOSE
106 G4cout << "@@@ G4ParticleHPInelasticData instantiated for particle "
107 << projectile->GetParticleName() << " data directory variable is " << dataDirVariable
108 << " pointing to " << dirName << G4endl;
109 }
110#endif
111
112 SetMinKinEnergy(0 * CLHEP::MeV);
113 SetMaxKinEnergy(20 * CLHEP::MeV);
114
115 theCrossSections = nullptr;
116 theProjectile = projectile;
117
118 theHPData = nullptr;
119 instanceOfWorker = false;
121 theHPData = new G4ParticleHPData(theProjectile);
122 }
123 else {
124 instanceOfWorker = true;
125 }
126 element_cache = nullptr;
127 material_cache = nullptr;
128 ke_cache = 0.0;
129 xs_cache = 0.0;
130}
131
133{
134 if (theCrossSections != nullptr && !instanceOfWorker) {
135 theCrossSections->clearAndDestroy();
136 delete theCrossSections;
137 theCrossSections = nullptr;
138 }
139 if (theHPData != nullptr && !instanceOfWorker) {
140 delete theHPData;
141 theHPData = nullptr;
142 }
143}
144
146 G4int /*A*/, const G4Element* /*elm*/,
147 const G4Material* /*mat*/)
148{
149 G4double eKin = dp->GetKineticEnergy();
150 return eKin <= GetMaxKinEnergy() && eKin >= GetMinKinEnergy()
151 && dp->GetDefinition() == theProjectile;
152}
153
155 G4int /*A*/, const G4Isotope* /*iso*/,
156 const G4Element* element,
157 const G4Material* material)
158{
159 if (dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache)
160 return xs_cache;
161
162 ke_cache = dp->GetKineticEnergy();
163 element_cache = element;
164 material_cache = material;
165 G4double xs = GetCrossSection(dp, element, material->GetTemperature());
166 xs_cache = xs;
167 return xs;
168}
169
171{
173 theCrossSections = G4ParticleHPManager::GetInstance()->GetInelasticCrossSections(&projectile);
174 return;
175 }
176 if (theHPData == nullptr)
177 theHPData = G4ParticleHPData::Instance(const_cast<G4ParticleDefinition*>(&projectile));
178
179 std::size_t numberOfElements = G4Element::GetNumberOfElements();
180 if (theCrossSections == nullptr)
181 theCrossSections = new G4PhysicsTable(numberOfElements);
182 else
183 theCrossSections->clearAndDestroy();
184
185 // make a PhysicsVector for each element
186
187 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
188 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
189 for (std::size_t i = 0; i < numberOfElements; ++i) {
190 G4PhysicsVector* physVec = theHPData->MakePhysicsVector((*theElementTable)[i], this);
191 theCrossSections->push_back(physVec);
192 }
194}
195
197{
198 if (&projectile != theProjectile)
199 throw G4HadronicException(__FILE__, __LINE__,
200 "Attempt to use ParticleHP data for a wrong projectile!!!");
201
202#ifdef G4VERBOSE
204
205 // Dump element based cross section
206 // range 10e-5 eV to 20 MeV
207 // 10 point per decade
208 // in barn
209
210 G4cout << G4endl;
211 G4cout << G4endl;
212 G4cout << "Inelastic Cross Section of Neutron HP" << G4endl;
213 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
214 G4cout << G4endl;
215 G4cout << "Name of Element" << G4endl;
216 G4cout << "Energy[eV] XS[barn]" << G4endl;
217 G4cout << G4endl;
218
219 std::size_t numberOfElements = G4Element::GetNumberOfElements();
220 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
221 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
222
223 for (std::size_t i = 0; i < numberOfElements; ++i) {
224 G4cout << (*theElementTable)[i]->GetName() << G4endl;
225
226 G4int ie = 0;
227
228 for (ie = 0; ie < 130; ie++) {
229 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * CLHEP::eV;
230 G4bool outOfRange = false;
231
232 if (eKinetic < 20 * CLHEP::MeV) {
233 G4cout << eKinetic / CLHEP::eV << " "
234 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / CLHEP::barn
235 << G4endl;
236 }
237 }
238 G4cout << G4endl;
239 }
240#endif
241}
242
244 const G4Element* anE, G4double aT)
245{
246 G4double result = 0;
247 G4bool outOfRange;
248 std::size_t index = anE->GetIndex();
249
250 // prepare neutron
251 G4double eKinetic = projectile->GetKineticEnergy();
252
253 if (G4ParticleHPManager::GetInstance()->GetNeglectDoppler()) {
254 // NEGLECT_DOPPLER
255 G4double factor = 1.0;
256 if (eKinetic < aT * CLHEP::k_Boltzmann) {
257 // below 0.1 eV neutrons
258 // Have to do some, but now just igonre.
259 // Will take care after performance check.
260 // factor = factor * targetV;
261 }
262 return ((*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange)) * factor;
263 }
264
265 G4ReactionProduct theNeutron(projectile->GetDefinition());
266 theNeutron.SetMomentum(projectile->GetMomentum());
267 theNeutron.SetKineticEnergy(eKinetic);
268
269 // prepare thermal nucleus
270 G4Nucleus aNuc;
271 G4double eps = 0.0001;
272 G4double theA = anE->GetN();
273 G4double theZ = anE->GetZ();
274 G4double eleMass;
275 eleMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA + eps),
276 static_cast<G4int>(theZ + eps));
277
278 G4ReactionProduct boosted;
279 G4double aXsection;
280
281 // MC integration loop
282 G4int counter = 0;
283 G4int failCount = 0;
284 G4double buffer = 0;
285 G4int size = G4int(std::max(10., aT / 60 * CLHEP::kelvin));
286 G4ThreeVector neutronVelocity = 1. / theProjectile->GetPDGMass() * theNeutron.GetMomentum();
287 G4double neutronVMag = neutronVelocity.mag();
288
289#ifndef G4PHP_DOPPLER_LOOP_ONCE
290 while (counter == 0
291 || std::abs(buffer - result / std::max(1, counter))
292 > 0.01 * buffer) // Loop checking, 11.05.2015, T. Koi
293 {
294 if (counter != 0) buffer = result / counter;
295 while (counter < size) // Loop checking, 11.05.2015, T. Koi
296 {
297 ++counter;
298#endif
299 G4ReactionProduct aThermalNuc =
300 aNuc.GetThermalNucleus(eleMass / G4Neutron::Neutron()->GetPDGMass(), aT);
301 boosted.Lorentz(theNeutron, aThermalNuc);
302 G4double theEkin = boosted.GetKineticEnergy();
303 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
304 if (aXsection < 0) {
305 if (failCount < 1000) {
306 ++failCount;
307#ifndef G4PHP_DOPPLER_LOOP_ONCE
308 --counter;
309 continue;
310#endif
311 }
312 else {
313 aXsection = 0;
314 }
315 }
316 // velocity correction.
317 G4ThreeVector targetVelocity = 1. / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
318 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
319 result += aXsection;
320 }
321#ifndef G4PHP_DOPPLER_LOOP_ONCE
322 size += size;
323 }
324 result /= counter;
325#endif
326
327 return result;
328}
329
334
339
341{
342 outFile << "Extension of High Precision cross section for inelastic reaction of proton, "
343 "deuteron, triton, He3 and alpha below 20MeV\n";
344}
std::vector< G4Element * > G4ElementTable
const char * G4FindDataDir(const char *)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double mag() const
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
G4double GetZ() const
Definition G4Element.hh:119
static size_t GetNumberOfElements()
Definition G4Element.cc:393
size_t GetIndex() const
Definition G4Element.hh:159
G4double GetN() const
Definition G4Element.hh:123
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition G4He3.cc:90
G4double GetTemperature() const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition G4Nucleus.cc:248
const G4String & GetParticleName() const
G4PhysicsVector * MakePhysicsVector(const G4Element *thE, G4ParticleHPFissionData *theP)
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
void CrossSectionDescription(std::ostream &) const override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4ParticleHPInelasticData(G4ParticleDefinition *projectile=G4Neutron::Neutron())
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *) override
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *) override
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
void DumpPhysicsTable(const G4ParticleDefinition &) override
G4PhysicsTable * GetInelasticCrossSections(const G4ParticleDefinition *part) const
void RegisterInelasticCrossSections(const G4ParticleDefinition *part, G4PhysicsTable *ptr)
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
void clearAndDestroy()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230
static G4Proton * Proton()
Definition G4Proton.cc:90
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetKineticEnergy(const G4double en)
G4double GetMass() const
static G4Triton * Triton()
Definition G4Triton.cc:90
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetName(const G4String &nam)
G4bool IsWorkerThread()
G4bool IsMasterThread()
#define G4ThreadLocal
Definition tls.hh:77