Geant4  10.00.p02
G4NeutronHPElasticData.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 //
38 #include "G4NeutronHPManager.hh"
39 #include "G4PhysicalConstants.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4Neutron.hh"
42 #include "G4ElementTable.hh"
43 #include "G4NeutronHPData.hh"
44 #include "G4NeutronHPManager.hh"
45 
47 :G4VCrossSectionDataSet("NeutronHPElasticXS")
48 {
49  SetMinKinEnergy( 0*MeV );
50  SetMaxKinEnergy( 20*MeV );
51 
52  ke_cache = 0.0;
53  xs_cache = 0.0;
54  element_cache = NULL;
55  material_cache = NULL;
56 
57  theCrossSections = 0;
58  onFlightDB = true;
59 // BuildPhysicsTable( *G4Neutron::Neutron() );
60 }
61 
63 {
64  //This should be now avoided in destructor because these are
65  //handled by allocators
66  //if ( theCrossSections != 0 ) theCrossSections->clearAndDestroy();
67  delete theCrossSections;
68 }
69 
71  G4int /*Z*/ , G4int /*A*/ ,
72  const G4Element* /*elm*/ ,
73  const G4Material* /*mat*/ )
74 {
75 
76  G4double eKin = dp->GetKineticEnergy();
77  if ( eKin > GetMaxKinEnergy()
78  || eKin < GetMinKinEnergy()
79  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
80 
81  return true;
82 }
83 
85  G4int /*Z*/ , G4int /*A*/ ,
86  const G4Isotope* /*iso*/ ,
87  const G4Element* element ,
88  const G4Material* material )
89 {
90  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
91 
92  ke_cache = dp->GetKineticEnergy();
93  element_cache = element;
94  material_cache = material;
95  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
96  xs_cache = xs;
97  return xs;
98 }
99 
100 /*
101 G4bool G4NeutronHPElasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
102 {
103  G4bool result = true;
104  G4double eKin = aP->GetKineticEnergy();
105  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
106  return result;
107 }
108 */
109 
111 {
112 
113  if(&aP!=G4Neutron::Neutron())
114  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
115 
116 //080428
117  if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
118  {
119  G4cout << "Find environment variable of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
120  G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of elastic scattering of neutrons (<20MeV)." << G4endl;
121  onFlightDB = false;
122  }
123 
124  size_t numberOfElements = G4Element::GetNumberOfElements();
125 // TKDB
126  //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
127  if ( theCrossSections == NULL )
128  theCrossSections = new G4PhysicsTable( numberOfElements );
129  else
131 
132  // make a PhysicsVector for each element
133 
134  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
135  for( size_t i=0; i<numberOfElements; ++i )
136  {
138  Instance()->MakePhysicsVector((*theElementTable)[i], this);
139  theCrossSections->push_back(physVec);
140  }
141 }
142 
144 {
145  if(&aP!=G4Neutron::Neutron())
146  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
147 
148 //
149 // Dump element based cross section
150 // range 10e-5 eV to 20 MeV
151 // 10 point per decade
152 // in barn
153 //
154 
155  G4cout << G4endl;
156  G4cout << G4endl;
157  G4cout << "Elastic Cross Section of Neutron HP"<< G4endl;
158  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
159  G4cout << G4endl;
160  G4cout << "Name of Element" << G4endl;
161  G4cout << "Energy[eV] XS[barn]" << G4endl;
162  G4cout << G4endl;
163 
164  size_t numberOfElements = G4Element::GetNumberOfElements();
165  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
166 
167  for ( size_t i = 0 ; i < numberOfElements ; ++i )
168  {
169 
170  G4cout << (*theElementTable)[i]->GetName() << G4endl;
171 
172  G4int ie = 0;
173 
174  for ( ie = 0 ; ie < 130 ; ie++ )
175  {
176  G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
177  G4bool outOfRange = false;
178 
179  if ( eKinetic < 20*MeV )
180  {
181  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
182  }
183 
184  }
185 
186  G4cout << G4endl;
187  }
188 
189 
190 // G4cout << "G4NeutronHPElasticData::DumpPhysicsTable still to be implemented"<<G4endl;
191 }
192 
193 #include "G4Nucleus.hh"
194 #include "G4NucleiProperties.hh"
195 #include "G4Neutron.hh"
196 #include "G4Electron.hh"
197 
200 {
201  G4double result = 0;
202  G4bool outOfRange;
203  G4int index = anE->GetIndex();
204 
205  // prepare neutron
206  G4double eKinetic = aP->GetKineticEnergy();
207 
208  // T. K.
209 // if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
210 //080428
211  if ( !onFlightDB )
212  {
213  G4double factor = 1.0;
214  if ( eKinetic < aT * k_Boltzmann )
215  {
216  // below 0.1 eV neutrons
217  // Have to do some, but now just igonre.
218  // Will take care after performance check.
219  // factor = factor * targetV;
220  }
221  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
222  }
223 
224  G4ReactionProduct theNeutron( aP->GetDefinition() );
225  theNeutron.SetMomentum( aP->GetMomentum() );
226  theNeutron.SetKineticEnergy( eKinetic );
227 
228  // prepare thermal nucleus
229  G4Nucleus aNuc;
230  G4double eps = 0.0001;
231  G4double theA = anE->GetN();
232  G4double theZ = anE->GetZ();
233  G4double eleMass;
234 
235 
236  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
238 
239  G4ReactionProduct boosted;
240  G4double aXsection;
241 
242  // MC integration loop
243  G4int counter = 0;
244  G4double buffer = 0;
245  G4int size = G4int(std::max(10., aT/60*kelvin));
246  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
247  G4double neutronVMag = neutronVelocity.mag();
248 
249  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
250  {
251  if(counter) buffer = result/counter;
252  while (counter<size)
253  {
254  counter ++;
255  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
256  boosted.Lorentz(theNeutron, aThermalNuc);
257  G4double theEkin = boosted.GetKineticEnergy();
258  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
259  // velocity correction.
260  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
261  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
262  result += aXsection;
263  }
264  size += size;
265  }
266  result /= counter;
267 /*
268  // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
269  G4cout << " result " << result << " "
270  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
271  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
272 */
273  return result;
274 }
275 
278 {
280 }
281 
283 SetVerboseLevel( G4int newValue )
284 {
286 }
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4NeutronHPFissionData *theP)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
static const double MeV
Definition: G4SIunits.hh:193
static G4double GetNuclearMass(const G4double A, const G4double Z)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetN() const
Definition: G4Element.hh:134
static G4NeutronHPManager * GetInstance()
void SetMomentum(const G4double x, const G4double y, const G4double z)
void push_back(G4PhysicsVector *)
G4double GetZ() const
Definition: G4Element.hh:131
#define buffer
Definition: xmlparse.cc:611
static const G4double eps
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocal
Definition: tls.hh:52
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
int G4int
Definition: G4Types.hh:78
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
void DumpPhysicsTable(const G4ParticleDefinition &)
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
void SetMinKinEnergy(G4double value)
bool G4bool
Definition: G4Types.hh:79
size_t GetIndex() const
Definition: G4Element.hh:181
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NeutronHPData * Instance()
static const double kelvin
Definition: G4SIunits.hh:260
G4double GetKineticEnergy() const
void SetVerboseLevel(G4int i)
static const double eV
Definition: G4SIunits.hh:194
G4double GetPDGMass() const
void SetMaxKinEnergy(G4double value)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const G4Material * material_cache
G4ThreeVector GetMomentum() const
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:95
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
void BuildPhysicsTable(const G4ParticleDefinition &)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4double GetMass() const
void clearAndDestroy()
G4ThreeVector GetMomentum() const