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

#include <G4NeutronInelasticXS.hh>

+ Inheritance diagram for G4NeutronInelasticXS:

Public Member Functions

 G4NeutronInelasticXS ()
 
virtual ~G4NeutronInelasticXS ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, 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 G4NeutronInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronInelasticXS()

G4NeutronInelasticXS::G4NeutronInelasticXS ( )

Definition at line 60 of file G4NeutronInelasticXS.cc.

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

◆ ~G4NeutronInelasticXS()

G4NeutronInelasticXS::~G4NeutronInelasticXS ( )
virtual

Definition at line 76 of file G4NeutronInelasticXS.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 229 of file G4NeutronInelasticXS.cc.

230{
231 if(isInitialized) { return; }
232 if(verboseLevel > 0){
233 G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
234 << p.GetParticleName() << G4endl;
235 }
236 if(p.GetParticleName() != "neutron") {
237 throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
238 return;
239 }
240 isInitialized = true;
241
242 // check environment variable
243 // Build the complete string identifying the file with the data set
244 char* path = getenv("G4NEUTRONXSDATA");
245 if (!path){
246 throw G4HadronicException(__FILE__, __LINE__,
247 "G4NEUTRONXSDATA environment variable not defined");
248 return;
249 }
250
251 G4DynamicParticle* dynParticle =
253
254 // Access to elements
255 const G4ElementTable* theElmTable = G4Element::GetElementTable();
256 size_t numOfElm = G4Element::GetNumberOfElements();
257 if(numOfElm > 0) {
258 for(size_t i=0; i<numOfElm; ++i) {
259 G4int Z = G4int(((*theElmTable)[i])->GetZ());
260 if(Z < 1) { Z = 1; }
261 else if(Z > maxZ) { Z = maxZ; }
262 //G4cout << "Z= " << Z << G4endl;
263 // Initialisation
264 if(!data[Z]) { Initialise(Z, dynParticle, path); }
265 }
266 }
267 delete dynParticle;
268}
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 G4NeutronInelasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 84 of file G4NeutronInelasticXS.cc.

85{
86 outFile << "G4NeutronInelasticXS calculates the neutron inelastic 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 << "G4NeutronInelasticXS is valid for energies up to 20 MeV, for\n"
91 << "nuclei through U.\n";
92}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 110 of file G4NeutronInelasticXS.cc.

112{
113 G4double xs = 0.0;
114 G4double ekin = aParticle->GetKineticEnergy();
115
116 if(Z < 1 || Z > maxZ) { return xs; }
117 G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
118 G4PhysicsVector* pv = data[Z];
119 // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin << " Z= " << Z << G4endl;
120
121 // element was not initialised
122 if(!pv) {
123 Initialise(Z);
124 pv = data[Z];
125 if(!pv) { return xs; }
126 }
127
128 G4double e1 = pv->Energy(0);
129 if(ekin <= e1) { return xs; }
130
131 G4int n = pv->GetVectorLength() - 1;
132 G4double e2 = pv->Energy(n);
133 if(ekin <= e2) {
134 xs = pv->Value(ekin);
135 } else if(1 == Z) {
136 fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
137 xs = coeff[1]*fNucleon->GetInelasticHadronNucleonXsc();
138 } else {
139 ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
140 xs = coeff[Z]*ggXsection->GetInelasticGlauberGribovXsc();
141 }
142
143 if(verboseLevel > 0) {
144 G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
145 }
146 return xs;
147}
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 GetInelasticHadronNucleonXsc()
static G4NistManager * Instance()
G4double Value(G4double theEnergy)
size_t GetVectorLength() const
G4double Energy(size_t index) const

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 95 of file G4NeutronInelasticXS.cc.

97{
98 return true;
99}

◆ IsIsoApplicable()

G4bool G4NeutronInelasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element ,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 102 of file G4NeutronInelasticXS.cc.

105{
106 return false;
107}

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