Geant4  10.02
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 #include "G4Threading.hh"
48 #include "G4Pow.hh"
49 
51 :G4VCrossSectionDataSet("NeutronHPCaptureXS")
52 {
53  SetMinKinEnergy( 0*MeV );
54  SetMaxKinEnergy( 20*MeV );
55 
56  theCrossSections = 0;
57  onFlightDB = true;
58 
59  //BuildPhysicsTable(*G4Neutron::Neutron());
60 }
61 
63 {
64  if ( theCrossSections != NULL ) {
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  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
91  return xs;
92 }
93 
94 /*
95 G4bool G4ParticleHPCaptureData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
96 {
97  G4bool result = true;
98  G4double eKin = aP->GetKineticEnergy();
99  if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
100  return result;
101 }
102 */
103 
105 {
106  if(&aP!=G4Neutron::Neutron())
107  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
108 
109 //080428
110  if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
111  {
112  G4cout << "Find a flag of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
113  G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of capture reaction of neutrons (<20MeV)." << G4endl;
114  onFlightDB = false;
115  }
116 
117  if ( G4Threading::IsWorkerThread() ) {
119  return;
120  }
121 
122  size_t numberOfElements = G4Element::GetNumberOfElements();
123  // G4cout << "CALLED G4ParticleHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl;
124  // TKDB
125  //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
126  if ( theCrossSections == NULL )
127  theCrossSections = new G4PhysicsTable( numberOfElements );
128  else
130 
131  // make a PhysicsVector for each element
132 
133  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
134  for( size_t i=0; i<numberOfElements; ++i )
135  {
136  if(getenv("CaptureDataIndexDebug"))
137  {
138  G4int index_debug = ((*theElementTable)[i])->GetIndex();
139  G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl;
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 << "Capture 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  G4int ie = 0;
179 
180  for ( ie = 0 ; ie < 130 ; ie++ )
181  {
182  G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
183  G4bool outOfRange = false;
184 
185  if ( eKinetic < 20*MeV )
186  {
187  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
188  }
189 
190  }
191 
192  G4cout << G4endl;
193  }
194 
195 
196 // G4cout << "G4ParticleHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
197 }
198 
199 #include "G4NucleiProperties.hh"
200 
203 {
204  G4double result = 0;
205  G4bool outOfRange;
206  G4int index = anE->GetIndex();
207 
208  // prepare neutron
209  G4double eKinetic = aP->GetKineticEnergy();
210 
211  if ( !onFlightDB )
212  {
213  //NEGLECT_DOPPLER
214  G4double factor = 1.0;
215  if ( eKinetic < aT * k_Boltzmann )
216  {
217  // below 0.1 eV neutrons
218  // Have to do some, but now just igonre.
219  // Will take care after performance check.
220  // factor = factor * targetV;
221  }
222  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
223  }
224 
225  G4ReactionProduct theNeutron( aP->GetDefinition() );
226  theNeutron.SetMomentum( aP->GetMomentum() );
227  theNeutron.SetKineticEnergy( eKinetic );
228 
229  // prepare thermal nucleus
230  G4Nucleus aNuc;
231  G4double eps = 0.0001;
232  G4double theA = anE->GetN();
233  G4double theZ = anE->GetZ();
234  G4double eleMass;
235  eleMass = G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) / G4Neutron::Neutron()->GetPDGMass();
236 
237  G4ReactionProduct boosted;
238  G4double aXsection;
239 
240  // MC integration loop
241  G4int counter = 0;
242  G4double buffer = 0;
243  G4int size = G4int(std::max(10., aT/60*kelvin));
244  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
245  G4double neutronVMag = neutronVelocity.mag();
246 
247  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
248  {
249  if(counter) buffer = result/counter;
250  while (counter<size) // Loop checking, 11.05.2015, T. Koi
251  {
252  counter ++;
253  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
254  boosted.Lorentz(theNeutron, aThermalNuc);
255  G4double theEkin = boosted.GetKineticEnergy();
256  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
257  // velocity correction, or luminosity factor...
258  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
259  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
260  result += aXsection;
261  }
262  size += size;
263  }
264  result /= counter;
265 /*
266  // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
267  G4cout << " result " << result << " "
268  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
269  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
270 */
271  return result;
272 }
273 
275 {
277 }
279 {
281 }
282 void G4ParticleHPCaptureData::CrossSectionDescription(std::ostream& outFile) const
283 {
284  outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for radiative capture reaction of neutrons below 20MeV\n" ;
285 }
static G4ParticleHPManager * GetInstance()
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
static const double MeV
Definition: G4SIunits.hh:211
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 *)
void RegisterCaptureCrossSections(G4PhysicsTable *val)
G4double GetZ() const
Definition: G4Element.hh:131
#define buffer
Definition: xmlparse.cc:628
static const G4double eps
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocal
Definition: tls.hh:89
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:143
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
G4PhysicsTable * GetCaptureCrossSections()
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:278
G4bool IsWorkerThread()
Definition: G4Threading.cc:129
G4double GetKineticEnergy() const
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
static const double eV
Definition: G4SIunits.hh:212
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
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
G4double GetMass() const
virtual void CrossSectionDescription(std::ostream &) const
void clearAndDestroy()
G4ThreeVector GetMomentum() const
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)