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