Geant4  10.01.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 "G4SystemOfUnits.hh"
40 #include "G4Neutron.hh"
41 #include "G4ElementTable.hh"
42 #include "G4ParticleHPData.hh"
43 #include "G4ParticleHPManager.hh"
44 
46 :G4VCrossSectionDataSet("NeutronHPFissionXS")
47 {
48  SetMinKinEnergy( 0*MeV );
49  SetMaxKinEnergy( 20*MeV );
50 
51  ke_cache = 0.0;
52  xs_cache = 0.0;
53  element_cache = NULL;
54  material_cache = NULL;
55 
56  theCrossSections = 0;
57  //BuildPhysicsTable(*G4Neutron::Neutron());
58 }
59 
61 {
62  if ( theCrossSections != NULL ) {
64  delete theCrossSections;
65  theCrossSections = NULL;
66  }
67 }
68 
70  G4int /*Z*/ , G4int /*A*/ ,
71  const G4Element* /*elm*/ ,
72  const G4Material* /*mat*/ )
73 {
74  G4double eKin = dp->GetKineticEnergy();
75  if ( eKin > GetMaxKinEnergy()
76  || eKin < GetMinKinEnergy()
77  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
78 
79  return true;
80 }
81 
83  G4int /*Z*/ , G4int /*A*/ ,
84  const G4Isotope* /*iso*/ ,
85  const G4Element* element ,
86  const G4Material* material )
87 {
88  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
89 
90  ke_cache = dp->GetKineticEnergy();
91  element_cache = element;
92  material_cache = material;
93  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
94  xs_cache = xs;
95  return xs;
96 }
97 
98 /*
99 G4bool G4ParticleHPFissionData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
100 {
101  G4bool result = true;
102  G4double eKin = aP->GetKineticEnergy();
103  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
104  return result;
105 }
106 */
107 
109 {
110  if(&aP!=G4Neutron::Neutron())
111  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
112  size_t numberOfElements = G4Element::GetNumberOfElements();
113  //theCrossSections = new G4PhysicsTable( numberOfElements );
114  // TKDB
115  //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
116  if ( theCrossSections == NULL )
117  theCrossSections = new G4PhysicsTable( numberOfElements );
118  else
120 
121  // make a PhysicsVector for each element
122 
123  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
124  for( size_t i=0; i<numberOfElements; ++i )
125  {
127  Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
128  theCrossSections->push_back(physVec);
129  }
130 }
131 
133 {
134  if(&aP!=G4Neutron::Neutron())
135  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
136 
137 //
138 // Dump element based cross section
139 // range 10e-5 eV to 20 MeV
140 // 10 point per decade
141 // in barn
142 //
143 
144  G4cout << G4endl;
145  G4cout << G4endl;
146  G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
147  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
148  G4cout << G4endl;
149  G4cout << "Name of Element" << G4endl;
150  G4cout << "Energy[eV] XS[barn]" << G4endl;
151  G4cout << G4endl;
152 
153  size_t numberOfElements = G4Element::GetNumberOfElements();
154  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
155 
156  for ( size_t i = 0 ; i < numberOfElements ; ++i )
157  {
158 
159  G4cout << (*theElementTable)[i]->GetName() << G4endl;
160 
161  if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
162  {
163  G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
164  G4cout << G4endl;
165  continue;
166  }
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  //G4cout << "G4ParticleHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
186 }
187 
188 #include "G4NucleiProperties.hh"
189 
192 {
193  G4double result = 0;
194  if(anE->GetZ()<90) return result;
195  G4bool outOfRange;
196  G4int index = anE->GetIndex();
197 
198 // 100729 TK add safety
199 if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
200 
201  // prepare neutron
202  G4double eKinetic = aP->GetKineticEnergy();
203  G4ReactionProduct theNeutronRP( aP->GetDefinition() );
204  theNeutronRP.SetMomentum( aP->GetMomentum() );
205  theNeutronRP.SetKineticEnergy( eKinetic );
206 
207  // prepare thermal nucleus
208  G4Nucleus aNuc;
209  G4double eps = 0.0001;
210  G4double theA = anE->GetN();
211  G4double theZ = anE->GetZ();
212  G4double eleMass;
213  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
215 
216  G4ReactionProduct boosted;
217  G4double aXsection;
218 
219  // MC integration loop
220  G4int counter = 0;
221  G4double buffer = 0;
222  G4int size = G4int(std::max(10., aT/60*kelvin));
223  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutronRP.GetMomentum();
224  G4double neutronVMag = neutronVelocity.mag();
225 
226  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
227  {
228  if(counter) buffer = result/counter;
229  while (counter<size)
230  {
231  counter ++;
232  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
233  boosted.Lorentz(theNeutronRP, aThermalNuc);
234  G4double theEkin = boosted.GetKineticEnergy();
235  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
236  // velocity correction.
237  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
238  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
239  result += aXsection;
240  }
241  size += size;
242  }
243  result /= counter;
244  return result;
245 }
246 
248 {
250 }
252 {
254 }
static G4ParticleHPManager * GetInstance()
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
void SetMomentum(const G4double x, const G4double y, const G4double z)
void push_back(G4PhysicsVector *)
G4double GetZ() const
Definition: G4Element.hh:131
void DumpPhysicsTable(const G4ParticleDefinition &)
#define buffer
Definition: xmlparse.cc:611
static const G4double eps
G4ParticleDefinition * GetDefinition() const
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
#define G4ThreadLocal
Definition: tls.hh:89
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
int G4int
Definition: G4Types.hh:78
void BuildPhysicsTable(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 const double kelvin
Definition: G4SIunits.hh:260
G4double GetKineticEnergy() const
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
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
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
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
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)