Geant4  10.02.p03
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 //
40 #include "G4ParticleHPManager.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4Neutron.hh"
44 #include "G4ElementTable.hh"
45 #include "G4ParticleHPData.hh"
46 #include "G4ParticleHPManager.hh"
47 #include "G4Pow.hh"
48 
50 :G4VCrossSectionDataSet("NeutronHPElasticXS")
51 {
52  SetMinKinEnergy( 0*MeV );
53  SetMaxKinEnergy( 20*MeV );
54 
55  theCrossSections = 0;
56  onFlightDB = true;
57  instanceOfWorker = false;
59  instanceOfWorker = true;
60  }
61 // BuildPhysicsTable( *G4Neutron::Neutron() );
62 }
63 
65 {
66  if ( theCrossSections != NULL && instanceOfWorker != true ) {
68  delete theCrossSections;
69  theCrossSections = NULL;
70  }
71 }
72 
74  G4int /*Z*/ , G4int /*A*/ ,
75  const G4Element* /*elm*/ ,
76  const G4Material* /*mat*/ )
77 {
78 
79  G4double eKin = dp->GetKineticEnergy();
80  if ( eKin > GetMaxKinEnergy()
81  || eKin < GetMinKinEnergy()
82  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
83 
84  return true;
85 }
86 
88  G4int /*Z*/ , G4int /*A*/ ,
89  const G4Isotope* /*iso*/ ,
90  const G4Element* element ,
91  const G4Material* material )
92 {
93  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
94 
95  ke_cache = dp->GetKineticEnergy();
96  element_cache = element;
98  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
99  xs_cache = xs;
100  return xs;
101 }
102 
103 /*
104 G4bool G4ParticleHPElasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
105 {
106  G4bool result = true;
107  G4double eKin = aP->GetKineticEnergy();
108  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
109  return result;
110 }
111 */
112 
114 {
115 
116  if(&aP!=G4Neutron::Neutron())
117  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
118 
119 //080428
120  if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
121  {
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 elastic scattering of neutrons (<20MeV)." << G4endl;
124  onFlightDB = false;
125  }
126 
127  if ( G4Threading::IsWorkerThread() ) {
129  return;
130  }
131 
132  size_t numberOfElements = G4Element::GetNumberOfElements();
133 // TKDB
134  //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
135  if ( theCrossSections == NULL )
136  theCrossSections = new G4PhysicsTable( numberOfElements );
137  else
139 
140  // make a PhysicsVector for each element
141 
142  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
143  for( size_t i=0; i<numberOfElements; ++i )
144  {
146  Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
147  theCrossSections->push_back(physVec);
148  }
149 
151 }
152 
154 {
155  if(&aP!=G4Neutron::Neutron())
156  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
157 
158 //
159 // Dump element based cross section
160 // range 10e-5 eV to 20 MeV
161 // 10 point per decade
162 // in barn
163 //
164 
165  G4cout << G4endl;
166  G4cout << G4endl;
167  G4cout << "Elastic Cross Section of Neutron HP"<< G4endl;
168  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
169  G4cout << G4endl;
170  G4cout << "Name of Element" << G4endl;
171  G4cout << "Energy[eV] XS[barn]" << G4endl;
172  G4cout << G4endl;
173 
174  size_t numberOfElements = G4Element::GetNumberOfElements();
175  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
176 
177  for ( size_t i = 0 ; i < numberOfElements ; ++i )
178  {
179 
180  G4cout << (*theElementTable)[i]->GetName() << G4endl;
181 
182  G4int ie = 0;
183 
184  for ( ie = 0 ; ie < 130 ; ie++ )
185  {
186  G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
187  G4bool outOfRange = false;
188 
189  if ( eKinetic < 20*MeV )
190  {
191  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
192  }
193 
194  }
195 
196  G4cout << G4endl;
197  }
198 
199 
200 // G4cout << "G4ParticleHPElasticData::DumpPhysicsTable still to be implemented"<<G4endl;
201 }
202 
203 #include "G4Nucleus.hh"
204 #include "G4NucleiProperties.hh"
205 #include "G4Neutron.hh"
206 #include "G4Electron.hh"
207 
210 {
211  G4double result = 0;
212  G4bool outOfRange;
213  G4int index = anE->GetIndex();
214 
215  // prepare neutron
216  G4double eKinetic = aP->GetKineticEnergy();
217 
218  if ( !onFlightDB )
219  {
220  //NEGLECT_DOPPLER_B.
221  G4double factor = 1.0;
222  if ( eKinetic < aT * k_Boltzmann )
223  {
224  // below 0.1 eV neutrons
225  // Have to do some, but now just igonre.
226  // Will take care after performance check.
227  // factor = factor * targetV;
228  }
229  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
230  }
231 
232  G4ReactionProduct theNeutron( aP->GetDefinition() );
233  theNeutron.SetMomentum( aP->GetMomentum() );
234  theNeutron.SetKineticEnergy( eKinetic );
235 
236  // prepare thermal nucleus
237  G4Nucleus aNuc;
238  G4double eps = 0.0001;
239  G4double theA = anE->GetN();
240  G4double theZ = anE->GetZ();
241  G4double eleMass;
242 
243 
244  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
246 
247  G4ReactionProduct boosted;
248  G4double aXsection;
249 
250  // MC integration loop
251  G4int counter = 0;
252  G4double buffer = 0;
253  G4int size = G4int(std::max(10., aT/60*kelvin));
254  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
255  G4double neutronVMag = neutronVelocity.mag();
256 
257  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
258  {
259  if(counter) buffer = result/counter;
260  while (counter<size) // Loop checking, 11.05.2015, T. Koi
261  {
262  counter ++;
263  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
264  boosted.Lorentz(theNeutron, aThermalNuc);
265  G4double theEkin = boosted.GetKineticEnergy();
266  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
267  // velocity correction.
268  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
269  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
270  result += aXsection;
271  }
272  size += size;
273  }
274  result /= counter;
275 /*
276  // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
277  G4cout << " result " << result << " "
278  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
279  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
280 */
281  return result;
282 }
283 
286 {
288 }
289 
291 SetVerboseLevel( G4int newValue )
292 {
294 }
295 void G4ParticleHPElasticData::CrossSectionDescription(std::ostream& outFile) const
296 {
297  outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for elastic reaction of neutrons below 20MeV\n" ;
298 }
static G4ParticleHPManager * GetInstance()
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double GetMass() const
static const double MeV
Definition: G4SIunits.hh:211
static G4double GetNuclearMass(const G4double A, const G4double Z)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
Int_t index
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:143
void BuildPhysicsTable(const G4ParticleDefinition &)
size_t GetIndex() const
Definition: G4Element.hh:181
void SetMomentum(const G4double x, const G4double y, const G4double z)
void push_back(G4PhysicsVector *)
void DumpPhysicsTable(const G4ParticleDefinition &)
#define buffer
Definition: xmlparse.cc:628
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
static const G4double eps
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
G4double GetN() const
Definition: G4Element.hh:134
string material
Definition: eplot.py:19
void RegisterElasticCrossSections(G4PhysicsTable *val)
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4GLOB_DLL std::ostream G4cout
float k_Boltzmann
Definition: hepunit.py:299
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
G4PhysicsTable * GetElasticCrossSections()
void SetMinKinEnergy(G4double value)
bool G4bool
Definition: G4Types.hh:79
double mag() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static const double kelvin
Definition: G4SIunits.hh:278
G4bool IsWorkerThread()
Definition: G4Threading.cc:135
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
static const double eV
Definition: G4SIunits.hh:212
static const G4double factor
void SetMaxKinEnergy(G4double value)
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
G4double GetKineticEnergy() const
static const double barn
Definition: G4SIunits.hh:104
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
G4double GetZ() const
Definition: G4Element.hh:131
virtual void CrossSectionDescription(std::ostream &) const
void clearAndDestroy()
G4ThreeVector GetMomentum() const