Geant4-11
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//
40#include "G4SystemOfUnits.hh"
41#include "G4Neutron.hh"
42#include "G4ElementTable.hh"
43#include "G4ParticleHPData.hh"
46#include "G4Pow.hh"
47
49:G4VCrossSectionDataSet("NeutronHPFissionXS")
50{
51 SetMinKinEnergy( 0*MeV );
52 SetMaxKinEnergy( 20*MeV );
53
55 instanceOfWorker = false;
57 instanceOfWorker = true;
58 }
59 element_cache = NULL;
60 material_cache = NULL;
61 ke_cache = 0.0;
62 xs_cache = 0.0;
63 //BuildPhysicsTable(*G4Neutron::Neutron());
64}
65
67{
68 if ( theCrossSections != NULL && instanceOfWorker != true ) {
70 delete theCrossSections;
71 theCrossSections = NULL;
72 }
73}
74
76 G4int /*Z*/ , G4int /*A*/ ,
77 const G4Element* /*elm*/ ,
78 const G4Material* /*mat*/ )
79{
80 G4double eKin = dp->GetKineticEnergy();
81 if ( eKin > GetMaxKinEnergy()
82 || eKin < GetMinKinEnergy()
83 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
84
85 return true;
86}
87
89 G4int /*Z*/ , G4int /*A*/ ,
90 const G4Isotope* /*iso*/ ,
91 const G4Element* element ,
92 const G4Material* material )
93{
94 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
95
97 element_cache = element;
99 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
100 xs_cache = xs;
101 return xs;
102}
103
104/*
105G4bool G4ParticleHPFissionData::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__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
118
121 return;
122 }
123
124 size_t numberOfElements = G4Element::GetNumberOfElements();
125 //theCrossSections = new G4PhysicsTable( numberOfElements );
126 // TKDB
127 //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
128 if ( theCrossSections == NULL )
129 theCrossSections = new G4PhysicsTable( numberOfElements );
130 else
132
133 // make a PhysicsVector for each element
134
135 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
136 for( size_t i=0; i<numberOfElements; ++i )
137 {
139 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
140 theCrossSections->push_back(physVec);
141 }
142
144}
145
147{
148 if(&aP!=G4Neutron::Neutron())
149 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
150
151 #ifdef G4VERBOSE
152 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
153
154//
155// Dump element based cross section
156// range 10e-5 eV to 20 MeV
157// 10 point per decade
158// in barn
159//
160
161 G4cout << G4endl;
162 G4cout << G4endl;
163 G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
164 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
165 G4cout << G4endl;
166 G4cout << "Name of Element" << G4endl;
167 G4cout << "Energy[eV] XS[barn]" << G4endl;
168 G4cout << G4endl;
169
170 size_t numberOfElements = G4Element::GetNumberOfElements();
171 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
172
173 for ( size_t i = 0 ; i < numberOfElements ; ++i )
174 {
175
176 G4cout << (*theElementTable)[i]->GetName() << G4endl;
177
178 if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
179 {
180 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
181 G4cout << G4endl;
182 continue;
183 }
184
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 }
198
199 G4cout << G4endl;
200 }
201
202 //G4cout << "G4ParticleHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
203 #endif
204}
205
206#include "G4NucleiProperties.hh"
207
209GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
210{
211 G4double result = 0;
212 if(anE->GetZ()<88) return result;
213 G4bool outOfRange;
214 G4int index = anE->GetIndex();
215
216// 100729 TK add safety
217if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
218
219 // prepare neutron
220 G4double eKinetic = aP->GetKineticEnergy();
221 G4ReactionProduct theNeutronRP( aP->GetDefinition() );
222 theNeutronRP.SetMomentum( aP->GetMomentum() );
223 theNeutronRP.SetKineticEnergy( eKinetic );
224
225 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
226 {
227 //NEGLECT_DOPPLER
228 G4double factor = 1.0;
229 if ( eKinetic < aT * k_Boltzmann ) {
230 // below 0.1 eV neutrons
231 // Have to do some, but now just igonre.
232 // Will take care after performance check.
233 // factor = factor * targetV;
234 }
235 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
236 }
237
238 // prepare thermal nucleus
239 G4Nucleus aNuc;
240 G4double eps = 0.0001;
241 G4double theA = anE->GetN();
242 G4double theZ = anE->GetZ();
243 G4double eleMass;
244 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
246
247 G4ReactionProduct boosted;
248 G4double aXsection;
249
250 // MC integration loop
251 G4int counter = 0;
252 G4double buffer = 0;
253 G4int size = G4int(std::max(10., aT/60*kelvin));
254 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutronRP.GetMomentum();
255 G4double neutronVMag = neutronVelocity.mag();
256
257 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
258 {
259 if(counter) buffer = result/counter;
260 while (counter<size) // Loop checking, 11.05.2015, T. Koi
261 {
262 counter ++;
263 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
264 boosted.Lorentz(theNeutronRP, aThermalNuc);
265 G4double theEkin = boosted.GetKineticEnergy();
266 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
267 // velocity correction.
268 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
269 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
270 result += aXsection;
271 }
272 size += size;
273 }
274 result /= counter;
275 return result;
276}
277
279{
281}
283{
285}
287{
288 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n" ;
289}
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)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
void DumpPhysicsTable(const G4ParticleDefinition &)
virtual void CrossSectionDescription(std::ostream &) const
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4PhysicsTable * GetFissionCrossSections()
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)
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