Geant4  10.02.p03
G4ParticleHPFissionData.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. Ko
35 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
36 //
38 #include "G4ParticleHPManager.hh"
39 #include "G4PhysicalConstants.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "G4Neutron.hh"
42 #include "G4ElementTable.hh"
43 #include "G4ParticleHPData.hh"
44 #include "G4ParticleHPManager.hh"
45 #include "G4Pow.hh"
46 
48 :G4VCrossSectionDataSet("NeutronHPFissionXS")
49 {
50  SetMinKinEnergy( 0*MeV );
51  SetMaxKinEnergy( 20*MeV );
52 
53  theCrossSections = 0;
54  onFlightDB = true;
55  instanceOfWorker = false;
57  instanceOfWorker = true;
58  }
59  //BuildPhysicsTable(*G4Neutron::Neutron());
60 }
61 
63 {
64  if ( theCrossSections != NULL && instanceOfWorker != true ) {
66  delete theCrossSections;
67  theCrossSections = NULL;
68  }
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;
95  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
96  xs_cache = xs;
97  return xs;
98 }
99 
100 /*
101 G4bool G4ParticleHPFissionData::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 ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() ) {
114  G4cout << "Find a flag of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
115  G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of fission reaction of neutrons (<20MeV)." << G4endl;
116  onFlightDB = false;
117  }
118 
119  if(&aP!=G4Neutron::Neutron())
120  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
121 
122  if ( G4Threading::IsWorkerThread() ) {
124  return;
125  }
126 
127  size_t numberOfElements = G4Element::GetNumberOfElements();
128  //theCrossSections = new G4PhysicsTable( numberOfElements );
129  // TKDB
130  //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
131  if ( theCrossSections == NULL )
132  theCrossSections = new G4PhysicsTable( numberOfElements );
133  else
135 
136  // make a PhysicsVector for each element
137 
138  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
139  for( size_t i=0; i<numberOfElements; ++i )
140  {
142  Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
143  theCrossSections->push_back(physVec);
144  }
145 
147 }
148 
150 {
151  if(&aP!=G4Neutron::Neutron())
152  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
153 
154 //
155 // Dump element based cross section
156 // range 10e-5 eV to 20 MeV
157 // 10 point per decade
158 // in barn
159 //
160 
161  G4cout << G4endl;
162  G4cout << G4endl;
163  G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
164  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
165  G4cout << G4endl;
166  G4cout << "Name of Element" << G4endl;
167  G4cout << "Energy[eV] XS[barn]" << G4endl;
168  G4cout << G4endl;
169 
170  size_t numberOfElements = G4Element::GetNumberOfElements();
171  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
172 
173  for ( size_t i = 0 ; i < numberOfElements ; ++i )
174  {
175 
176  G4cout << (*theElementTable)[i]->GetName() << G4endl;
177 
178  if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
179  {
180  G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
181  G4cout << G4endl;
182  continue;
183  }
184 
185  G4int ie = 0;
186 
187  for ( ie = 0 ; ie < 130 ; ie++ )
188  {
189  G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
190  G4bool outOfRange = false;
191 
192  if ( eKinetic < 20*MeV )
193  {
194  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
195  }
196 
197  }
198 
199  G4cout << G4endl;
200  }
201 
202  //G4cout << "G4ParticleHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
203 }
204 
205 #include "G4NucleiProperties.hh"
206 
209 {
210  G4double result = 0;
211  if(anE->GetZ()<88) return result;
212  G4bool outOfRange;
213  G4int index = anE->GetIndex();
214 
215 // 100729 TK add safety
216 if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
217 
218  // prepare neutron
219  G4double eKinetic = aP->GetKineticEnergy();
220  G4ReactionProduct theNeutronRP( aP->GetDefinition() );
221  theNeutronRP.SetMomentum( aP->GetMomentum() );
222  theNeutronRP.SetKineticEnergy( eKinetic );
223 
224  if ( !onFlightDB ) {
225  //NEGLECT_DOPPLER
226  G4double factor = 1.0;
227  if ( eKinetic < aT * k_Boltzmann ) {
228  // below 0.1 eV neutrons
229  // Have to do some, but now just igonre.
230  // Will take care after performance check.
231  // factor = factor * targetV;
232  }
233  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
234  }
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  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
244 
245  G4ReactionProduct boosted;
246  G4double aXsection;
247 
248  // MC integration loop
249  G4int counter = 0;
250  G4double buffer = 0;
251  G4int size = G4int(std::max(10., aT/60*kelvin));
252  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutronRP.GetMomentum();
253  G4double neutronVMag = neutronVelocity.mag();
254 
255  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
256  {
257  if(counter) buffer = result/counter;
258  while (counter<size) // Loop checking, 11.05.2015, T. Koi
259  {
260  counter ++;
261  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
262  boosted.Lorentz(theNeutronRP, aThermalNuc);
263  G4double theEkin = boosted.GetKineticEnergy();
264  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
265  // velocity correction.
266  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
267  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
268  result += aXsection;
269  }
270  size += size;
271  }
272  result /= counter;
273  return result;
274 }
275 
277 {
279 }
281 {
283 }
284 void G4ParticleHPFissionData::CrossSectionDescription(std::ostream& outFile) const
285 {
286  outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n" ;
287 }
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 RegisterFissionCrossSections(G4PhysicsTable *val)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
Int_t index
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:143
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
static const G4double eps
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
#define G4ThreadLocal
Definition: tls.hh:89
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * GetFissionCrossSections()
G4double GetN() const
Definition: G4Element.hh:134
string material
Definition: eplot.py:19
void BuildPhysicsTable(const G4ParticleDefinition &)
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
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
void clearAndDestroy()
G4ThreeVector GetMomentum() const
virtual void CrossSectionDescription(std::ostream &) const
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)