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