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