Geant4  10.01.p02
G4NeutronHPCaptureData.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 #include "G4Threading.hh"
46 
48 :G4VCrossSectionDataSet("NeutronHPCaptureXS")
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  theCrossSections = 0;
59  onFlightDB = true;
60 
61  //BuildPhysicsTable(*G4Neutron::Neutron());
62 }
63 
65 {
66  if ( theCrossSections != 0 ) {
68  delete theCrossSections;
69  theCrossSections = 0;
70  }
71 }
72 
74  G4int /*Z*/ , G4int /*A*/ ,
75  const G4Element* /*elm*/ ,
76  const G4Material* /*mat*/ )
77 {
78  G4double eKin = dp->GetKineticEnergy();
79  if ( eKin > GetMaxKinEnergy()
80  || eKin < GetMinKinEnergy()
81  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
82 
83  return true;
84 }
85 
87  G4int /*Z*/ , G4int /*A*/ ,
88  const G4Isotope* /*iso*/ ,
89  const G4Element* element ,
90  const G4Material* material )
91 {
92  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
93 
94  ke_cache = dp->GetKineticEnergy();
95  element_cache = element;
96  material_cache = material;
97  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
98  xs_cache = xs;
99  return xs;
100 }
101 
102 /*
103 G4bool G4NeutronHPCaptureData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
104 {
105  G4bool result = true;
106  G4double eKin = aP->GetKineticEnergy();
107  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
108  return result;
109 }
110 */
111 
113 {
114  if(&aP!=G4Neutron::Neutron())
115  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
116 
117 //080428
118  //if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
119  if ( G4NeutronHPManager::GetInstance()->GetNeglectDoppler() )
120  {
121  G4cout << "Find a flag 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  if ( G4Threading::IsWorkerThread() ) {
128  return;
129  }
130 
131  size_t numberOfElements = G4Element::GetNumberOfElements();
132  // G4cout << "CALLED G4NeutronHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl;
133  // TKDB
134  //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
135  if ( theCrossSections == NULL )
136  theCrossSections = new G4PhysicsTable( numberOfElements );
137  else
139 
140  // make a PhysicsVector for each element
141 
142  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
143  for( size_t i=0; i<numberOfElements; ++i )
144  {
145  if(getenv("CaptureDataIndexDebug"))
146  {
147  G4int index_debug = ((*theElementTable)[i])->GetIndex();
148  G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl;
149  }
151  Instance()->MakePhysicsVector((*theElementTable)[i], this);
152  theCrossSections->push_back(physVec);
153  }
154 
156 }
157 
159 {
160  if(&aP!=G4Neutron::Neutron())
161  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
162 
163 //
164 // Dump element based cross section
165 // range 10e-5 eV to 20 MeV
166 // 10 point per decade
167 // in barn
168 //
169 
170  G4cout << G4endl;
171  G4cout << G4endl;
172  G4cout << "Capture Cross Section of Neutron HP"<< G4endl;
173  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
174  G4cout << G4endl;
175  G4cout << "Name of Element" << G4endl;
176  G4cout << "Energy[eV] XS[barn]" << G4endl;
177  G4cout << G4endl;
178 
179  size_t numberOfElements = G4Element::GetNumberOfElements();
180  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
181 
182  for ( size_t i = 0 ; i < numberOfElements ; ++i )
183  {
184 
185  G4cout << (*theElementTable)[i]->GetName() << G4endl;
186 
187  G4int ie = 0;
188 
189  for ( ie = 0 ; ie < 130 ; ie++ )
190  {
191  G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
192  G4bool outOfRange = false;
193 
194  if ( eKinetic < 20*MeV )
195  {
196  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
197  }
198 
199  }
200 
201  G4cout << G4endl;
202  }
203 
204 
205 // G4cout << "G4NeutronHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
206 }
207 
208 #include "G4NucleiProperties.hh"
209 
212 {
213  G4double result = 0;
214  G4bool outOfRange;
215  G4int index = anE->GetIndex();
216 
217  // prepare neutron
218  G4double eKinetic = aP->GetKineticEnergy();
219 
220  if ( !onFlightDB )
221  {
222  //NEGLECT_DOPPLER
223  G4double factor = 1.0;
224  if ( eKinetic < aT * k_Boltzmann )
225  {
226  // below 0.1 eV neutrons
227  // Have to do some, but now just igonre.
228  // Will take care after performance check.
229  // factor = factor * targetV;
230  }
231  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
232  }
233 
234  G4ReactionProduct theNeutron( aP->GetDefinition() );
235  theNeutron.SetMomentum( aP->GetMomentum() );
236  theNeutron.SetKineticEnergy( eKinetic );
237 
238  // prepare thermal nucleus
239  G4Nucleus aNuc;
240  G4double eps = 0.0001;
241  G4double theA = anE->GetN();
242  G4double theZ = anE->GetZ();
243  G4double eleMass;
244  eleMass = G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) / G4Neutron::Neutron()->GetPDGMass();
245 
246  G4ReactionProduct boosted;
247  G4double aXsection;
248 
249  // MC integration loop
250  G4int counter = 0;
251  G4double buffer = 0;
252  G4int size = G4int(std::max(10., aT/60*kelvin));
253  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
254  G4double neutronVMag = neutronVelocity.mag();
255 
256  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
257  {
258  if(counter) buffer = result/counter;
259  while (counter<size)
260  {
261  counter ++;
262  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
263  boosted.Lorentz(theNeutron, aThermalNuc);
264  G4double theEkin = boosted.GetKineticEnergy();
265  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
266  // velocity correction, or luminosity factor...
267  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
268  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
269  result += aXsection;
270  }
271  size += size;
272  }
273  result /= counter;
274 /*
275  // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
276  G4cout << " result " << result << " "
277  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
278  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
279 */
280  return result;
281 }
282 
284 {
286 }
288 {
290 }
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
const G4Material * material_cache
G4double GetN() const
Definition: G4Element.hh:134
void BuildPhysicsTable(const G4ParticleDefinition &)
static G4NeutronHPManager * GetInstance()
void RegisterCaptureCrossSections(G4PhysicsTable *val)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void push_back(G4PhysicsVector *)
G4double GetZ() const
Definition: G4Element.hh:131
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
#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
G4PhysicsTable * GetCaptureCrossSections()
void DumpPhysicsTable(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 G4NeutronHPData * Instance()
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
static const double kelvin
Definition: G4SIunits.hh:260
G4bool IsWorkerThread()
Definition: G4Threading.cc:128
G4double GetKineticEnergy() const
void SetVerboseLevel(G4int i)
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
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
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