G4NeutronHPorLElastic Class Reference

#include <G4NeutronHPorLElastic.hh>

Inheritance diagram for G4NeutronHPorLElastic:

G4HadronicInteraction

Public Member Functions

 G4NeutronHPorLElastic ()
 ~G4NeutronHPorLElastic ()
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
virtual const std::pair< G4double,
G4double
GetFatalEnergyCheckLevels () const
G4int GetNiso ()
void DoNotSuspend ()
G4bool IsThisElementOK (G4String)
G4VCrossSectionDataSetGiveXSectionDataSet ()

Detailed Description

Definition at line 56 of file G4NeutronHPorLElastic.hh.


Constructor & Destructor Documentation

G4NeutronHPorLElastic::G4NeutronHPorLElastic (  ) 

Definition at line 42 of file G4NeutronHPorLElastic.cc.

References G4cout, G4endl, G4Element::GetElementTable(), G4Element::GetNumberOfElements(), G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

00043   :G4HadronicInteraction("NeutronHPorLElastic")
00044 {
00045    overrideSuspension = false;
00046    G4NeutronHPElasticFS * theFS = new G4NeutronHPElasticFS;
00047    if(!getenv("G4NEUTRONHPDATA")) 
00048        throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
00049    dirName = getenv("G4NEUTRONHPDATA");
00050    G4String tString = "/Elastic/";
00051    dirName = dirName + tString;
00052 //    G4cout <<"G4NeutronHPorLElastic::G4NeutronHPorLElastic testit "<<dirName<<G4endl;
00053    numEle = G4Element::GetNumberOfElements();
00054    theElastic = new G4NeutronHPChannel[numEle];
00055    unavailable_elements.clear();
00056    for (G4int i=0; i<numEle; i++)
00057    {
00058       theElastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
00059       //while(!theElastic[i].Register(theFS));
00060       try { while(!theElastic[i].Register(theFS)) ; }
00061       catch ( G4HadronicException )
00062       {
00063           unavailable_elements.insert ( (*(G4Element::GetElementTable()))[i]->GetName() ); 
00064       }
00065    }
00066    delete theFS;
00067    SetMinEnergy(0.*eV);
00068    SetMaxEnergy(20.*MeV);
00069    if ( unavailable_elements.size() > 0 ) 
00070    {
00071       std::set< G4String>::iterator it;
00072       G4cout << "HP Elastic data are not available for thess  elements "<< G4endl;
00073       for ( it = unavailable_elements.begin() ; it != unavailable_elements.end() ; it++ )
00074          G4cout << *it << G4endl;
00075       G4cout << "Low Energy Parameterization Models will be used."<< G4endl;
00076    }
00077 
00078    createXSectionDataSet();
00079 }

G4NeutronHPorLElastic::~G4NeutronHPorLElastic (  ) 

Definition at line 81 of file G4NeutronHPorLElastic.cc.

00082 {
00083    delete [] theElastic;
00084    delete theDataSet; 
00085 }


Member Function Documentation

G4HadFinalState * G4NeutronHPorLElastic::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
) [virtual]

Implements G4HadronicInteraction.

Definition at line 89 of file G4NeutronHPorLElastic.cc.

References G4NeutronHPChannel::ApplyYourself(), G4UniformRand, G4Material::GetElement(), G4Element::GetIndex(), G4HadProjectile::GetMaterial(), G4Material::GetNumberOfElements(), G4Material::GetTemperature(), G4NeutronHPThermalBoost::GetThermalEnergy(), G4Material::GetVecNbOfAtomsPerVolume(), G4NeutronHPChannel::GetXsec(), isAlive, CLHEP::detail::n, and G4HadFinalState::SetStatusChange().

Referenced by G4NeutronHPorLElasticModel::ApplyYourself().

00090 {
00091    const G4Material * theMaterial = aTrack.GetMaterial();
00092    G4int n = theMaterial->GetNumberOfElements();
00093    G4int index = theMaterial->GetElement(0)->GetIndex();
00094    if(n!=1)
00095    {
00096       G4int i;
00097       xSec = new G4double[n];
00098       G4double sum=0;
00099       const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
00100       G4double rWeight;    
00101       G4NeutronHPThermalBoost aThermalE;
00102       for (i=0; i<n; i++)
00103       {
00104         index = theMaterial->GetElement(i)->GetIndex();
00105         rWeight = NumAtomsPerVolume[i];
00106         G4double x = aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i), theMaterial->GetTemperature());
00107 
00108         //xSec[i] = theElastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
00109         //                                                                   theMaterial->GetElement(i),
00110         //                                                                   theMaterial->GetTemperature()));
00111         xSec[i] = theElastic[index].GetXsec(x);
00112 
00113         xSec[i] *= rWeight;
00114         sum+=xSec[i];
00115       }
00116       G4double random = G4UniformRand();
00117       G4double running = 0;
00118       for (i=0; i<n; i++)
00119       {
00120         running += xSec[i];
00121         index = theMaterial->GetElement(i)->GetIndex();
00122         if(random<=running/sum) break;
00123       }
00124       delete [] xSec;
00125       // it is element-wise initialised.
00126     }
00127     G4HadFinalState* finalState = theElastic[index].ApplyYourself(aTrack);
00128     if (overrideSuspension) finalState->SetStatusChange(isAlive);
00129     return finalState;
00130 }

void G4NeutronHPorLElastic::DoNotSuspend (  )  [inline]

Definition at line 70 of file G4NeutronHPorLElastic.hh.

00070 {overrideSuspension = true;}

const std::pair< G4double, G4double > G4NeutronHPorLElastic::GetFatalEnergyCheckLevels (  )  const [virtual]

Reimplemented from G4HadronicInteraction.

Definition at line 149 of file G4NeutronHPorLElastic.cc.

References DBL_MAX.

00150 {
00151    //return std::pair<G4double, G4double>(10*perCent,10*GeV);
00152    return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
00153 }

G4int G4NeutronHPorLElastic::GetNiso (  )  [inline]

Definition at line 68 of file G4NeutronHPorLElastic.hh.

References G4NeutronHPChannel::GetNiso().

00068 {return theElastic[0].GetNiso();}

G4VCrossSectionDataSet* G4NeutronHPorLElastic::GiveXSectionDataSet (  )  [inline]

Definition at line 86 of file G4NeutronHPorLElastic.hh.

Referenced by G4NeutronHPorLElasticModel::GiveHPXSectionDataSet().

00086 { return theDataSet; }; 

G4bool G4NeutronHPorLElastic::IsThisElementOK ( G4String   ) 

Definition at line 134 of file G4NeutronHPorLElastic.cc.

References FALSE, and TRUE.

Referenced by G4NeutronHPorLElasticModel::ApplyYourself().

00135 {
00136    if ( unavailable_elements.find( name ) == unavailable_elements.end() ) 
00137       return TRUE;
00138    else 
00139       return FALSE; 
00140 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:40 2013 for Geant4 by  doxygen 1.4.7