Geant4  10.00.p02
G4NeutronHPFissionData.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 // 070618 fix memory leaking by T. Koi
31 // 071002 enable cross section dump by T. Koi
32 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
33 // 081124 Protect invalid read which caused run time errors by T. Koi
34 // 100729 Add safty for 0 lenght cross sections by T. Koi
35 
37 #include "G4NeutronHPManager.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "G4Neutron.hh"
40 #include "G4ElementTable.hh"
41 #include "G4NeutronHPData.hh"
42 #include "G4NeutronHPManager.hh"
43 
45 :G4VCrossSectionDataSet("NeutronHPFissionXS")
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  //BuildPhysicsTable(*G4Neutron::Neutron());
57 }
58 
60 {
61  //This should be now avoided in destructor because these are
62  //handled by allocators
63  //if ( theCrossSections != NULL ) theCrossSections->clearAndDestroy();
64  delete theCrossSections;
65 }
66 
68  G4int /*Z*/ , G4int /*A*/ ,
69  const G4Element* /*elm*/ ,
70  const G4Material* /*mat*/ )
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 G4NeutronHPFissionData::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  if(&aP!=G4Neutron::Neutron())
109  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
110  size_t numberOfElements = G4Element::GetNumberOfElements();
111  //theCrossSections = new G4PhysicsTable( numberOfElements );
112  // TKDB
113  //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
114  if ( theCrossSections == NULL )
115  theCrossSections = new G4PhysicsTable( numberOfElements );
116  else
118 
119  // make a PhysicsVector for each element
120 
121  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
122  for( size_t i=0; i<numberOfElements; ++i )
123  {
125  Instance()->MakePhysicsVector((*theElementTable)[i], this);
126  theCrossSections->push_back(physVec);
127  }
128 }
129 
131 {
132  if(&aP!=G4Neutron::Neutron())
133  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
134 
135 //
136 // Dump element based cross section
137 // range 10e-5 eV to 20 MeV
138 // 10 point per decade
139 // in barn
140 //
141 
142  G4cout << G4endl;
143  G4cout << G4endl;
144  G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
145  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
146  G4cout << G4endl;
147  G4cout << "Name of Element" << G4endl;
148  G4cout << "Energy[eV] XS[barn]" << G4endl;
149  G4cout << G4endl;
150 
151  size_t numberOfElements = G4Element::GetNumberOfElements();
152  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
153 
154  for ( size_t i = 0 ; i < numberOfElements ; ++i )
155  {
156 
157  G4cout << (*theElementTable)[i]->GetName() << G4endl;
158 
159  if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
160  {
161  G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
162  G4cout << G4endl;
163  continue;
164  }
165 
166  G4int ie = 0;
167 
168  for ( ie = 0 ; ie < 130 ; ie++ )
169  {
170  G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
171  G4bool outOfRange = false;
172 
173  if ( eKinetic < 20*MeV )
174  {
175  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
176  }
177 
178  }
179 
180  G4cout << G4endl;
181  }
182 
183  //G4cout << "G4NeutronHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
184 }
185 
186 #include "G4NucleiProperties.hh"
187 
190 {
191  G4double result = 0;
192  if(anE->GetZ()<90) return result;
193  G4bool outOfRange;
194  G4int index = anE->GetIndex();
195 
196 // 100729 TK add safety
197 if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
198 
199  // prepare neutron
200  G4double eKinetic = aP->GetKineticEnergy();
201  G4ReactionProduct theNeutron( aP->GetDefinition() );
202  theNeutron.SetMomentum( aP->GetMomentum() );
203  theNeutron.SetKineticEnergy( eKinetic );
204 
205  // prepare thermal nucleus
206  G4Nucleus aNuc;
207  G4double eps = 0.0001;
208  G4double theA = anE->GetN();
209  G4double theZ = anE->GetZ();
210  G4double eleMass;
211  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
213 
214  G4ReactionProduct boosted;
215  G4double aXsection;
216 
217  // MC integration loop
218  G4int counter = 0;
219  G4double buffer = 0;
220  G4int size = G4int(std::max(10., aT/60*kelvin));
221  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
222  G4double neutronVMag = neutronVelocity.mag();
223 
224  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
225  {
226  if(counter) buffer = result/counter;
227  while (counter<size)
228  {
229  counter ++;
230  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
231  boosted.Lorentz(theNeutron, aThermalNuc);
232  G4double theEkin = boosted.GetKineticEnergy();
233  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
234  // velocity correction.
235  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
236  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
237  result += aXsection;
238  }
239  size += size;
240  }
241  result /= counter;
242  return result;
243 }
244 
246 {
248 }
250 {
252 }
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
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
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
G4GLOB_DLL std::ostream G4cout
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
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
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NeutronHPData * Instance()
static const double kelvin
Definition: G4SIunits.hh:260
G4double GetKineticEnergy() const
const G4Material * material_cache
void BuildPhysicsTable(const G4ParticleDefinition &)
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
void DumpPhysicsTable(const G4ParticleDefinition &)
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