Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPFissionData.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. Koi
35
37#include "G4SystemOfUnits.hh"
38#include "G4Neutron.hh"
39#include "G4ElementTable.hh"
40#include "G4NeutronHPData.hh"
41
43:G4VCrossSectionDataSet("NeutronHPFissionXS")
44{
45 SetMinKinEnergy( 0*MeV );
46 SetMaxKinEnergy( 20*MeV );
47
48 ke_cache = 0.0;
49 xs_cache = 0.0;
50 element_cache = NULL;
51 material_cache = NULL;
52
53 theCrossSections = 0;
55}
56
58{
59 if ( theCrossSections != NULL ) theCrossSections->clearAndDestroy();
60 delete theCrossSections;
61}
62
64 G4int /*Z*/ , G4int /*A*/ ,
65 const G4Element* /*elm*/ ,
66 const G4Material* /*mat*/ )
67{
68 G4double eKin = dp->GetKineticEnergy();
69 if ( eKin > GetMaxKinEnergy()
70 || eKin < GetMinKinEnergy()
71 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
72
73 return true;
74}
75
77 G4int /*Z*/ , G4int /*A*/ ,
78 const G4Isotope* /*iso*/ ,
79 const G4Element* element ,
80 const G4Material* material )
81{
82 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
83
84 ke_cache = dp->GetKineticEnergy();
85 element_cache = element;
86 material_cache = material;
87 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
88 xs_cache = xs;
89 return xs;
90}
91
92/*
93G4bool G4NeutronHPFissionData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
94{
95 G4bool result = true;
96 G4double eKin = aP->GetKineticEnergy();
97 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
98 return result;
99}
100*/
101
103{
104 if(&aP!=G4Neutron::Neutron())
105 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
106 size_t numberOfElements = G4Element::GetNumberOfElements();
107 //theCrossSections = new G4PhysicsTable( numberOfElements );
108 // TKDB
109 //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
110 if ( theCrossSections == NULL )
111 theCrossSections = new G4PhysicsTable( numberOfElements );
112 else
113 theCrossSections->clearAndDestroy();
114
115 // make a PhysicsVector for each element
116
117 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
118 for( size_t i=0; i<numberOfElements; ++i )
119 {
121 Instance()->MakePhysicsVector((*theElementTable)[i], this);
122 theCrossSections->push_back(physVec);
123 }
124}
125
127{
128 if(&aP!=G4Neutron::Neutron())
129 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
130
131//
132// Dump element based cross section
133// range 10e-5 eV to 20 MeV
134// 10 point per decade
135// in barn
136//
137
138 G4cout << G4endl;
139 G4cout << G4endl;
140 G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
141 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
142 G4cout << G4endl;
143 G4cout << "Name of Element" << G4endl;
144 G4cout << "Energy[eV] XS[barn]" << G4endl;
145 G4cout << G4endl;
146
147 size_t numberOfElements = G4Element::GetNumberOfElements();
148 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
149
150 for ( size_t i = 0 ; i < numberOfElements ; ++i )
151 {
152
153 G4cout << (*theElementTable)[i]->GetName() << G4endl;
154
155 if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
156 {
157 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
158 G4cout << G4endl;
159 continue;
160 }
161
162 G4int ie = 0;
163
164 for ( ie = 0 ; ie < 130 ; ie++ )
165 {
166 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
167 G4bool outOfRange = false;
168
169 if ( eKinetic < 20*MeV )
170 {
171 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
172 }
173
174 }
175
176 G4cout << G4endl;
177 }
178
179 //G4cout << "G4NeutronHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
180}
181
182#include "G4NucleiProperties.hh"
183
185GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
186{
187 G4double result = 0;
188 if(anE->GetZ()<90) return result;
189 G4bool outOfRange;
190 G4int index = anE->GetIndex();
191
192// 100729 TK add safety
193if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
194
195 // prepare neutron
196 G4double eKinetic = aP->GetKineticEnergy();
197 G4ReactionProduct theNeutron( aP->GetDefinition() );
198 theNeutron.SetMomentum( aP->GetMomentum() );
199 theNeutron.SetKineticEnergy( eKinetic );
200
201 // prepare thermal nucleus
202 G4Nucleus aNuc;
203 G4double eps = 0.0001;
204 G4double theA = anE->GetN();
205 G4double theZ = anE->GetZ();
206 G4double eleMass;
207 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
209
210 G4ReactionProduct boosted;
211 G4double aXsection;
212
213 // MC integration loop
214 G4int counter = 0;
215 G4double buffer = 0;
216 G4int size = G4int(std::max(10., aT/60*kelvin));
217 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
218 G4double neutronVMag = neutronVelocity.mag();
219
220 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
221 {
222 if(counter) buffer = result/counter;
223 while (counter<size)
224 {
225 counter ++;
226 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
227 boosted.Lorentz(theNeutron, aThermalNuc);
228 G4double theEkin = boosted.GetKineticEnergy();
229 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
230 // velocity correction.
231 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
232 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
233 result += aXsection;
234 }
235 size += size;
236 }
237 result /= counter;
238 return result;
239}
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()
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 &)
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