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