Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPCaptureData.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("NeutronHPCaptureXS")
46{
47 SetMinKinEnergy( 0*MeV );
48 SetMaxKinEnergy( 20*MeV );
49
50 ke_cache = 0.0;
51 xs_cache = 0.0;
52 element_cache = NULL;
53 material_cache = NULL;
54
55 theCrossSections = 0;
56 onFlightDB = true;
57
59}
60
62{
63 if ( theCrossSections != 0 ) theCrossSections->clearAndDestroy();
64
65 delete theCrossSections;
66}
67
69 G4int /*Z*/ , G4int /*A*/ ,
70 const G4Element* /*elm*/ ,
71 const G4Material* /*mat*/ )
72{
73 G4double eKin = dp->GetKineticEnergy();
74 if ( eKin > GetMaxKinEnergy()
75 || eKin < GetMinKinEnergy()
76 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
77
78 return true;
79}
80
82 G4int /*Z*/ , G4int /*A*/ ,
83 const G4Isotope* /*iso*/ ,
84 const G4Element* element ,
85 const G4Material* material )
86{
87 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
88
89 ke_cache = dp->GetKineticEnergy();
90 element_cache = element;
91 material_cache = material;
92 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
93 xs_cache = xs;
94 return xs;
95}
96
97/*
98G4bool G4NeutronHPCaptureData::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 capture reaction of neutrons (<20MeV)." << G4endl;
117 onFlightDB = false;
118 }
119
120 size_t numberOfElements = G4Element::GetNumberOfElements();
121 // G4cout << "CALLED G4NeutronHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl;
122 // TKDB
123 //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
124 if ( theCrossSections == NULL )
125 theCrossSections = new G4PhysicsTable( numberOfElements );
126 else
127 theCrossSections->clearAndDestroy();
128
129 // make a PhysicsVector for each element
130
131 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
132 for( size_t i=0; i<numberOfElements; ++i )
133 {
134 if(getenv("CaptureDataIndexDebug"))
135 {
136 G4int index_debug = ((*theElementTable)[i])->GetIndex();
137 G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl;
138 }
140 Instance()->MakePhysicsVector((*theElementTable)[i], this);
141 theCrossSections->push_back(physVec);
142 }
143}
144
146{
147 if(&aP!=G4Neutron::Neutron())
148 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
149
150//
151// Dump element based cross section
152// range 10e-5 eV to 20 MeV
153// 10 point per decade
154// in barn
155//
156
157 G4cout << G4endl;
158 G4cout << G4endl;
159 G4cout << "Capture Cross Section of Neutron HP"<< G4endl;
160 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
161 G4cout << G4endl;
162 G4cout << "Name of Element" << G4endl;
163 G4cout << "Energy[eV] XS[barn]" << G4endl;
164 G4cout << G4endl;
165
166 size_t numberOfElements = G4Element::GetNumberOfElements();
167 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
168
169 for ( size_t i = 0 ; i < numberOfElements ; ++i )
170 {
171
172 G4cout << (*theElementTable)[i]->GetName() << G4endl;
173
174 G4int ie = 0;
175
176 for ( ie = 0 ; ie < 130 ; ie++ )
177 {
178 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
179 G4bool outOfRange = false;
180
181 if ( eKinetic < 20*MeV )
182 {
183 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
184 }
185
186 }
187
188 G4cout << G4endl;
189 }
190
191
192// G4cout << "G4NeutronHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
193}
194
195#include "G4NucleiProperties.hh"
196
198GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
199{
200 G4double result = 0;
201 G4bool outOfRange;
202 G4int index = anE->GetIndex();
203
204 // prepare neutron
205 G4double eKinetic = aP->GetKineticEnergy();
206
207//if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
208//080428
209 if ( !onFlightDB )
210 {
211 G4double factor = 1.0;
212 if ( eKinetic < aT * k_Boltzmann )
213 {
214 // below 0.1 eV neutrons
215 // Have to do some, but now just igonre.
216 // Will take care after performance check.
217 // factor = factor * targetV;
218 }
219 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
220 }
221
222 G4ReactionProduct theNeutron( aP->GetDefinition() );
223 theNeutron.SetMomentum( aP->GetMomentum() );
224 theNeutron.SetKineticEnergy( eKinetic );
225
226 // prepare thermal nucleus
227 G4Nucleus aNuc;
228 G4double eps = 0.0001;
229 G4double theA = anE->GetN();
230 G4double theZ = anE->GetZ();
231 G4double eleMass;
232 eleMass = G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) / G4Neutron::Neutron()->GetPDGMass();
233
234 G4ReactionProduct boosted;
235 G4double aXsection;
236
237 // MC integration loop
238 G4int counter = 0;
239 G4double buffer = 0;
240 G4int size = G4int(std::max(10., aT/60*kelvin));
241 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
242 G4double neutronVMag = neutronVelocity.mag();
243
244 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
245 {
246 if(counter) buffer = result/counter;
247 while (counter<size)
248 {
249 counter ++;
250 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
251 boosted.Lorentz(theNeutron, aThermalNuc);
252 G4double theEkin = boosted.GetKineticEnergy();
253 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
254 // velocity correction, or luminosity factor...
255 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
256 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
257 result += aXsection;
258 }
259 size += size;
260 }
261 result /= counter;
262/*
263 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
264 G4cout << " result " << result << " "
265 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
266 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
267*/
268 return result;
269}
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
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4NeutronHPFissionData *theP)
static G4NeutronHPData * Instance()
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