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