Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPThermalScatteringData.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 // Thermal Neutron Scattering
27 // Koi, Tatsumi (SCCS/SLAC)
28 //
29 // Class Description
30 // Cross Sections for a high precision (based on evaluated data
31 // libraries) description of themal neutron scattering below 4 eV;
32 // Based on Thermal neutron scattering files
33 // from the evaluated nuclear data files ENDF/B-VI, Release2
34 // To be used in your physics list in case you need this physics.
35 // In this case you want to register an object of this class with
36 // the corresponding process.
37 // Class Description - End
38 
39 // 15-Nov-06 First implementation is done by T. Koi (SLAC/SCCS)
40 // 070625 implement clearCurrentXSData to fix memory leaking by T. Koi
41 
42 #include <list>
43 #include <algorithm>
44 
46 #include "G4NeutronHPManager.hh"
47 
48 #include "G4SystemOfUnits.hh"
49 #include "G4Neutron.hh"
50 #include "G4ElementTable.hh"
51 //#include "G4NeutronHPData.hh"
52 
54 :G4VCrossSectionDataSet("NeutronHPThermalScatteringData")
55 {
56 // Upper limit of neutron energy
57  emax = 4*eV;
58  SetMinKinEnergy( 0*MeV );
59  SetMaxKinEnergy( emax );
60 
61  ke_cache = 0.0;
62  xs_cache = 0.0;
63  element_cache = NULL;
64  material_cache = NULL;
65 
66  indexOfThermalElement.clear();
67 
69 
70  //BuildPhysicsTable( *G4Neutron::Neutron() );
71 }
72 
74 {
75 
76  clearCurrentXSData();
77 
78  delete names;
79 }
80 
82  G4int /*Z*/ , G4int /*A*/ ,
83  const G4Element* element ,
84  const G4Material* material )
85 {
86  G4double eKin = dp->GetKineticEnergy();
87  if ( eKin > 4.0*eV //GetMaxKinEnergy()
88  || eKin < 0 //GetMinKinEnergy()
89  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
90 
91  if ( dic.find( std::pair < const G4Material* , const G4Element* > ( (G4Material*)NULL , element ) ) != dic.end()
92  || dic.find( std::pair < const G4Material* , const G4Element* > ( material , element ) ) != dic.end() ) return true;
93 
94  return false;
95 
96 // return IsApplicable( dp , element );
97 /*
98  G4double eKin = dp->GetKineticEnergy();
99  if ( eKin > 4.0*eV //GetMaxKinEnergy()
100  || eKin < 0 //GetMinKinEnergy()
101  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
102  return true;
103 */
104 }
105 
107  G4int /*Z*/ , G4int /*A*/ ,
108  const G4Isotope* /*iso*/ ,
109  const G4Element* element ,
110  const G4Material* material )
111 {
112  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
113 
114  ke_cache = dp->GetKineticEnergy();
115  element_cache = element;
116  material_cache = material;
117  //G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
118  G4double xs = GetCrossSection( dp , element , material );
119  xs_cache = xs;
120  return xs;
121  //return GetCrossSection( dp , element , material->GetTemperature() );
122 }
123 
124 void G4NeutronHPThermalScatteringData::clearCurrentXSData()
125 {
126  std::map< G4int , std::map< G4double , G4NeutronHPVector* >* >::iterator it;
127  std::map< G4double , G4NeutronHPVector* >::iterator itt;
128 
129  for ( it = coherent.begin() ; it != coherent.end() ; it++ )
130  {
131  if ( it->second != NULL )
132  {
133  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
134  {
135  delete itt->second;
136  }
137  }
138  delete it->second;
139  }
140 
141  for ( it = incoherent.begin() ; it != incoherent.end() ; it++ )
142  {
143  if ( it->second != NULL )
144  {
145  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
146  {
147  delete itt->second;
148  }
149  }
150  delete it->second;
151  }
152 
153  for ( it = inelastic.begin() ; it != inelastic.end() ; it++ )
154  {
155  if ( it->second != NULL )
156  {
157  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
158  {
159  delete itt->second;
160  }
161  }
162  delete it->second;
163  }
164 
165  coherent.clear();
166  incoherent.clear();
167  inelastic.clear();
168 
169 }
170 
171 
172 
174 {
175  G4bool result = false;
176 
177  G4double eKin = aP->GetKineticEnergy();
178  // Check energy
179  if ( eKin < emax )
180  {
181  // Check Particle Species
182  if ( aP->GetDefinition() == G4Neutron::Neutron() )
183  {
184  // anEle is one of Thermal elements
185  G4int ie = (G4int) anEle->GetIndex();
186  std::vector < G4int >::iterator it;
187  for ( it = indexOfThermalElement.begin() ; it != indexOfThermalElement.end() ; it++ )
188  {
189  if ( ie == *it ) return true;
190  }
191  }
192  }
193 
194 /*
195  if ( names->IsThisThermalElement ( anEle->GetName() ) )
196  {
197  // Check energy and projectile species
198  G4double eKin = aP->GetKineticEnergy();
199  if ( eKin < emax && aP->GetDefinition() == G4Neutron::Neutron() ) result = true;
200  }
201 */
202  return result;
203 }
204 
205 
207 {
208 
209  if ( &aP != G4Neutron::Neutron() )
210  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
211 
212  //std::map < std::pair < G4Material* , const G4Element* > , G4int > dic;
213  dic.clear();
214  clearCurrentXSData();
215  std::map < G4String , G4int > co_dic;
216 
217  //Searching Nist Materials
218  static G4ThreadLocal G4MaterialTable* theMaterialTable = 0 ; if (!theMaterialTable) theMaterialTable= G4Material::GetMaterialTable();
219  size_t numberOfMaterials = G4Material::GetNumberOfMaterials();
220  for ( size_t i = 0 ; i < numberOfMaterials ; i++ )
221  {
222  G4Material* material = (*theMaterialTable)[i];
223  size_t numberOfElements = material->GetNumberOfElements();
224  for ( size_t j = 0 ; j < numberOfElements ; j++ )
225  {
226  const G4Element* element = material->GetElement(j);
227  if ( names->IsThisThermalElement ( material->GetName() , element->GetName() ) )
228  {
229  G4int ts_ID_of_this_geometry;
230  G4String ts_ndl_name = names->GetTS_NDL_Name( material->GetName() , element->GetName() );
231  if ( co_dic.find ( ts_ndl_name ) != co_dic.end() )
232  {
233  ts_ID_of_this_geometry = co_dic.find ( ts_ndl_name ) -> second;
234  }
235  else
236  {
237  ts_ID_of_this_geometry = co_dic.size();
238  co_dic.insert ( std::pair< G4String , G4int >( ts_ndl_name , ts_ID_of_this_geometry ) );
239  }
240 
241  //G4cout << "Neutron HP Thermal Scattering Data : Registering a material-element pair of "
242  // << material->GetName() << " " << element->GetName()
243  // << " as internal thermal scattering id of " << ts_ID_of_this_geometry << "." << G4endl;
244 
245  dic.insert( std::pair < std::pair < G4Material* , const G4Element* > , G4int > ( std::pair < G4Material* , const G4Element* > ( material , element ) , ts_ID_of_this_geometry ) );
246  }
247  }
248  }
249 
250  //Searching TS Elements
251  static G4ThreadLocal G4ElementTable* theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
252  size_t numberOfElements = G4Element::GetNumberOfElements();
253  //size_t numberOfThermalElements = 0;
254  for ( size_t i = 0 ; i < numberOfElements ; i++ )
255  {
256  const G4Element* element = (*theElementTable)[i];
257  if ( names->IsThisThermalElement ( element->GetName() ) )
258  {
259  if ( names->IsThisThermalElement ( element->GetName() ) )
260  {
261  G4int ts_ID_of_this_geometry;
262  G4String ts_ndl_name = names->GetTS_NDL_Name( element->GetName() );
263  if ( co_dic.find ( ts_ndl_name ) != co_dic.end() )
264  {
265  ts_ID_of_this_geometry = co_dic.find ( ts_ndl_name ) -> second;
266  }
267  else
268  {
269  ts_ID_of_this_geometry = co_dic.size();
270  co_dic.insert ( std::pair< G4String , G4int >( ts_ndl_name , ts_ID_of_this_geometry ) );
271  }
272 
273  //G4cout << "Neutron HP Thermal Scattering: Registering an element of "
274  // << material->GetName() << " " << element->GetName()
275  // << " as internal thermal scattering id of " << ts_ID_of_this_geometry << "." << G4endl;
276 
277  dic.insert( std::pair < std::pair < const G4Material* , const G4Element* > , G4int > ( std::pair < const G4Material* , const G4Element* > ( (G4Material*)NULL , element ) , ts_ID_of_this_geometry ) );
278  }
279  }
280  }
281 
282  G4cout << G4endl;
283  G4cout << "Neutron HP Thermal Scattering Data: Following material-element pairs and/or elements are registered." << G4endl;
284  for ( std::map < std::pair < const G4Material* , const G4Element* > , G4int >::iterator it = dic.begin() ; it != dic.end() ; it++ )
285  {
286  if ( it->first.first != NULL )
287  {
288  G4cout << "Material " << it->first.first->GetName() << " - Element " << it->first.second->GetName() << ", internal thermal scattering id " << it->second << G4endl;
289  }
290  else
291  {
292  G4cout << "Element " << it->first.second->GetName() << ", internal thermal scattering id " << it->second << G4endl;
293  }
294  }
295  G4cout << G4endl;
296 
297 
298  //G4cout << "Neutron HP Thermal Scattering Data: Following NDL thermal scattering files are assigned to the internal thermal scattering id." << G4endl;
299  //for ( std::map < G4String , G4int >::iterator it = co_dic.begin() ; it != co_dic.end() ; it++ )
300  //{
301  // G4cout << "NDL file name " << it->first << ", internal thermal scattering id " << it->second << G4endl;
302  //}
303 
304 
305  // Read Cross Section Data files
306 
307  G4String dirName;
308  if ( !getenv( "G4NEUTRONHPDATA" ) )
309  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
310  G4String baseName = getenv( "G4NEUTRONHPDATA" );
311 
312  dirName = baseName + "/ThermalScattering";
313 
314  G4String ndl_filename;
315  G4String full_name;
316 
317  for ( std::map < G4String , G4int >::iterator it = co_dic.begin() ; it != co_dic.end() ; it++ )
318  {
319  ndl_filename = it->first;
320  G4int ts_ID = it->second;
321 
322  // Coherent
323  full_name = dirName + "/Coherent/CrossSection/" + ndl_filename;
324  std::map< G4double , G4NeutronHPVector* >* coh_amapTemp_EnergyCross = readData( full_name );
325  coherent.insert ( std::pair < G4int , std::map< G4double , G4NeutronHPVector* >* > ( ts_ID , coh_amapTemp_EnergyCross ) );
326 
327  // Incoherent
328  full_name = dirName + "/Incoherent/CrossSection/" + ndl_filename;
329  std::map< G4double , G4NeutronHPVector* >* incoh_amapTemp_EnergyCross = readData( full_name );
330  incoherent.insert ( std::pair < G4int , std::map< G4double , G4NeutronHPVector* >* > ( ts_ID , incoh_amapTemp_EnergyCross ) );
331 
332  // Inelastic
333  full_name = dirName + "/Inelastic/CrossSection/" + ndl_filename;
334  std::map< G4double , G4NeutronHPVector* >* inela_amapTemp_EnergyCross = readData( full_name );
335  inelastic.insert ( std::pair < G4int , std::map< G4double , G4NeutronHPVector* >* > ( ts_ID , inela_amapTemp_EnergyCross ) );
336 
337  }
338 
339 }
340 
341 
342 
343 std::map< G4double , G4NeutronHPVector* >* G4NeutronHPThermalScatteringData::readData ( G4String full_name )
344 {
345 
346  std::map< G4double , G4NeutronHPVector* >* aData = new std::map< G4double , G4NeutronHPVector* >;
347 
348  //std::ifstream theChannel( full_name.c_str() );
349  std::istringstream theChannel;
350  G4NeutronHPManager::GetInstance()->GetDataStream(full_name,theChannel);
351 
352  //G4cout << "G4NeutronHPThermalScatteringData " << name << G4endl;
353 
354  G4int dummy;
355  while ( theChannel >> dummy ) // MF
356  {
357  theChannel >> dummy; // MT
358  G4double temp;
359  theChannel >> temp;
360  G4NeutronHPVector* anEnergyCross = new G4NeutronHPVector;
361  G4int nData;
362  theChannel >> nData;
363  anEnergyCross->Init ( theChannel , nData , eV , barn );
364  aData->insert ( std::pair < G4double , G4NeutronHPVector* > ( temp , anEnergyCross ) );
365  }
366  //theChannel.close();
367 
368  return aData;
369 
370 }
371 
372 
373 
375 {
376  if( &aP != G4Neutron::Neutron() )
377  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
378 // G4cout << "G4NeutronHPThermalScatteringData::DumpPhysicsTable still to be implemented"<<G4endl;
379 }
380 
381 //#include "G4Nucleus.hh"
382 //#include "G4NucleiPropertiesTable.hh"
383 //#include "G4Neutron.hh"
384 //#include "G4Electron.hh"
385 
386 
387 
388 /*
389 G4double G4NeutronHPThermalScatteringData::GetCrossSection( const G4DynamicParticle* aP , const G4Element*anE , G4double aT )
390 {
391 
392  G4double result = 0;
393  const G4Material* aM = NULL;
394 
395  G4int iele = anE->GetIndex();
396 
397  if ( dic.find( std::pair < const G4Material* , const G4Element* > ( (G4Material*)NULL , anE ) ) != dic.end() )
398  {
399  iele = dic.find( std::pair < const G4Material* , const G4Element* > ( (G4Material*)NULL , anE ) )->second;
400  }
401  else if ( dic.find( std::pair < const G4Material* , const G4Element* > ( aM , anE ) ) != dic.end() )
402  {
403  iele = dic.find( std::pair < const G4Material* , const G4Element* > ( aM , anE ) )->second;
404  }
405  else
406  {
407  return result;
408  }
409 
410  G4double Xcoh = GetX ( aP , aT , coherent.find(iele)->second );
411  G4double Xincoh = GetX ( aP , aT , incoherent.find(iele)->second );
412  G4double Xinela = GetX ( aP , aT , inelastic.find(iele)->second );
413 
414  result = Xcoh + Xincoh + Xinela;
415 
416  //G4cout << "G4NeutronHPThermalScatteringData::GetCrossSection Tot= " << result/barn << " Coherent= " << Xcoh/barn << " Incoherent= " << Xincoh/barn << " Inelastic= " << Xinela/barn << G4endl;
417 
418  return result;
419 
420 }
421 */
422 
424 {
425  G4double result = 0;
426 
427  G4int ts_id =getTS_ID( aM , anE );
428 
429  if ( ts_id == -1 ) return result;
430 
431  G4double aT = aM->GetTemperature();
432 
433  G4double Xcoh = GetX ( aP , aT , coherent.find(ts_id)->second );
434  G4double Xincoh = GetX ( aP , aT , incoherent.find(ts_id)->second );
435  G4double Xinela = GetX ( aP , aT , inelastic.find(ts_id)->second );
436 
437  result = Xcoh + Xincoh + Xinela;
438 
439  //G4cout << "G4NeutronHPThermalScatteringData::GetCrossSection Tot= " << result/barn << " Coherent= " << Xcoh/barn << " Incoherent= " << Xincoh/barn << " Inelastic= " << Xinela/barn << G4endl;
440 
441  return result;
442 }
443 
444 
446 {
447  G4double result = 0;
448  G4int ts_id = getTS_ID( aM , anE );
449  G4double aT = aM->GetTemperature();
450  result = GetX ( aP , aT , inelastic.find( ts_id )->second );
451  return result;
452 }
453 
455 {
456  G4double result = 0;
457  G4int ts_id = getTS_ID( aM , anE );
458  G4double aT = aM->GetTemperature();
459  result = GetX ( aP , aT , coherent.find( ts_id )->second );
460  return result;
461 }
462 
464 {
465  G4double result = 0;
466  G4int ts_id = getTS_ID( aM , anE );
467  G4double aT = aM->GetTemperature();
468  result = GetX ( aP , aT , incoherent.find( ts_id )->second );
469  return result;
470 }
471 
472 
473 
474 G4int G4NeutronHPThermalScatteringData::getTS_ID ( const G4Material* material , const G4Element* element )
475 {
476  G4int result = -1;
477  if ( dic.find( std::pair < const G4Material* , const G4Element* > ( (G4Material*)NULL , element ) ) != dic.end() )
478  return dic.find( std::pair < const G4Material* , const G4Element* > ( (G4Material*)NULL , element ) )->second;
479  if ( dic.find( std::pair < const G4Material* , const G4Element* > ( material , element ) ) != dic.end() )
480  return dic.find( std::pair < const G4Material* , const G4Element* > ( material , element ) )->second;
481  return result;
482 }
483 
484 
485 
486 
487 G4double G4NeutronHPThermalScatteringData::GetX ( const G4DynamicParticle* aP, G4double aT , std::map < G4double , G4NeutronHPVector* >* amapTemp_EnergyCross )
488 {
489 
490  G4double result = 0;
491  if ( amapTemp_EnergyCross->size() == 0 ) return result;
492 
493 
494  G4double eKinetic = aP->GetKineticEnergy();
495 
496  if ( amapTemp_EnergyCross->size() == 1 ) {
497  if ( std::fabs ( aT - amapTemp_EnergyCross->begin()->first ) / amapTemp_EnergyCross->begin()->first > 0.1 ) {
498  G4cout << "G4NeutronHPThermalScatteringData:: The temperature of material ("
499  << aT/kelvin << "K) is different more than 10% from temperature of thermal scattering file expected ("
500  << amapTemp_EnergyCross->begin()->first << "K). Result may not be reliable."
501  << G4endl;
502  }
503  result = amapTemp_EnergyCross->begin()->second->GetXsec ( eKinetic );
504  return result;
505  }
506 
507  std::map< G4double , G4NeutronHPVector* >::iterator it;
508  for ( it = amapTemp_EnergyCross->begin() ; it != amapTemp_EnergyCross->end() ; it++ ) {
509  if ( aT < it->first ) break;
510  }
511  if ( it == amapTemp_EnergyCross->begin() && it != amapTemp_EnergyCross->end() ) it++; // lower than the first
512  if ( it != amapTemp_EnergyCross->begin() && it == amapTemp_EnergyCross->end() ) it--; // upper than the last
513 
514  G4double TH = it->first;
515  G4double XH = it->second->GetXsec ( eKinetic );
516 
517  //G4cout << "G4NeutronHPThermalScatteringData::GetX TH " << TH << " E " << eKinetic << " XH " << XH << G4endl;
518 
519  if ( it != amapTemp_EnergyCross->begin() ) it--;
520  G4double TL = it->first;
521  G4double XL = it->second->GetXsec ( eKinetic );
522 
523  //G4cout << "G4NeutronHPThermalScatteringData::GetX TL " << TL << " E " << eKinetic << " XL " << XL << G4endl;
524 
525  if ( TH == TL )
526  throw G4HadronicException(__FILE__, __LINE__, "Thermal Scattering Data Error!");
527 
528  G4double T = aT;
529  G4double X = ( XH - XL ) / ( TH - TL ) * ( T - TL ) + XL;
530  result = X;
531 
532  return result;
533 }
534 
535 
537 {
538  names->AddThermalElement( nameG4Element , filename );
539 }
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
G4double GetKineticEnergy() const
G4int first(char) const
static G4NeutronHPManager * GetInstance()
G4double GetCoherentCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
const G4String & GetName() const
Definition: G4Material.hh:176
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
std::vector< G4Material * > G4MaterialTable
G4bool IsApplicable(const G4DynamicParticle *, const G4Element *)
void Init(std::istream &aDataFile, G4int total, G4double ux=1., G4double uy=1.)
G4ParticleDefinition * GetDefinition() const
void GetDataStream(G4String, std::istringstream &iss)
#define G4ThreadLocal
Definition: tls.hh:52
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
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 size_t GetNumberOfMaterials()
Definition: G4Material.cc:571
G4double GetIncoherentCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
void SetMaxKinEnergy(G4double value)
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double GetInelasticCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
G4String GetTS_NDL_Name(G4String nameG4Element)
G4double GetTemperature() const
Definition: G4Material.hh:180
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395