Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronElasticXS Class Reference

#include <G4NeutronElasticXS.hh>

+ Inheritance diagram for G4NeutronElasticXS:

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 IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
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 GetElementCrossSection (const G4DynamicParticle *, G4int Z, 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 BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
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::G4NeutronElasticXS ( )

Definition at line 59 of file G4NeutronElasticXS.cc.

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

◆ ~G4NeutronElasticXS()

G4NeutronElasticXS::~G4NeutronElasticXS ( )
virtual

Definition at line 75 of file G4NeutronElasticXS.cc.

76{
77 delete fNucleon;
78 for(G4int i=0; i<=maxZ; ++i) {
79 delete data[i];
80 }
81}
int G4int
Definition: G4Types.hh:66

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 143 of file G4NeutronElasticXS.cc.

144{
145 if(isInitialized) { return; }
146 if(verboseLevel > 0){
147 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
148 << p.GetParticleName() << G4endl;
149 }
150 if(p.GetParticleName() != "neutron") {
151 throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
152 return;
153 }
154 isInitialized = true;
155
156 // check environment variable
157 // Build the complete string identifying the file with the data set
158 char* path = getenv("G4NEUTRONXSDATA");
159 if (!path){
160 throw G4HadronicException(__FILE__, __LINE__,
161 "G4NEUTRONXSDATA environment variable not defined");
162 return;
163 }
164
165 G4DynamicParticle* dynParticle =
167
168 // Access to elements
169 const G4ElementTable* theElmTable = G4Element::GetElementTable();
170 size_t numOfElm = G4Element::GetNumberOfElements();
171 if(numOfElm > 0) {
172 for(size_t i=0; i<numOfElm; ++i) {
173 G4int Z = G4int(((*theElmTable)[i])->GetZ());
174 if(Z < 1) { Z = 1; }
175 else if(Z > maxZ) { Z = maxZ; }
176 //G4cout << "Z= " << Z << G4endl;
177 // Initialisation
178 if(!data[Z]) { Initialise(Z, dynParticle, path); }
179 }
180 }
181 delete dynParticle;
182}
std::vector< G4Element * > G4ElementTable
CLHEP::Hep3Vector G4ThreeVector
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4String & GetParticleName() const

◆ CrossSectionDescription()

void G4NeutronElasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 83 of file G4NeutronElasticXS.cc.

84{
85 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
86 << "cross section on nuclei using data from the high precision\n"
87 << "neutron database. These data are simplified and smoothed over\n"
88 << "the resonance region in order to reduce CPU time.\n"
89 << "G4NeutronElasticXS is valid for energies up to 20 MeV, for all\n"
90 << "targets through U.\n";
91}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 101 of file G4NeutronElasticXS.cc.

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

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 94 of file G4NeutronElasticXS.cc.

96{
97 return true;
98}

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