Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 "G4NucleiProperties.hh"
47#include "G4Pow.hh"
48
50:G4VCrossSectionDataSet("NeutronHPFissionXS")
51{
52 SetMinKinEnergy( 0*MeV );
53 SetMaxKinEnergy( 20*MeV );
54
55 theCrossSections = 0;
56 instanceOfWorker = false;
58 instanceOfWorker = true;
59 }
60 element_cache = nullptr;
61 material_cache = nullptr;
62 ke_cache = 0.0;
63 xs_cache = 0.0;
64 //BuildPhysicsTable(*G4Neutron::Neutron());
65}
66
68{
69 if ( theCrossSections != nullptr && instanceOfWorker != true ) {
70 theCrossSections->clearAndDestroy();
71 delete theCrossSections;
72 theCrossSections = nullptr;
73 }
74}
75
77 G4int /*Z*/ , G4int /*A*/ ,
78 const G4Element* /*elm*/ ,
79 const G4Material* /*mat*/ )
80{
81 G4double eKin = dp->GetKineticEnergy();
82 if ( eKin > GetMaxKinEnergy()
83 || eKin < GetMinKinEnergy()
84 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
85
86 return true;
87}
88
90 G4int /*Z*/ , G4int /*A*/ ,
91 const G4Isotope* /*iso*/ ,
92 const G4Element* element ,
93 const G4Material* material )
94{
95 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
96
97 ke_cache = dp->GetKineticEnergy();
98 element_cache = element;
99 material_cache = material;
100 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
101 xs_cache = xs;
102 return xs;
103}
104
106{
107 if(&aP!=G4Neutron::Neutron())
108 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
109
112 return;
113 }
114
115 std::size_t numberOfElements = G4Element::GetNumberOfElements();
116 if ( theCrossSections == nullptr )
117 theCrossSections = new G4PhysicsTable( numberOfElements );
118 else
119 theCrossSections->clearAndDestroy();
120
121 // make a PhysicsVector for each element
122
123 static G4ThreadLocal G4ElementTable *theElementTable = nullptr ;
124 if (!theElementTable) theElementTable= G4Element::GetElementTable();
125 for( std::size_t i=0; i<numberOfElements; ++i )
126 {
128 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
129 theCrossSections->push_back(physVec);
130 }
131
133}
134
136{
137 if(&aP!=G4Neutron::Neutron())
138 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
139
140#ifdef G4VERBOSE
141 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
142
143 //
144 // Dump element based cross section
145 // range 10e-5 eV to 20 MeV
146 // 10 point per decade
147 // in barn
148 //
149 G4cout << G4endl;
150 G4cout << G4endl;
151 G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
152 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
153 G4cout << G4endl;
154 G4cout << "Name of Element" << G4endl;
155 G4cout << "Energy[eV] XS[barn]" << G4endl;
156 G4cout << G4endl;
157
158 std::size_t numberOfElements = G4Element::GetNumberOfElements();
159 static G4ThreadLocal G4ElementTable *theElementTable = nullptr ;
160 if (!theElementTable) theElementTable= G4Element::GetElementTable();
161
162 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
163 {
164 G4cout << (*theElementTable)[i]->GetName() << G4endl;
165
166 if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
167 {
168 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
169 G4cout << G4endl;
170 continue;
171 }
172
173 for ( G4int ie = 0 ; ie < 130 ; ++ie )
174 {
175 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
176 G4bool outOfRange = false;
177
178 if ( eKinetic < 20*MeV )
179 {
180 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
181 }
182 }
183
184 G4cout << G4endl;
185 }
186#endif
187}
188
190GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
191{
192 G4double result = 0;
193 if(anE->GetZ()<88) return result;
194 G4bool outOfRange;
195 G4int index = (G4int)anE->GetIndex();
196
197 if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 )
198 return result;
199
200 // prepare neutron
201 G4double eKinetic = aP->GetKineticEnergy();
202 G4ReactionProduct theNeutronRP( aP->GetDefinition() );
203 theNeutronRP.SetMomentum( aP->GetMomentum() );
204 theNeutronRP.SetKineticEnergy( eKinetic );
205
206 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
207 {
208 //NEGLECT_DOPPLER
209 G4double factor = 1.0;
210 if ( eKinetic < aT * k_Boltzmann ) {
211 // below 0.1 eV neutrons
212 // Have to do some, but now just igonre.
213 // Will take care after performance check.
214 // factor = factor * targetV;
215 }
216 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
217 }
218
219 // prepare thermal nucleus
220 G4Nucleus aNuc;
221 G4double eps = 0.0001;
222 G4double theA = anE->GetN();
223 G4double theZ = anE->GetZ();
224 G4double eleMass;
225 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
227
228 G4ReactionProduct boosted;
229 G4double aXsection;
230
231 // MC integration loop
232 G4int counter = 0;
233 G4double buffer = 0;
234 G4int size = G4int(std::max(10., aT/60*kelvin));
235 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutronRP.GetMomentum();
236 G4double neutronVMag = neutronVelocity.mag();
237
238 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
239 {
240 if(counter) buffer = result/counter;
241 while (counter<size) // Loop checking, 11.05.2015, T. Koi
242 {
243 counter ++;
244 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
245 boosted.Lorentz(theNeutronRP, aThermalNuc);
246 G4double theEkin = boosted.GetKineticEnergy();
247 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
248 // velocity correction.
249 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
250 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
251 result += aXsection;
252 }
253 size += size;
254 }
255 result /= counter;
256 return result;
257}
258
260{
262}
263
265{
267}
268
270{
271 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF)\n"
272 << "for induced fission reaction of neutrons below 20MeV\n" ;
273}
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)
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)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
#define G4ThreadLocal
Definition: tls.hh:77