Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4NeutronHPThermalScattering Class Reference

#include <G4NeutronHPThermalScattering.hh>

Inheritance diagram for G4NeutronHPThermalScattering:
G4HadronicInteraction

Public Member Functions

 G4NeutronHPThermalScattering ()
 
 ~G4NeutronHPThermalScattering ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
void AddUserThermalScatteringFile (G4String, G4String)
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 64 of file G4NeutronHPThermalScattering.hh.

Constructor & Destructor Documentation

G4NeutronHPThermalScattering::G4NeutronHPThermalScattering ( )

Definition at line 51 of file G4NeutronHPThermalScattering.cc.

References G4NeutronHPThermalScatteringData::BuildPhysicsTable(), python.hepunit::eV, G4Material::GetMaterialTable(), G4Neutron::Neutron(), G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

52  :G4HadronicInteraction("NeutronHPThermalScattering")
53 {
54 
55  theHPElastic = new G4NeutronHPElastic();
56 
57  SetMinEnergy( 0.*eV );
58  SetMaxEnergy( 4*eV );
59  theXSection = new G4NeutronHPThermalScatteringData();
60  theXSection->BuildPhysicsTable( *(G4Neutron::Neutron()) );
61 
62  sizeOfMaterialTable = G4Material::GetMaterialTable()->size();
63  buildPhysicsTable();
64 
65 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void BuildPhysicsTable(const G4ParticleDefinition &)
void SetMaxEnergy(const G4double anEnergy)
G4NeutronHPThermalScattering::~G4NeutronHPThermalScattering ( )

Definition at line 69 of file G4NeutronHPThermalScattering.cc.

70 {
71 
72  for ( std::map < G4int , std::map < G4double , std::vector < E_isoAng* >* >* >::iterator it = incoherentFSs.begin() ; it != incoherentFSs.end() ; it++ )
73  {
74  std::map < G4double , std::vector < E_isoAng* >* >::iterator itt;
75  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
76  {
77  std::vector< E_isoAng* >::iterator ittt;
78  for ( ittt = itt->second->begin(); ittt != itt->second->end() ; ittt++ )
79  {
80  delete *ittt;
81  }
82  delete itt->second;
83  }
84  delete it->second;
85  }
86 
87  for ( std::map < G4int , std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >* >::iterator it = coherentFSs.begin() ; it != coherentFSs.end() ; it++ )
88  {
89  std::map < G4double , std::vector < std::pair< G4double , G4double >* >* >::iterator itt;
90  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
91  {
92  std::vector < std::pair< G4double , G4double >* >::iterator ittt;
93  for ( ittt = itt->second->begin(); ittt != itt->second->end() ; ittt++ )
94  {
95  delete *ittt;
96  }
97  delete itt->second;
98  }
99  delete it->second;
100  }
101 
102  for ( std::map < G4int , std::map < G4double , std::vector < E_P_E_isoAng* >* >* >::iterator it = inelasticFSs.begin() ; it != inelasticFSs.end() ; it++ )
103  {
104  std::map < G4double , std::vector < E_P_E_isoAng* >* >::iterator itt;
105  for ( itt = it->second->begin() ; itt != it->second->end() ; itt++ )
106  {
107  std::vector < E_P_E_isoAng* >::iterator ittt;
108  for ( ittt = itt->second->begin(); ittt != itt->second->end() ; ittt++ )
109  {
110  std::vector < E_isoAng* >::iterator it4;
111  for ( it4 = (*ittt)->vE_isoAngle.begin() ; it4 != (*ittt)->vE_isoAngle.end() ; it4++ )
112  {
113  delete *it4;
114  }
115  delete *ittt;
116  }
117  delete itt->second;
118  }
119  delete it->second;
120  }
121 
122  delete theHPElastic;
123  delete theXSection;
124 }
int G4int
Definition: G4Types.hh:78
double G4double
Definition: G4Types.hh:76

Member Function Documentation

void G4NeutronHPThermalScattering::AddUserThermalScatteringFile ( G4String  nameG4Element,
G4String  filename 
)

Definition at line 1037 of file G4NeutronHPThermalScattering.cc.

References G4NeutronHPThermalScatteringNames::AddThermalElement(), and G4NeutronHPThermalScatteringData::AddUserThermalScatteringFile().

1038 {
1039  names.AddThermalElement( nameG4Element , filename );
1040  theXSection->AddUserThermalScatteringFile( nameG4Element , filename );
1041  buildPhysicsTable();
1042 }
G4HadFinalState * G4NeutronHPThermalScattering::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 311 of file G4NeutronHPThermalScattering.cc.

References G4NeutronHPThermalScatteringData::BuildPhysicsTable(), E_isoAng::energy, python.hepunit::eV, G4UniformRand, G4HadProjectile::Get4Momentum(), G4NeutronHPThermalScatteringData::GetCoherentCrossSection(), G4NeutronHPThermalScatteringData::GetCrossSection(), G4HadProjectile::GetDefinition(), G4Material::GetElement(), G4NeutronHPThermalScatteringData::GetInelasticCrossSection(), G4HadProjectile::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4Material::GetMaterialTable(), G4Material::GetNumberOfElements(), G4Material::GetTemperature(), G4Element::GetZ(), G4Nucleus::GetZ_asInt(), E_isoAng::isoAngle, python.hepunit::kelvin, E_isoAng::n, n, python.hepunit::second, G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetMomentumChange(), and G4HadronicInteraction::theParticleChange.

312 {
313 
314  //Trick for dynamically generated materials
315  if ( sizeOfMaterialTable != G4Material::GetMaterialTable()->size() ) {
316  sizeOfMaterialTable = G4Material::GetMaterialTable()->size();
317  buildPhysicsTable();
318  theXSection->BuildPhysicsTable( *aTrack.GetDefinition() );
319  }
320 // Select Element > Reaction >
321 
322  const G4Material * theMaterial = aTrack.GetMaterial();
323  G4double aTemp = theMaterial->GetTemperature();
324  G4int n = theMaterial->GetNumberOfElements();
325  //static const G4ElementTable* theElementTable = G4Element::GetElementTable();
326 
327  G4bool findThermalElement = false;
328  G4int ielement;
329  const G4Element* theElement = NULL;
330  for ( G4int i = 0; i < n ; i++ )
331  {
332  theElement = theMaterial->GetElement(i);
333  //Select target element
334  if ( aNucleus.GetZ_asInt() == (G4int)(theElement->GetZ() + 0.5 ) )
335  {
336  //Check Applicability of Thermal Scattering
337  if ( getTS_ID( NULL , theElement ) != -1 )
338  {
339  ielement = getTS_ID( NULL , theElement );
340  findThermalElement = true;
341  break;
342  }
343  else if ( getTS_ID( theMaterial , theElement ) != -1 )
344  {
345  ielement = getTS_ID( theMaterial , theElement );
346  findThermalElement = true;
347  break;
348  }
349  }
350  }
351 
352  if ( findThermalElement == true )
353  {
354 
355 // Select Reaction (Inelastic, coherent, incoherent)
356 
357  G4ParticleDefinition* pd = const_cast< G4ParticleDefinition* >( aTrack.GetDefinition() );
358  G4DynamicParticle* dp = new G4DynamicParticle ( pd , aTrack.Get4Momentum() );
359  G4double total = theXSection->GetCrossSection( dp , theElement , theMaterial );
360  G4double inelastic = theXSection->GetInelasticCrossSection( dp , theElement , theMaterial );
361 
362 
363  G4double random = G4UniformRand();
364  if ( random <= inelastic/total )
365  {
366  // Inelastic
367 
368  // T_L and T_H
369  std::map < G4double , std::vector< E_P_E_isoAng* >* >::iterator it;
370  std::vector<G4double> v_temp;
371  v_temp.clear();
372  for ( it = inelasticFSs.find( ielement )->second->begin() ; it != inelasticFSs.find( ielement )->second->end() ; it++ )
373  {
374  v_temp.push_back( it->first );
375  }
376 
377 // T_L T_H
378  std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
379 //
380 // For T_L aNEP_EPM_TL and T_H aNEP_EPM_TH
381 //
382  std::vector< E_P_E_isoAng* >* vNEP_EPM_TL = 0;
383  std::vector< E_P_E_isoAng* >* vNEP_EPM_TH = 0;
384 
385  if ( tempLH.first != 0.0 && tempLH.second != 0.0 )
386  {
387  vNEP_EPM_TL = inelasticFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second;
388  vNEP_EPM_TH = inelasticFSs.find( ielement )->second->find ( tempLH.second/kelvin )->second;
389  }
390  else if ( tempLH.first == 0.0 )
391  {
392  std::map < G4double , std::vector< E_P_E_isoAng* >* >::iterator itm;
393  itm = inelasticFSs.find( ielement )->second->begin();
394  vNEP_EPM_TL = itm->second;
395  itm++;
396  vNEP_EPM_TH = itm->second;
397  tempLH.first = tempLH.second;
398  tempLH.second = itm->first;
399  }
400  else if ( tempLH.second == 0.0 )
401  {
402  std::map < G4double , std::vector< E_P_E_isoAng* >* >::iterator itm;
403  itm = inelasticFSs.find( ielement )->second->end();
404  itm--;
405  vNEP_EPM_TH = itm->second;
406  itm--;
407  vNEP_EPM_TL = itm->second;
408  tempLH.second = tempLH.first;
409  tempLH.first = itm->first;
410  }
411 
412  G4double rand_for_sE = G4UniformRand();
413 
414  std::pair< G4double , E_isoAng > TL = create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( rand_for_sE , aTrack.GetKineticEnergy() , vNEP_EPM_TL );
415  std::pair< G4double , E_isoAng > TH = create_sE_and_EPM_from_pE_and_vE_P_E_isoAng ( rand_for_sE , aTrack.GetKineticEnergy() , vNEP_EPM_TH );
416 
417  G4double sE;
418  sE = get_linear_interpolated ( aTemp , std::pair < G4double , G4double > ( tempLH.first , TL.first ) , std::pair < G4double , G4double > ( tempLH.second , TH.first ) );
419 
420  G4double mu=1.0;
421  E_isoAng anE_isoAng;
422  if ( TL.second.n == TH.second.n )
423  {
424  anE_isoAng.energy = sE;
425  anE_isoAng.n = TL.second.n;
426  for ( G4int i=0 ; i < anE_isoAng.n ; i++ )
427  {
428  G4double angle;
429  angle = get_linear_interpolated ( aTemp , std::pair< G4double , G4double > ( tempLH.first , TL.second.isoAngle[ i ] ) , std::pair< G4double , G4double > ( tempLH.second , TH.second.isoAngle[ i ] ) );
430  anE_isoAng.isoAngle.push_back( angle );
431  }
432  mu = getMu( &anE_isoAng );
433 
434  } else {
435  //TL.second.n != TH.second.n
436  G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data! Do not yet supported");
437  }
438 
439  //set
441  theParticleChange.SetMomentumChange( 0.0 , std::sqrt ( 1 - mu*mu ) , mu );
442 
443  }
444  //else if ( random <= ( inelastic + theXSection->GetCoherentCrossSection( dp , (*theElementTable)[ ielement ] , aTemp ) ) / total )
445  else if ( random <= ( inelastic + theXSection->GetCoherentCrossSection( dp , theElement , theMaterial ) ) / total )
446  {
447  // Coherent Elastic
448 
449  G4double E = aTrack.GetKineticEnergy();
450 
451  // T_L and T_H
452  std::map < G4double , std::vector< std::pair< G4double , G4double >* >* >::iterator it;
453  std::vector<G4double> v_temp;
454  v_temp.clear();
455  for ( it = coherentFSs.find( ielement )->second->begin() ; it != coherentFSs.find( ielement )->second->end() ; it++ )
456  {
457  v_temp.push_back( it->first );
458  }
459 
460 // T_L T_H
461  std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
462 //
463 //
464 // For T_L anEPM_TL and T_H anEPM_TH
465 //
466  std::vector< std::pair< G4double , G4double >* >* pvE_p_TL = NULL;
467  std::vector< std::pair< G4double , G4double >* >* pvE_p_TH = NULL;
468 
469  if ( tempLH.first != 0.0 && tempLH.second != 0.0 )
470  {
471  pvE_p_TL = coherentFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second;
472  pvE_p_TH = coherentFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second;
473  }
474  else if ( tempLH.first == 0.0 )
475  {
476  pvE_p_TL = coherentFSs.find( ielement )->second->find ( v_temp[ 0 ] )->second;
477  pvE_p_TH = coherentFSs.find( ielement )->second->find ( v_temp[ 1 ] )->second;
478  tempLH.first = tempLH.second;
479  tempLH.second = v_temp[ 1 ];
480  }
481  else if ( tempLH.second == 0.0 )
482  {
483  pvE_p_TL = coherentFSs.find( ielement )->second->find ( v_temp.back() )->second;
484  std::vector< G4double >::iterator itv;
485  itv = v_temp.end();
486  itv--;
487  itv--;
488  pvE_p_TL = coherentFSs.find( ielement )->second->find ( *itv )->second;
489  tempLH.second = tempLH.first;
490  tempLH.first = *itv;
491  }
492  else
493  {
494  //tempLH.first == 0.0 && tempLH.second
495  G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data! Unexpected temperature values in data");
496  }
497 
498  std::vector< G4double > vE_T;
499  std::vector< G4double > vp_T;
500 
501  G4int n1 = pvE_p_TL->size();
502  //G4int n2 = pvE_p_TH->size();
503 
504  for ( G4int i=1 ; i < n1 ; i++ )
505  {
506  if ( (*pvE_p_TL)[i]->first != (*pvE_p_TH)[i]->first ) G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data!");
507  vE_T.push_back ( (*pvE_p_TL)[i]->first );
508  vp_T.push_back ( get_linear_interpolated ( aTemp , std::pair< G4double , G4double > ( tempLH.first , (*pvE_p_TL)[i]->second ) , std::pair< G4double , G4double > ( tempLH.second , (*pvE_p_TL)[i]->second ) ) );
509  }
510 
511  G4int j = 0;
512  for ( G4int i = 1 ; i < n ; i++ )
513  {
514  if ( E/eV < vE_T[ i ] )
515  {
516  j = i-1;
517  break;
518  }
519  }
520 
521  G4double rand_for_mu = G4UniformRand();
522 
523  G4int k = 0;
524  for ( G4int i = 1 ; i < j ; i++ )
525  {
526  G4double Pi = vp_T[ i ] / vp_T[ j ];
527  if ( rand_for_mu < Pi )
528  {
529  k = i-1;
530  break;
531  }
532  }
533 
534  //G4double Ei = vE_T[ j ];
535  G4double Ei = vE_T[ k ];
536 
537  G4double mu = 1 - 2 * Ei / (E/eV) ;
538  //111102
539  if ( mu < -1.0 ) mu = -1.0;
540 
542  theParticleChange.SetMomentumChange( 0.0 , std::sqrt ( 1 - mu*mu ) , mu );
543 
544 
545  }
546  else
547  {
548  // InCoherent Elastic
549 
550  // T_L and T_H
551  std::map < G4double , std::vector < E_isoAng* >* >::iterator it;
552  std::vector<G4double> v_temp;
553  v_temp.clear();
554  for ( it = incoherentFSs.find( ielement )->second->begin() ; it != incoherentFSs.find( ielement )->second->end() ; it++ )
555  {
556  v_temp.push_back( it->first );
557  }
558 
559 // T_L T_H
560  std::pair < G4double , G4double > tempLH = find_LH ( aTemp , &v_temp );
561 
562 //
563 // For T_L anEPM_TL and T_H anEPM_TH
564 //
565 
566  E_isoAng anEPM_TL_E;
567  E_isoAng anEPM_TH_E;
568 
569  if ( tempLH.first != 0.0 && tempLH.second != 0.0 )
570  {
571  anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( tempLH.first/kelvin )->second );
572  anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( tempLH.second/kelvin )->second );
573  }
574  else if ( tempLH.first == 0.0 )
575  {
576  anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( v_temp[ 0 ] )->second );
577  anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( v_temp[ 1 ] )->second );
578  tempLH.first = tempLH.second;
579  tempLH.second = v_temp[ 1 ];
580  }
581  else if ( tempLH.second == 0.0 )
582  {
583  anEPM_TH_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( v_temp.back() )->second );
584  std::vector< G4double >::iterator itv;
585  itv = v_temp.end();
586  itv--;
587  itv--;
588  anEPM_TL_E = create_E_isoAng_from_energy ( aTrack.GetKineticEnergy() , incoherentFSs.find( ielement )->second->find ( *itv )->second );
589  tempLH.second = tempLH.first;
590  tempLH.first = *itv;
591  }
592 
593  // E_isoAng for aTemp and aTrack.GetKineticEnergy()
594  G4double mu=1.0;
595  E_isoAng anEPM_T_E;
596 
597  if ( anEPM_TL_E.n == anEPM_TH_E.n )
598  {
599  anEPM_T_E.n = anEPM_TL_E.n;
600  for ( G4int i=0 ; i < anEPM_TL_E.n ; i++ )
601  {
602  G4double angle;
603  angle = get_linear_interpolated ( aTemp , std::pair< G4double , G4double > ( tempLH.first , anEPM_TL_E.isoAngle[ i ] ) , std::pair< G4double , G4double > ( tempLH.second , anEPM_TH_E.isoAngle[ i ] ) );
604  anEPM_T_E.isoAngle.push_back( angle );
605  }
606  mu = getMu ( &anEPM_T_E );
607 
608  } else {
609  // anEPM_TL_E.n != anEPM_TH_E.n
610  G4HadronicException(__FILE__, __LINE__, "A problem is found in Thermal Scattering Data! Do not yet supported");
611  }
612 
613  // Set Final State
614  theParticleChange.SetEnergyChange( aTrack.GetKineticEnergy() ); // No energy change in Elastic
615  theParticleChange.SetMomentumChange( 0.0 , std::sqrt ( 1 - mu*mu ) , mu );
616 
617  }
618  delete dp;
619 
620  return &theParticleChange;
621 
622  }
623  else
624  {
625  // Not thermal element
626  // Neutron HP will handle
627  return theHPElastic -> ApplyYourself( aTrack, aNucleus );
628  }
629 
630 }
std::vector< G4double > isoAngle
G4double GetCoherentCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
G4double GetZ() const
Definition: G4Element.hh:131
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
const G4int n
const G4LorentzVector & Get4Momentum() const
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
void SetEnergyChange(G4double anEnergy)
G4double total(Particle const *const p1, Particle const *const p2)
G4double GetInelasticCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *)
G4double GetTemperature() const
Definition: G4Material.hh:180
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
void SetMomentumChange(const G4ThreeVector &aV)
const std::pair< G4double, G4double > G4NeutronHPThermalScattering::GetFatalEnergyCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 1031 of file G4NeutronHPThermalScattering.cc.

References DBL_MAX, and python.hepunit::perCent.

1032 {
1033  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
1034  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
1035 }
float perCent
Definition: hepunit.py:239
#define DBL_MAX
Definition: templates.hh:83

The documentation for this class was generated from the following files: