G4NeutronElasticXS Class Reference

#include <G4NeutronElasticXS.hh>

Inheritance diagram for G4NeutronElasticXS:

G4VCrossSectionDataSet

Public Member Functions

 G4NeutronElasticXS ()
virtual ~G4NeutronElasticXS ()
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *)
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
virtual void CrossSectionDescription (std::ostream &) const

Detailed Description

Definition at line 59 of file G4NeutronElasticXS.hh.


Constructor & Destructor Documentation

G4NeutronElasticXS::G4NeutronElasticXS (  ) 

Definition at line 59 of file G4NeutronElasticXS.cc.

References G4cout, G4endl, and G4VCrossSectionDataSet::verboseLevel.

00060  : G4VCrossSectionDataSet("G4NeutronElasticXS"),
00061    proton(G4Proton::Proton()), maxZ(92)
00062 {
00063   //  verboseLevel = 0;
00064   if(verboseLevel > 0){
00065     G4cout  << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < " 
00066             << maxZ + 1 << G4endl;
00067   }
00068   data.resize(maxZ+1, 0);
00069   coeff.resize(maxZ+1, 1.0);
00070   ggXsection = new G4GlauberGribovCrossSection();
00071   fNucleon = new G4HadronNucleonXsc();
00072   isInitialized = false;
00073 }

G4NeutronElasticXS::~G4NeutronElasticXS (  )  [virtual]

Definition at line 75 of file G4NeutronElasticXS.cc.

00076 {
00077   delete fNucleon;
00078   for(G4int i=0; i<=maxZ; ++i) {
00079     delete data[i];
00080   }
00081 }


Member Function Documentation

void G4NeutronElasticXS::BuildPhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 143 of file G4NeutronElasticXS.cc.

References G4cout, G4endl, G4Element::GetElementTable(), G4Element::GetNumberOfElements(), G4ParticleDefinition::GetParticleName(), G4Neutron::Neutron(), and G4VCrossSectionDataSet::verboseLevel.

00144 {
00145   if(isInitialized) { return; }
00146   if(verboseLevel > 0){
00147     G4cout << "G4NeutronElasticXS::BuildPhysicsTable for " 
00148            << p.GetParticleName() << G4endl;
00149   }
00150   if(p.GetParticleName() != "neutron") { 
00151     throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
00152     return; 
00153   }
00154   isInitialized = true;
00155 
00156   // check environment variable 
00157   // Build the complete string identifying the file with the data set
00158   char* path = getenv("G4NEUTRONXSDATA");
00159   if (!path){
00160     throw G4HadronicException(__FILE__, __LINE__, 
00161                               "G4NEUTRONXSDATA environment variable not defined");
00162     return;
00163   }
00164 
00165   G4DynamicParticle* dynParticle = 
00166     new G4DynamicParticle(G4Neutron::Neutron(),G4ThreeVector(1,0,0),1);
00167 
00168   // Access to elements
00169   const G4ElementTable* theElmTable = G4Element::GetElementTable();
00170   size_t numOfElm = G4Element::GetNumberOfElements();
00171   if(numOfElm > 0) {
00172     for(size_t i=0; i<numOfElm; ++i) {
00173       G4int Z = G4int(((*theElmTable)[i])->GetZ());
00174       if(Z < 1)         { Z = 1; }
00175       else if(Z > maxZ) { Z = maxZ; }
00176       //G4cout << "Z= " << Z << G4endl;
00177       // Initialisation 
00178       if(!data[Z]) { Initialise(Z, dynParticle, path); }
00179     }
00180   }
00181   delete dynParticle;
00182 }

void G4NeutronElasticXS::CrossSectionDescription ( std::ostream &   )  const [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 83 of file G4NeutronElasticXS.cc.

00084 {
00085   outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
00086           << "cross section on nuclei using data from the high precision\n"
00087           << "neutron database.  These data are simplified and smoothed over\n"
00088           << "the resonance region in order to reduce CPU time.\n"
00089           << "G4NeutronElasticXS is valid for energies up to 20 MeV, for all\n"
00090           << "targets through U.\n";  
00091 }

G4double G4NeutronElasticXS::GetElementCrossSection ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = 0 
) [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 101 of file G4NeutronElasticXS.cc.

References G4PhysicsVector::Energy(), G4cout, G4endl, G4GlauberGribovCrossSection::GetElasticGlauberGribovXsc(), G4HadronNucleonXsc::GetElasticHadronNucleonXsc(), G4HadronNucleonXsc::GetHadronNucleonXscPDG(), G4GlauberGribovCrossSection::GetIsoCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4PhysicsVector::GetVectorLength(), G4NistManager::Instance(), CLHEP::detail::n, G4PhysicsVector::Value(), and G4VCrossSectionDataSet::verboseLevel.

00103 {
00104   G4double xs = 0.0;
00105   G4double ekin = aParticle->GetKineticEnergy();
00106 
00107   if(Z < 1 || Z > maxZ) { return xs; }
00108 
00109   G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
00110   G4PhysicsVector* pv = data[Z];
00111   //  G4cout  << "G4NeutronElasticXS::GetCrossSection e= " << ekin << " Z= " << Z << G4endl;
00112 
00113   // element was not initialised
00114   if(!pv) {
00115     Initialise(Z);
00116     pv = data[Z];
00117     if(!pv) { return xs; }
00118   }
00119 
00120   G4double e1 = pv->Energy(0);
00121   if(ekin <= e1) { return (*pv)[0]; }
00122 
00123   G4int n = pv->GetVectorLength() - 1;
00124   G4double e2 = pv->Energy(n);
00125 
00126   if(ekin <= e2) { 
00127     xs = pv->Value(ekin); 
00128   } else if(1 == Z) { 
00129     fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
00130     xs = coeff[1]*fNucleon->GetElasticHadronNucleonXsc();
00131   } else {          
00132     ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
00133     xs = coeff[Z]*ggXsection->GetElasticGlauberGribovXsc();
00134   }
00135 
00136   if(verboseLevel > 0){
00137     G4cout  << "ekin= " << ekin << ",  XSinel= " << xs << G4endl;
00138   }
00139   return xs;
00140 }

G4bool G4NeutronElasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
) [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 94 of file G4NeutronElasticXS.cc.

00096 {
00097   return true;
00098 }


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