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