Geant4 10.7.0
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 "G4Pow.hh"
47
49:G4VCrossSectionDataSet("NeutronHPFissionXS")
50{
51 SetMinKinEnergy( 0*MeV );
52 SetMaxKinEnergy( 20*MeV );
53
54 theCrossSections = 0;
55 onFlightDB = true;
56 instanceOfWorker = false;
58 instanceOfWorker = true;
59 }
60 element_cache = NULL;
61 material_cache = NULL;
62 ke_cache = 0.0;
63 xs_cache = 0.0;
64 //BuildPhysicsTable(*G4Neutron::Neutron());
65}
66
68{
69 if ( theCrossSections != NULL && instanceOfWorker != true ) {
70 theCrossSections->clearAndDestroy();
71 delete theCrossSections;
72 theCrossSections = NULL;
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
105/*
106G4bool G4ParticleHPFissionData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
107{
108 G4bool result = true;
109 G4double eKin = aP->GetKineticEnergy();
110 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
111 return result;
112}
113*/
114
116{
117
118 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() ) {
119 onFlightDB = false;
120 #ifdef G4VERBOSE
122 G4cout << "Find a flag of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
123 G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of fission reaction of neutrons (<20MeV)." << G4endl;
124 }
125 #endif
126 }
127
128 if(&aP!=G4Neutron::Neutron())
129 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
130
133 return;
134 }
135
136 size_t numberOfElements = G4Element::GetNumberOfElements();
137 //theCrossSections = new G4PhysicsTable( numberOfElements );
138 // TKDB
139 //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
140 if ( theCrossSections == NULL )
141 theCrossSections = new G4PhysicsTable( numberOfElements );
142 else
143 theCrossSections->clearAndDestroy();
144
145 // make a PhysicsVector for each element
146
147 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
148 for( size_t i=0; i<numberOfElements; ++i )
149 {
151 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
152 theCrossSections->push_back(physVec);
153 }
154
156}
157
159{
160 if(&aP!=G4Neutron::Neutron())
161 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
162
163 #ifdef G4VERBOSE
164 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
165
166//
167// Dump element based cross section
168// range 10e-5 eV to 20 MeV
169// 10 point per decade
170// in barn
171//
172
173 G4cout << G4endl;
174 G4cout << G4endl;
175 G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
176 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
177 G4cout << G4endl;
178 G4cout << "Name of Element" << G4endl;
179 G4cout << "Energy[eV] XS[barn]" << G4endl;
180 G4cout << G4endl;
181
182 size_t numberOfElements = G4Element::GetNumberOfElements();
183 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
184
185 for ( size_t i = 0 ; i < numberOfElements ; ++i )
186 {
187
188 G4cout << (*theElementTable)[i]->GetName() << G4endl;
189
190 if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
191 {
192 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
193 G4cout << G4endl;
194 continue;
195 }
196
197 G4int ie = 0;
198
199 for ( ie = 0 ; ie < 130 ; ie++ )
200 {
201 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
202 G4bool outOfRange = false;
203
204 if ( eKinetic < 20*MeV )
205 {
206 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
207 }
208
209 }
210
211 G4cout << G4endl;
212 }
213
214 //G4cout << "G4ParticleHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
215 #endif
216}
217
218#include "G4NucleiProperties.hh"
219
221GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
222{
223 G4double result = 0;
224 if(anE->GetZ()<88) return result;
225 G4bool outOfRange;
226 G4int index = anE->GetIndex();
227
228// 100729 TK add safety
229if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
230
231 // prepare neutron
232 G4double eKinetic = aP->GetKineticEnergy();
233 G4ReactionProduct theNeutronRP( aP->GetDefinition() );
234 theNeutronRP.SetMomentum( aP->GetMomentum() );
235 theNeutronRP.SetKineticEnergy( eKinetic );
236
237 if ( !onFlightDB ) {
238 //NEGLECT_DOPPLER
239 G4double factor = 1.0;
240 if ( eKinetic < aT * k_Boltzmann ) {
241 // below 0.1 eV neutrons
242 // Have to do some, but now just igonre.
243 // Will take care after performance check.
244 // factor = factor * targetV;
245 }
246 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
247 }
248
249 // prepare thermal nucleus
250 G4Nucleus aNuc;
251 G4double eps = 0.0001;
252 G4double theA = anE->GetN();
253 G4double theZ = anE->GetZ();
254 G4double eleMass;
255 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
257
258 G4ReactionProduct boosted;
259 G4double aXsection;
260
261 // MC integration loop
262 G4int counter = 0;
263 G4double buffer = 0;
264 G4int size = G4int(std::max(10., aT/60*kelvin));
265 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutronRP.GetMomentum();
266 G4double neutronVMag = neutronVelocity.mag();
267
268 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
269 {
270 if(counter) buffer = result/counter;
271 while (counter<size) // Loop checking, 11.05.2015, T. Koi
272 {
273 counter ++;
274 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
275 boosted.Lorentz(theNeutronRP, aThermalNuc);
276 G4double theEkin = boosted.GetKineticEnergy();
277 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
278 // velocity correction.
279 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
280 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
281 result += aXsection;
282 }
283 size += size;
284 }
285 result /= counter;
286 return result;
287}
288
290{
292}
294{
296}
298{
299 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n" ;
300}
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:397
G4double GetZ() const
Definition: G4Element.hh:130
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
size_t GetIndex() const
Definition: G4Element.hh:181
G4double GetN() const
Definition: G4Element.hh:134
static G4HadronicParameters * Instance()
G4double GetTemperature() const
Definition: G4Material.hh:180
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:143
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
#define buffer
Definition: xmlparse.cc:628