Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
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
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 59 of file G4NeutronElasticXS.hh.

Constructor & Destructor Documentation

G4NeutronElasticXS::G4NeutronElasticXS ( )

Definition at line 58 of file G4NeutronElasticXS.cc.

References G4cout, G4endl, and G4VCrossSectionDataSet::verboseLevel.

59  : G4VCrossSectionDataSet("G4NeutronElasticXS"),
60  proton(G4Proton::Proton()), maxZ(92)
61 {
62  // verboseLevel = 0;
63  if(verboseLevel > 0){
64  G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
65  << maxZ + 1 << G4endl;
66  }
67  data.resize(maxZ+1, 0);
68  coeff.resize(maxZ+1, 1.0);
69  ggXsection = new G4GlauberGribovCrossSection();
70  fNucleon = new G4HadronNucleonXsc();
71  isInitialized = false;
72 }
G4VCrossSectionDataSet(const G4String &nam="")
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition: G4Proton.cc:93
#define G4endl
Definition: G4ios.hh:61
const XML_Char const XML_Char * data
G4NeutronElasticXS::~G4NeutronElasticXS ( )
virtual

Definition at line 74 of file G4NeutronElasticXS.cc.

75 {
76  delete fNucleon;
77  /*
78  for(G4int i=0; i<=maxZ; ++i) {
79  delete data[i];
80  }
81  */
82 }

Member Function Documentation

void G4NeutronElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 145 of file G4NeutronElasticXS.cc.

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

146 {
147  if(isInitialized) { return; }
148  if(verboseLevel > 0){
149  G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
150  << p.GetParticleName() << G4endl;
151  }
152  if(p.GetParticleName() != "neutron") {
154  ed << p.GetParticleName() << " is a wrong particle type -"
155  << " only neutron is allowed";
156  G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
157  FatalException, ed, "");
158  return;
159  }
160  isInitialized = true;
161 
162  // check environment variable
163  // Build the complete string identifying the file with the data set
164  char* path = getenv("G4NEUTRONXSDATA");
165 
166  G4DynamicParticle* dynParticle =
168 
169  // Access to elements
170  const G4ElementTable* theElmTable = G4Element::GetElementTable();
171  size_t numOfElm = G4Element::GetNumberOfElements();
172  if(numOfElm > 0) {
173  for(size_t i=0; i<numOfElm; ++i) {
174  G4int Z = G4int(((*theElmTable)[i])->GetZ());
175  if(Z < 1) { Z = 1; }
176  else if(Z > maxZ) { Z = maxZ; }
177  //G4cout << "Z= " << Z << G4endl;
178  // Initialisation
179  if(!data[Z]) { Initialise(Z, dynParticle, path); }
180  }
181  }
182  delete dynParticle;
183 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
const XML_Char const XML_Char * data
void G4NeutronElasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 84 of file G4NeutronElasticXS.cc.

85 {
86  outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
87  << "cross section on nuclei using data from the high precision\n"
88  << "neutron database. These data are simplified and smoothed over\n"
89  << "the resonance region in order to reduce CPU time.\n"
90  << "G4NeutronElasticXS is valid for energies up to 20 MeV, for all\n"
91  << "targets through U.\n";
92 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double G4NeutronElasticXS::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 102 of file G4NeutronElasticXS.cc.

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

104 {
105  G4double xs = 0.0;
106  G4double ekin = aParticle->GetKineticEnergy();
107 
108  if(Z < 1 || Z > maxZ) { return xs; }
109 
110  G4int Amean = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
111  G4PhysicsVector* pv = data[Z];
112  // G4cout << "G4NeutronElasticXS::GetCrossSection e= " << ekin
113  // << " Z= " << Z << G4endl;
114 
115  // element was not initialised
116  if(!pv) {
117  Initialise(Z);
118  pv = data[Z];
119  if(!pv) { return xs; }
120  }
121 
122  G4double e1 = pv->Energy(0);
123  if(ekin <= e1) { return (*pv)[0]; }
124 
125  G4int n = pv->GetVectorLength() - 1;
126  G4double e2 = pv->Energy(n);
127 
128  if(ekin <= e2) {
129  xs = pv->Value(ekin);
130  } else if(1 == Z) {
131  fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
132  xs = coeff[1]*fNucleon->GetElasticHadronNucleonXsc();
133  } else {
134  ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
135  xs = coeff[Z]*ggXsection->GetElasticGlauberGribovXsc();
136  }
137 
138  if(verboseLevel > 0){
139  G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
140  }
141  return xs;
142 }
G4double GetElasticHadronNucleonXsc()
G4double GetKineticEnergy() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
const G4int n
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const XML_Char const XML_Char * data
G4bool G4NeutronElasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 95 of file G4NeutronElasticXS.cc.

97 {
98  return true;
99 }

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