Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPElasticData.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// 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 "G4Electron.hh"
42#include "G4ElementTable.hh"
44#include "G4Neutron.hh"
45#include "G4NucleiProperties.hh"
46#include "G4Nucleus.hh"
47#include "G4ParticleHPData.hh"
50#include "G4Pow.hh"
51#include "G4SystemOfUnits.hh"
52
54{
55 SetMinKinEnergy(0 * MeV);
56 SetMaxKinEnergy(20 * MeV);
57
58 theCrossSections = nullptr;
59 instanceOfWorker = false;
61 instanceOfWorker = true;
62 }
63 element_cache = nullptr;
64 material_cache = nullptr;
65 ke_cache = 0.0;
66 xs_cache = 0.0;
67 // BuildPhysicsTable( *G4Neutron::Neutron() );
68}
69
71{
72 if (theCrossSections != nullptr && !instanceOfWorker) {
73 theCrossSections->clearAndDestroy();
74 delete theCrossSections;
75 theCrossSections = nullptr;
76 }
77}
78
80 G4int /*A*/, const G4Element* /*elm*/,
81 const G4Material* /*mat*/)
82{
83 G4double eKin = dp->GetKineticEnergy();
84 return eKin <= GetMaxKinEnergy() && eKin >= GetMinKinEnergy()
86}
87
89 G4int /*A*/, const G4Isotope* /*iso*/,
90 const G4Element* element,
91 const G4Material* material)
92{
93 if (dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache)
94 return xs_cache;
95
96 ke_cache = dp->GetKineticEnergy();
97 element_cache = element;
98 material_cache = material;
99 G4double xs = GetCrossSection(dp, element, material->GetTemperature());
100 xs_cache = xs;
101 return xs;
102}
103
104/*
105G4bool G4ParticleHPElasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
106{
107 G4bool result = true;
108 G4double eKin = aP->GetKineticEnergy();
109 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
110 return result;
111}
112*/
113
115{
116 if (&aP != G4Neutron::Neutron())
117 throw G4HadronicException(__FILE__, __LINE__,
118 "Attempt to use NeutronHP data for particles other than neutrons!!!");
119
122 return;
123 }
124
125 std::size_t numberOfElements = G4Element::GetNumberOfElements();
126 // TKDB
127 // if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
128 if (theCrossSections == nullptr)
129 theCrossSections = new G4PhysicsTable(numberOfElements);
130 else
131 theCrossSections->clearAndDestroy();
132
133 // make a PhysicsVector for each element
134
135 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
136 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
137 for (std::size_t i = 0; i < numberOfElements; ++i) {
139 ->MakePhysicsVector((*theElementTable)[i], this);
140 theCrossSections->push_back(physVec);
141 }
142
144}
145
147{
148 if (&aP != G4Neutron::Neutron())
149 throw G4HadronicException(__FILE__, __LINE__,
150 "Attempt to use NeutronHP data for particles other than neutrons!!!");
151
152#ifdef G4VERBOSE
154
155 //
156 // Dump element based cross section
157 // range 10e-5 eV to 20 MeV
158 // 10 point per decade
159 // in barn
160 //
161
162 G4cout << G4endl;
163 G4cout << G4endl;
164 G4cout << "Elastic Cross Section of Neutron HP" << G4endl;
165 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
166 G4cout << G4endl;
167 G4cout << "Name of Element" << G4endl;
168 G4cout << "Energy[eV] XS[barn]" << G4endl;
169 G4cout << G4endl;
170
171 std::size_t numberOfElements = G4Element::GetNumberOfElements();
172 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
173 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
174
175 for (std::size_t i = 0; i < numberOfElements; ++i) {
176 G4cout << (*theElementTable)[i]->GetName() << G4endl;
177 G4int ie = 0;
178
179 for (ie = 0; ie < 130; ++ie) {
180 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * eV;
181 G4bool outOfRange = false;
182
183 if (eKinetic < 20 * MeV) {
184 G4cout << eKinetic / eV << " "
185 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / barn << G4endl;
186 }
187 }
188 G4cout << G4endl;
189 }
190#endif
191}
192
194 G4double aT)
195{
196 G4double result = 0;
197 G4bool outOfRange;
198 auto index = (G4int)anE->GetIndex();
199
200 // prepare neutron
201 G4double eKinetic = aP->GetKineticEnergy();
202
203 if (G4ParticleHPManager::GetInstance()->GetNeglectDoppler()) {
204 // NEGLECT_DOPPLER_B.
205 G4double factor = 1.0;
206 if (eKinetic < aT * k_Boltzmann) {
207 // below 0.1 eV neutrons
208 // Have to do some, but now just igonre.
209 // Will take care after performance check.
210 // factor = factor * targetV;
211 }
212 return ((*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange)) * factor;
213 }
214
215 G4ReactionProduct theNeutron(aP->GetDefinition());
216 theNeutron.SetMomentum(aP->GetMomentum());
217 theNeutron.SetKineticEnergy(eKinetic);
218
219 // prepare thermal nucleus
220 G4Nucleus aNuc;
221 G4double eps = 0.0001;
222 G4double theA = anE->GetN();
223 G4double theZ = anE->GetZ();
224 G4double eleMass;
225
226 eleMass = (G4NucleiProperties::GetNuclearMass(G4int(theA + eps), G4int(theZ + eps)))
228
229 G4ReactionProduct boosted;
230 G4double aXsection;
231
232 // MC integration loop
233 G4int counter = 0;
234 G4double buffer = 0;
235 G4int size = G4int(std::max(10., aT / 60 * kelvin));
236 G4ThreeVector neutronVelocity =
237 1. / G4Neutron::Neutron()->GetPDGMass() * theNeutron.GetMomentum();
238 G4double neutronVMag = neutronVelocity.mag();
239
240 while (counter == 0
241 || std::abs(buffer - result / std::max(1, counter))
242 > 0.03 * buffer) // Loop checking, 11.05.2015, T. Koi
243 {
244 if (counter != 0) buffer = result / counter;
245 while (counter < size) // Loop checking, 11.05.2015, T. Koi
246 {
247 counter++;
248 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
249 boosted.Lorentz(theNeutron, aThermalNuc);
250 G4double theEkin = boosted.GetKineticEnergy();
251 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
252 // velocity correction.
253 G4ThreeVector targetVelocity = 1. / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
254 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
255 result += aXsection;
256 }
257 size += size;
258 }
259 result /= counter;
260 /*
261 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
262 G4cout << " result " << result << " "
263 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
264 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
265 */
266 return result;
267}
268
273
278
280{
281 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for elastic "
282 "reaction of neutrons below 20MeV\n";
283}
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)
void DumpPhysicsTable(const G4ParticleDefinition &) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
void CrossSectionDescription(std::ostream &) const override
void SetVerboseLevel(G4int) override
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)
G4PhysicsTable * GetElasticCrossSections() const
void RegisterElasticCrossSections(G4PhysicsTable *val)
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
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