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