Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4NeutronHPInelasticData.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//
39#include "G4SystemOfUnits.hh"
40#include "G4Neutron.hh"
41#include "G4ElementTable.hh"
42#include "G4NeutronHPData.hh"
43
45:G4VCrossSectionDataSet("NeutronHPInelasticXS")
46{
47
48 SetMinKinEnergy( 0*MeV );
49 SetMaxKinEnergy( 20*MeV );
50
51 ke_cache = 0.0;
52 xs_cache = 0.0;
53 element_cache = NULL;
54 material_cache = NULL;
55
56 onFlightDB = true;
57 theCrossSections = 0;
59}
60
62{
63 if ( theCrossSections != 0 ) theCrossSections->clearAndDestroy();
64 delete theCrossSections;
65}
66
68 G4int /*Z*/ , G4int /*A*/ ,
69 const G4Element* /*elm*/ ,
70 const G4Material* /*mat*/ )
71{
72 G4double eKin = dp->GetKineticEnergy();
73 if ( eKin > GetMaxKinEnergy()
74 || eKin < GetMinKinEnergy()
75 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
76
77 return true;
78}
79
81 G4int /*Z*/ , G4int /*A*/ ,
82 const G4Isotope* /*iso*/ ,
83 const G4Element* element ,
84 const G4Material* material )
85{
86 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
87
88 ke_cache = dp->GetKineticEnergy();
89 element_cache = element;
90 material_cache = material;
91 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
92 xs_cache = xs;
93 return xs;
94 //return GetCrossSection( dp , element , material->GetTemperature() );
95}
96
97/*
98G4bool G4NeutronHPInelasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
99{
100 G4bool result = true;
101 G4double eKin = aP->GetKineticEnergy();
102 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
103 return result;
104}
105*/
106
108{
109 if(&aP!=G4Neutron::Neutron())
110 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
111
112//080428
113 if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
114 {
115 G4cout << "Find environment variable of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
116 G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of inelastic scattering of neutrons (<20MeV)." << G4endl;
117 onFlightDB = false;
118 }
119
120 size_t numberOfElements = G4Element::GetNumberOfElements();
121// theCrossSections = new G4PhysicsTable( numberOfElements );
122// TKDB
123 //if ( theCrossSections == 0 )
124 //{ theCrossSections = new G4PhysicsTable( numberOfElements ); }
125 if ( theCrossSections == NULL )
126 theCrossSections = new G4PhysicsTable( numberOfElements );
127 else
128 theCrossSections->clearAndDestroy();
129
130 // make a PhysicsVector for each element
131
132 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
133 for( size_t i=0; i<numberOfElements; ++i )
134 {
136 Instance()->MakePhysicsVector((*theElementTable)[i], this);
137 theCrossSections->push_back(physVec);
138 }
139}
140
142{
143 if(&aP!=G4Neutron::Neutron())
144 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
145
146//
147// Dump element based cross section
148// range 10e-5 eV to 20 MeV
149// 10 point per decade
150// in barn
151//
152
153 G4cout << G4endl;
154 G4cout << G4endl;
155 G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
156 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
157 G4cout << G4endl;
158 G4cout << "Name of Element" << G4endl;
159 G4cout << "Energy[eV] XS[barn]" << G4endl;
160 G4cout << G4endl;
161
162 size_t numberOfElements = G4Element::GetNumberOfElements();
163 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
164
165 for ( size_t i = 0 ; i < numberOfElements ; ++i )
166 {
167
168 G4cout << (*theElementTable)[i]->GetName() << G4endl;
169
170 G4int ie = 0;
171
172 for ( ie = 0 ; ie < 130 ; ie++ )
173 {
174 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
175 G4bool outOfRange = false;
176
177 if ( eKinetic < 20*MeV )
178 {
179 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
180 }
181
182 }
183
184 G4cout << G4endl;
185 }
186
187 //G4cout << "G4NeutronHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
188}
189
190#include "G4NucleiProperties.hh"
191
193GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
194{
195 G4double result = 0;
196 G4bool outOfRange;
197 G4int index = anE->GetIndex();
198
199 // prepare neutron
200 G4double eKinetic = aP->GetKineticEnergy();
201
202 // T. K.
203//if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
204//080428
205 if ( !onFlightDB )
206 {
207 G4double factor = 1.0;
208 if ( eKinetic < aT * k_Boltzmann )
209 {
210 // below 0.1 eV neutrons
211 // Have to do some, but now just igonre.
212 // Will take care after performance check.
213 // factor = factor * targetV;
214 }
215 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
216 }
217
218 G4ReactionProduct theNeutron( aP->GetDefinition() );
219 theNeutron.SetMomentum( aP->GetMomentum() );
220 theNeutron.SetKineticEnergy( eKinetic );
221
222 // prepare thermal nucleus
223 G4Nucleus aNuc;
224 G4double eps = 0.0001;
225 G4double theA = anE->GetN();
226 G4double theZ = anE->GetZ();
227 G4double eleMass;
228 eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
230
231 G4ReactionProduct boosted;
232 G4double aXsection;
233
234 // MC integration loop
235 G4int counter = 0;
236 G4int failCount = 0;
237 G4double buffer = 0;
238 G4int size = G4int(std::max(10., aT/60*kelvin));
239 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
240 G4double neutronVMag = neutronVelocity.mag();
241
242 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
243 {
244 if(counter) buffer = result/counter;
245 while (counter<size)
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 if(aXsection <0)
253 {
254 if(failCount<1000)
255 {
256 failCount++;
257 counter--;
258 continue;
259 }
260 else
261 {
262 aXsection = 0;
263 }
264 }
265 // velocity correction.
266 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
267 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
268 result += aXsection;
269 }
270 size += size;
271 }
272 result /= counter;
273/*
274 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
275 G4cout << " result " << result << " "
276 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
277 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
278*/
279 return result;
280}
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
size_t GetIndex() const
Definition: G4Element.hh:182
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
G4double GetN() const
Definition: G4Element.hh:134
G4double GetTemperature() const
Definition: G4Material.hh:181
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4NeutronHPFissionData *theP)
static G4NeutronHPData * Instance()
void DumpPhysicsTable(const G4ParticleDefinition &)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
void push_back(G4PhysicsVector *)
void clearAndDestroy()
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)
#define buffer
Definition: xmlparse.cc:611