Geant4  10.00.p02
G4NeutronHPInelasticData.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("NeutronHPInelasticXS")
48 {
49 
50  SetMinKinEnergy( 0*MeV );
51  SetMaxKinEnergy( 20*MeV );
52 
53  ke_cache = 0.0;
54  xs_cache = 0.0;
55  element_cache = NULL;
56  material_cache = NULL;
57 
58  onFlightDB = true;
59  theCrossSections = 0;
60  //BuildPhysicsTable(*G4Neutron::Neutron());
61 }
62 
64 {
65  //This should now be avoided in destructors because content of
66  //theCrossSections is managed by allocators
67  //if ( theCrossSections != 0 ) theCrossSections->clearAndDestroy();
68  delete theCrossSections;
69 }
70 
72  G4int /*Z*/ , G4int /*A*/ ,
73  const G4Element* /*elm*/ ,
74  const G4Material* /*mat*/ )
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  //return GetCrossSection( dp , element , material->GetTemperature() );
99 }
100 
101 /*
102 G4bool G4NeutronHPInelasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
103 {
104  G4bool result = true;
105  G4double eKin = aP->GetKineticEnergy();
106  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
107  return result;
108 }
109 */
110 
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 inelastic scattering of neutrons (<20MeV)." << G4endl;
121  onFlightDB = false;
122  }
123 
124  size_t numberOfElements = G4Element::GetNumberOfElements();
125 // theCrossSections = new G4PhysicsTable( numberOfElements );
126 // TKDB
127  //if ( theCrossSections == 0 )
128  //{ theCrossSections = new G4PhysicsTable( numberOfElements ); }
129  if ( theCrossSections == NULL )
130  theCrossSections = new G4PhysicsTable( numberOfElements );
131  else
133 
134  // make a PhysicsVector for each element
135 
136  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
137  for( size_t i=0; i<numberOfElements; ++i )
138  {
140  Instance()->MakePhysicsVector((*theElementTable)[i], this);
141  theCrossSections->push_back(physVec);
142  }
143 }
144 
146 {
147  if(&aP!=G4Neutron::Neutron())
148  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
149 
150 //
151 // Dump element based cross section
152 // range 10e-5 eV to 20 MeV
153 // 10 point per decade
154 // in barn
155 //
156 
157  G4cout << G4endl;
158  G4cout << G4endl;
159  G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
160  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
161  G4cout << G4endl;
162  G4cout << "Name of Element" << G4endl;
163  G4cout << "Energy[eV] XS[barn]" << G4endl;
164  G4cout << G4endl;
165 
166  size_t numberOfElements = G4Element::GetNumberOfElements();
167  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
168 
169  for ( size_t i = 0 ; i < numberOfElements ; ++i )
170  {
171 
172  G4cout << (*theElementTable)[i]->GetName() << G4endl;
173 
174  G4int ie = 0;
175 
176  for ( ie = 0 ; ie < 130 ; ie++ )
177  {
178  G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
179  G4bool outOfRange = false;
180 
181  if ( eKinetic < 20*MeV )
182  {
183  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
184  }
185 
186  }
187 
188  G4cout << G4endl;
189  }
190 
191  //G4cout << "G4NeutronHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
192 }
193 
194 #include "G4NucleiProperties.hh"
195 
198 {
199  G4double result = 0;
200  G4bool outOfRange;
201  G4int index = anE->GetIndex();
202 
203  // prepare neutron
204  G4double eKinetic = aP->GetKineticEnergy();
205 
206  // T. K.
207 //if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
208 //080428
209  if ( !onFlightDB )
210  {
211  G4double factor = 1.0;
212  if ( eKinetic < aT * k_Boltzmann )
213  {
214  // below 0.1 eV neutrons
215  // Have to do some, but now just igonre.
216  // Will take care after performance check.
217  // factor = factor * targetV;
218  }
219  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
220  }
221 
222  G4ReactionProduct theNeutron( aP->GetDefinition() );
223  theNeutron.SetMomentum( aP->GetMomentum() );
224  theNeutron.SetKineticEnergy( eKinetic );
225 
226  // prepare thermal nucleus
227  G4Nucleus aNuc;
228  G4double eps = 0.0001;
229  G4double theA = anE->GetN();
230  G4double theZ = anE->GetZ();
231  G4double eleMass;
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  G4int failCount = 0;
241  G4double buffer = 0;
242  G4int size = G4int(std::max(10., aT/60*kelvin));
243  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
244  G4double neutronVMag = neutronVelocity.mag();
245 
246  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
247  {
248  if(counter) buffer = result/counter;
249  while (counter<size)
250  {
251  counter ++;
252  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
253  boosted.Lorentz(theNeutron, aThermalNuc);
254  G4double theEkin = boosted.GetKineticEnergy();
255  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
256  if(aXsection <0)
257  {
258  if(failCount<1000)
259  {
260  failCount++;
261  counter--;
262  continue;
263  }
264  else
265  {
266  aXsection = 0;
267  }
268  }
269  // velocity correction.
270  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
271  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
272  result += aXsection;
273  }
274  size += size;
275  }
276  result /= counter;
277 /*
278  // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
279  G4cout << " result " << result << " "
280  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
281  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
282 */
283  return result;
284 }
285 
287 {
289 }
291 {
293 }
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4NeutronHPFissionData *theP)
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
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
void SetMinKinEnergy(G4double value)
bool G4bool
Definition: G4Types.hh:79
void DumpPhysicsTable(const G4ParticleDefinition &)
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
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
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
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double GetMass() const
void clearAndDestroy()
G4ThreeVector GetMomentum() const