Geant4 11.1.1
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 "G4NucleiProperties.hh"
50#include "G4Pow.hh"
51
53:G4VCrossSectionDataSet("NeutronHPCaptureXS")
54{
55 SetMinKinEnergy( 0*MeV );
56 SetMaxKinEnergy( 20*MeV );
57
58 theCrossSections = 0;
59
60 instanceOfWorker = false;
62 instanceOfWorker = true;
63 }
64
65 element_cache = nullptr;
66 material_cache = nullptr;
67 ke_cache = 0.0;
68 xs_cache = 0.0;
69
70 //BuildPhysicsTable(*G4Neutron::Neutron());
71}
72
74{
75 if ( theCrossSections != nullptr && instanceOfWorker != true ) {
76 theCrossSections->clearAndDestroy();
77 delete theCrossSections;
78 theCrossSections = nullptr;
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
128 return;
129 }
130
131 std::size_t numberOfElements = G4Element::GetNumberOfElements();
132 // G4cout << "CALLED G4ParticleHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl;
133 // TKDB
134 //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
135 if ( theCrossSections == nullptr )
136 theCrossSections = new G4PhysicsTable( numberOfElements );
137 else
138 theCrossSections->clearAndDestroy();
139
140 // make a PhysicsVector for each element
141
142 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
143 for( std::size_t i=0; i<numberOfElements; ++i )
144 {
145 #ifdef G4VERBOSE
146 if(std::getenv("CaptureDataIndexDebug"))
147 {
148 std::size_t index_debug = ((*theElementTable)[i])->GetIndex();
149 if ( G4HadronicParameters::Instance()->GetVerboseLevel() > 0 ) G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl;
150 }
151 #endif
153 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
154 theCrossSections->push_back(physVec);
155 }
156
158}
159
161{
162 if(&aP!=G4Neutron::Neutron())
163 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
164
165 #ifdef G4VERBOSE
166 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
167
168//
169// Dump element based cross section
170// range 10e-5 eV to 20 MeV
171// 10 point per decade
172// in barn
173//
174
175 G4cout << G4endl;
176 G4cout << G4endl;
177 G4cout << "Capture Cross Section of Neutron HP"<< G4endl;
178 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
179 G4cout << G4endl;
180 G4cout << "Name of Element" << G4endl;
181 G4cout << "Energy[eV] XS[barn]" << G4endl;
182 G4cout << G4endl;
183
184 std::size_t numberOfElements = G4Element::GetNumberOfElements();
185 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
186 if (!theElementTable) theElementTable= G4Element::GetElementTable();
187
188 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
189 {
190
191 G4cout << (*theElementTable)[i]->GetName() << G4endl;
192
193 G4int ie = 0;
194
195 for ( ie = 0 ; ie < 130 ; ++ie )
196 {
197 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
198 G4bool outOfRange = false;
199
200 if ( eKinetic < 20*MeV )
201 {
202 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
203 }
204
205 }
206
207 G4cout << G4endl;
208 }
209
210 //G4cout << "G4ParticleHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
211 #endif
212}
213
215GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
216{
217 G4double result = 0;
218 G4bool outOfRange;
219 G4int index = (G4int)anE->GetIndex();
220
221 // prepare neutron
222 G4double eKinetic = aP->GetKineticEnergy();
223
225 {
226 //NEGLECT_DOPPLER
227 G4double factor = 1.0;
228 if ( eKinetic < aT * k_Boltzmann )
229 {
230 // below 0.1 eV neutrons
231 // Have to do some, but now just igonre.
232 // Will take care after performance check.
233 // factor = factor * targetV;
234 }
235 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
236 }
237
238 G4ReactionProduct theNeutron( aP->GetDefinition() );
239 theNeutron.SetMomentum( aP->GetMomentum() );
240 theNeutron.SetKineticEnergy( eKinetic );
241
242 // prepare thermal nucleus
243 G4Nucleus aNuc;
244 G4double eps = 0.0001;
245 G4double theA = anE->GetN();
246 G4double theZ = anE->GetZ();
247 G4double eleMass;
248 eleMass = G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) / G4Neutron::Neutron()->GetPDGMass();
249
250 G4ReactionProduct boosted;
251 G4double aXsection;
252
253 // MC integration loop
254 G4int counter = 0;
255 G4double buffer = 0;
256 G4int size = G4int(std::max(10., aT/60*kelvin));
257 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
258 G4double neutronVMag = neutronVelocity.mag();
259
260 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
261 {
262 if(counter) buffer = result/counter;
263 while (counter<size) // Loop checking, 11.05.2015, T. Koi
264 {
265 ++counter;
266 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
267 boosted.Lorentz(theNeutron, aThermalNuc);
268 G4double theEkin = boosted.GetKineticEnergy();
269 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
270 // velocity correction, or luminosity factor...
271 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
272 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
273 result += aXsection;
274 }
275 size += size;
276 }
277 result /= counter;
278/*
279 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
280 G4cout << " result " << result << " "
281 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
282 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
283*/
284 return result;
285}
286
288{
290}
291
293{
295}
296
298{
299 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for radiative capture 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: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
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