Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPElasticData.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//
37// P. Arce, June-2014 Conversion neutron_hp to particle_hp
38//
42#include "G4SystemOfUnits.hh"
43#include "G4Neutron.hh"
44#include "G4ElementTable.hh"
45#include "G4ParticleHPData.hh"
48#include "G4Pow.hh"
49
51:G4VCrossSectionDataSet("NeutronHPElasticXS")
52{
53 SetMinKinEnergy( 0*MeV );
54 SetMaxKinEnergy( 20*MeV );
55
56 theCrossSections = 0;
57 onFlightDB = true;
58 instanceOfWorker = false;
60 instanceOfWorker = true;
61 }
62 element_cache = NULL;
63 material_cache = NULL;
64 ke_cache = 0.0;
65 xs_cache = 0.0;
66// BuildPhysicsTable( *G4Neutron::Neutron() );
67}
68
70{
71 if ( theCrossSections != NULL && instanceOfWorker != true ) {
72 theCrossSections->clearAndDestroy();
73 delete theCrossSections;
74 theCrossSections = NULL;
75 }
76}
77
79 G4int /*Z*/ , G4int /*A*/ ,
80 const G4Element* /*elm*/ ,
81 const G4Material* /*mat*/ )
82{
83
84 G4double eKin = dp->GetKineticEnergy();
85 if ( eKin > GetMaxKinEnergy()
86 || eKin < GetMinKinEnergy()
87 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
88
89 return true;
90}
91
93 G4int /*Z*/ , G4int /*A*/ ,
94 const G4Isotope* /*iso*/ ,
95 const G4Element* element ,
96 const G4Material* material )
97{
98 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
99
100 ke_cache = dp->GetKineticEnergy();
101 element_cache = element;
102 material_cache = material;
103 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
104 xs_cache = xs;
105 return xs;
106}
107
108/*
109G4bool G4ParticleHPElasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
110{
111 G4bool result = true;
112 G4double eKin = aP->GetKineticEnergy();
113 if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
114 return result;
115}
116*/
117
119{
120
121 if(&aP!=G4Neutron::Neutron())
122 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
123
124//080428
125 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
126 {
127 onFlightDB = false;
128 #ifdef G4VERBOSE
130 G4cout << "Find a flag of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
131 G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of elastic scattering of neutrons (<20MeV)." << G4endl;
132 }
133 #endif
134 }
135
138 return;
139 }
140
141 size_t numberOfElements = G4Element::GetNumberOfElements();
142// TKDB
143 //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
144 if ( theCrossSections == NULL )
145 theCrossSections = new G4PhysicsTable( numberOfElements );
146 else
147 theCrossSections->clearAndDestroy();
148
149 // make a PhysicsVector for each element
150
151 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
152 for( size_t i=0; i<numberOfElements; ++i )
153 {
155 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
156 theCrossSections->push_back(physVec);
157 }
158
160}
161
163{
164 if(&aP!=G4Neutron::Neutron())
165 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
166
167 #ifdef G4VERBOSE
168 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
169
170//
171// Dump element based cross section
172// range 10e-5 eV to 20 MeV
173// 10 point per decade
174// in barn
175//
176
177 G4cout << G4endl;
178 G4cout << G4endl;
179 G4cout << "Elastic Cross Section of Neutron HP"<< G4endl;
180 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
181 G4cout << G4endl;
182 G4cout << "Name of Element" << G4endl;
183 G4cout << "Energy[eV] XS[barn]" << G4endl;
184 G4cout << G4endl;
185
186 size_t numberOfElements = G4Element::GetNumberOfElements();
187 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
188
189 for ( size_t i = 0 ; i < numberOfElements ; ++i )
190 {
191
192 G4cout << (*theElementTable)[i]->GetName() << G4endl;
193
194 G4int ie = 0;
195
196 for ( ie = 0 ; ie < 130 ; ie++ )
197 {
198 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
199 G4bool outOfRange = false;
200
201 if ( eKinetic < 20*MeV )
202 {
203 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
204 }
205
206 }
207
208 G4cout << G4endl;
209 }
210
211 //G4cout << "G4ParticleHPElasticData::DumpPhysicsTable still to be implemented"<<G4endl;
212 #endif
213}
214
215#include "G4Nucleus.hh"
216#include "G4NucleiProperties.hh"
217#include "G4Neutron.hh"
218#include "G4Electron.hh"
219
221GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
222{
223 G4double result = 0;
224 G4bool outOfRange;
225 G4int index = anE->GetIndex();
226
227 // prepare neutron
228 G4double eKinetic = aP->GetKineticEnergy();
229
230 if ( !onFlightDB )
231 {
232 //NEGLECT_DOPPLER_B.
233 G4double factor = 1.0;
234 if ( eKinetic < aT * k_Boltzmann )
235 {
236 // below 0.1 eV neutrons
237 // Have to do some, but now just igonre.
238 // Will take care after performance check.
239 // factor = factor * targetV;
240 }
241 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
242 }
243
244 G4ReactionProduct theNeutron( aP->GetDefinition() );
245 theNeutron.SetMomentum( aP->GetMomentum() );
246 theNeutron.SetKineticEnergy( eKinetic );
247
248 // prepare thermal nucleus
249 G4Nucleus aNuc;
250 G4double eps = 0.0001;
251 G4double theA = anE->GetN();
252 G4double theZ = anE->GetZ();
253 G4double eleMass;
254
255
256 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
258
259 G4ReactionProduct boosted;
260 G4double aXsection;
261
262 // MC integration loop
263 G4int counter = 0;
264 G4double buffer = 0;
265 G4int size = G4int(std::max(10., aT/60*kelvin));
266 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
267 G4double neutronVMag = neutronVelocity.mag();
268
269 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
270 {
271 if(counter) buffer = result/counter;
272 while (counter<size) // Loop checking, 11.05.2015, T. Koi
273 {
274 counter ++;
275 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
276 boosted.Lorentz(theNeutron, aThermalNuc);
277 G4double theEkin = boosted.GetKineticEnergy();
278 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
279 // velocity correction.
280 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
281 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
282 result += aXsection;
283 }
284 size += size;
285 }
286 result /= counter;
287/*
288 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
289 G4cout << " result " << result << " "
290 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
291 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
292*/
293 return result;
294}
295
297GetVerboseLevel() const
298{
300}
301
303SetVerboseLevel( G4int newValue )
304{
306}
308{
309 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for elastic reaction of neutrons below 20MeV\n" ;
310}
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)
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 &)
virtual void CrossSectionDescription(std::ostream &) const
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4PhysicsTable * GetElasticCrossSections()
void RegisterElasticCrossSections(G4PhysicsTable *val)
static G4ParticleHPManager * GetInstance()
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