Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPFissionData.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// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 070618 fix memory leaking by T. Koi
31// 071002 enable cross section dump by T. Koi
32// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
33// 081124 Protect invalid read which caused run time errors by T. Koi
34// 100729 Add safty for 0 lenght cross sections by T. Ko
35// P. Arce, June-2014 Conversion neutron_hp to particle_hp
36//
38
39#include "G4ElementTable.hh"
41#include "G4Neutron.hh"
42#include "G4NucleiProperties.hh"
43#include "G4ParticleHPData.hh"
46#include "G4Pow.hh"
47#include "G4SystemOfUnits.hh"
48
50{
51 SetMinKinEnergy(0 * MeV);
52 SetMaxKinEnergy(20 * MeV);
53
54 theCrossSections = nullptr;
55 instanceOfWorker = false;
57 instanceOfWorker = true;
58 }
59 element_cache = nullptr;
60 material_cache = nullptr;
61 ke_cache = 0.0;
62 xs_cache = 0.0;
63 // BuildPhysicsTable(*G4Neutron::Neutron());
64}
65
67{
68 if (theCrossSections != nullptr && !instanceOfWorker) {
69 theCrossSections->clearAndDestroy();
70 delete theCrossSections;
71 theCrossSections = nullptr;
72 }
73}
74
76 G4int /*A*/, const G4Element* /*elm*/,
77 const G4Material* /*mat*/)
78{
79 G4double eKin = dp->GetKineticEnergy();
80 return eKin <= GetMaxKinEnergy() && eKin >= GetMinKinEnergy()
82}
83
85 G4int /*A*/, const G4Isotope* /*iso*/,
86 const G4Element* element,
87 const G4Material* material)
88{
89 if (dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache)
90 return xs_cache;
91
92 ke_cache = dp->GetKineticEnergy();
93 element_cache = element;
94 material_cache = material;
95 G4double xs = GetCrossSection(dp, element, material->GetTemperature());
96 xs_cache = xs;
97 return xs;
98}
99
101{
102 if (&aP != G4Neutron::Neutron())
103 throw G4HadronicException(__FILE__, __LINE__,
104 "Attempt to use NeutronHP data for particles other than neutrons!!!");
105
108 return;
109 }
110
111 std::size_t numberOfElements = G4Element::GetNumberOfElements();
112 if (theCrossSections == nullptr)
113 theCrossSections = new G4PhysicsTable(numberOfElements);
114 else
115 theCrossSections->clearAndDestroy();
116
117 // make a PhysicsVector for each element
118
119 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
120 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
121 for (std::size_t i = 0; i < numberOfElements; ++i) {
123 ->MakePhysicsVector((*theElementTable)[i], this);
124 theCrossSections->push_back(physVec);
125 }
126
128}
129
131{
132 if (&aP != G4Neutron::Neutron())
133 throw G4HadronicException(__FILE__, __LINE__,
134 "Attempt to use NeutronHP data for particles other than neutrons!!!");
135
136#ifdef G4VERBOSE
138
139 //
140 // Dump element based cross section
141 // range 10e-5 eV to 20 MeV
142 // 10 point per decade
143 // in barn
144 //
145 G4cout << G4endl;
146 G4cout << G4endl;
147 G4cout << "Fission Cross Section of Neutron HP" << G4endl;
148 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
149 G4cout << G4endl;
150 G4cout << "Name of Element" << G4endl;
151 G4cout << "Energy[eV] XS[barn]" << G4endl;
152 G4cout << G4endl;
153
154 std::size_t numberOfElements = G4Element::GetNumberOfElements();
155 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
156 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
157
158 for (std::size_t i = 0; i < numberOfElements; ++i) {
159 G4cout << (*theElementTable)[i]->GetName() << G4endl;
160
161 if ((*((*theCrossSections)(i))).GetVectorLength() == 0) {
162 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
163 G4cout << G4endl;
164 continue;
165 }
166
167 for (G4int ie = 0; ie < 130; ++ie) {
168 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * eV;
169 G4bool outOfRange = false;
170
171 if (eKinetic < 20 * MeV) {
172 G4cout << eKinetic / eV << " "
173 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / barn << G4endl;
174 }
175 }
176
177 G4cout << G4endl;
178 }
179#endif
180}
181
183 G4double aT)
184{
185 G4double result = 0;
186 if (anE->GetZ() < 88) return result;
187 G4bool outOfRange;
188 auto index = (G4int)anE->GetIndex();
189
190 if (((*theCrossSections)(index))->GetVectorLength() == 0) return result;
191
192 // prepare neutron
193 G4double eKinetic = aP->GetKineticEnergy();
194 G4ReactionProduct theNeutronRP(aP->GetDefinition());
195 theNeutronRP.SetMomentum(aP->GetMomentum());
196 theNeutronRP.SetKineticEnergy(eKinetic);
197
198 if (G4ParticleHPManager::GetInstance()->GetNeglectDoppler()) {
199 // NEGLECT_DOPPLER
200 G4double factor = 1.0;
201 if (eKinetic < aT * k_Boltzmann) {
202 // below 0.1 eV neutrons
203 // Have to do some, but now just igonre.
204 // Will take care after performance check.
205 // factor = factor * targetV;
206 }
207 return ((*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange)) * factor;
208 }
209
210 // prepare thermal nucleus
211 G4Nucleus aNuc;
212 G4double eps = 0.0001;
213 G4double theA = anE->GetN();
214 G4double theZ = anE->GetZ();
215 G4double eleMass;
216 eleMass = (G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA + eps),
217 static_cast<G4int>(theZ + eps)))
219
220 G4ReactionProduct boosted;
221 G4double aXsection;
222
223 // MC integration loop
224 G4int counter = 0;
225 G4double buffer = 0;
226 G4int size = G4int(std::max(10., aT / 60 * kelvin));
227 G4ThreeVector neutronVelocity =
228 1. / G4Neutron::Neutron()->GetPDGMass() * theNeutronRP.GetMomentum();
229 G4double neutronVMag = neutronVelocity.mag();
230
231 while (counter == 0
232 || std::abs(buffer - result / std::max(1, counter))
233 > 0.01 * buffer) // Loop checking, 11.05.2015, T. Koi
234 {
235 if (counter != 0) buffer = result / counter;
236 while (counter < size) // Loop checking, 11.05.2015, T. Koi
237 {
238 counter++;
239 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
240 boosted.Lorentz(theNeutronRP, aThermalNuc);
241 G4double theEkin = boosted.GetKineticEnergy();
242 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
243 // velocity correction.
244 G4ThreeVector targetVelocity = 1. / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
245 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
246 result += aXsection;
247 }
248 size += size;
249 }
250 result /= counter;
251 return result;
252}
253
258
263
265{
266 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF)\n"
267 << "for induced fission reaction of neutrons below 20MeV\n";
268}
std::vector< G4Element * > G4ElementTable
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
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()
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
G4PhysicsVector * MakePhysicsVector(const G4Element *thE, G4ParticleHPFissionData *theP)
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *) override
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
void CrossSectionDescription(std::ostream &) const override
void SetVerboseLevel(G4int) override
void DumpPhysicsTable(const G4ParticleDefinition &) override
G4PhysicsTable * GetFissionCrossSections() const
static G4ParticleHPManager * GetInstance()
void RegisterFissionCrossSections(G4PhysicsTable *val)
void push_back(G4PhysicsVector *)
void clearAndDestroy()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230
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
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()
#define G4ThreadLocal
Definition tls.hh:77