#include <G4NeutronInelasticXS.hh>
|
| 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 G4double | GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) |
|
virtual G4Isotope * | SelectIsotope (const G4Element *, G4double kinEnergy) |
|
virtual void | BuildPhysicsTable (const G4ParticleDefinition &) |
|
virtual void | CrossSectionDescription (std::ostream &) const |
|
| G4VCrossSectionDataSet (const G4String &nam="") |
|
virtual | ~G4VCrossSectionDataSet () |
|
G4double | GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0) |
|
G4double | ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0) |
|
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 G4String & | GetName () const |
|
Definition at line 62 of file G4NeutronInelasticXS.hh.
G4NeutronInelasticXS::G4NeutronInelasticXS |
( |
| ) |
|
Definition at line 93 of file G4NeutronInelasticXS.cc.
99 G4cout <<
"G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
G4VCrossSectionDataSet(const G4String &nam="")
static const char * Default_Name()
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
G4NeutronInelasticXS::~G4NeutronInelasticXS |
( |
| ) |
|
|
virtual |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 279 of file G4NeutronInelasticXS.cc.
282 G4cout <<
"G4NeutronInelasticXS::BuildPhysicsTable for "
288 <<
" only neutron is allowed";
289 G4Exception(
"G4NeutronInelasticXS::BuildPhysicsTable(..)",
"had012",
297 data->SetName(
"NeutronInelastic");
306 char* path = getenv(
"G4NEUTRONXSDATA");
315 for(
size_t i=0; i<numOfElm; ++i) {
321 if(!(
data->GetElementData(Z))) {
322 Initialise(Z, dynParticle, path);
std::ostringstream G4ExceptionDescription
CLHEP::Hep3Vector G4ThreeVector
const G4String & GetParticleName() const
const XML_Char const XML_Char * data
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
static G4Neutron * Neutron()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
void G4NeutronInelasticXS::CrossSectionDescription |
( |
std::ostream & |
outFile | ) |
const |
|
virtual |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 116 of file G4NeutronInelasticXS.cc.
118 outFile <<
"G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
119 <<
"cross section on nuclei using data from the high precision\n"
120 <<
"neutron database. These data are simplified and smoothed over\n"
121 <<
"the resonance region in order to reduce CPU time.\n"
122 <<
"G4NeutronInelasticXS is valid for energies up to 20 MeV, for\n"
123 <<
"nuclei through U.\n";
static const char* G4NeutronInelasticXS::Default_Name |
( |
| ) |
|
|
inlinestatic |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 141 of file G4NeutronInelasticXS.cc.
149 else if(
Z < 1) {
Z = 1; }
160 pv =
data->GetElementData(
Z);
161 if(!pv) {
return xs; }
165 if(ekin <= e1) {
return xs; }
170 xs = pv->
Value(ekin);
180 G4cout <<
"ekin= " << ekin <<
", XSinel= " << xs <<
G4endl;
G4double GetMaxEnergy() const
G4double GetKineticEnergy() const
G4double GetInelasticGlauberGribovXsc()
static G4NistManager * Instance()
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
const XML_Char const XML_Char * data
G4GLOB_DLL std::ostream G4cout
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetInelasticHadronNucleonXsc()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 224 of file G4NeutronInelasticXS.cc.
245 for (j = 0; j<
nIso; ++j) {
246 sum += abundVector[j];
248 iso = (*isoVector)[j];
255 if(!
data->GetElementData(Z)) { Initialise(Z); }
256 size_t nn = temp.size();
257 if(nn < nIso) { temp.resize(nIso, 0.); }
259 for (j=0; j<
nIso; ++j) {
262 sum += abundVector[j]*IsoCrossSection(kinEnergy, Z,
263 (*isoVector)[j]->GetN());
267 for (j = 0; j<
nIso; ++j) {
269 iso = (*isoVector)[j];
size_t GetNumberOfIsotopes() const
std::vector< G4Isotope * > G4IsotopeVector
const XML_Char const XML_Char * data
G4double * GetRelativeAbundanceVector() const
G4IsotopeVector * GetIsotopeVector() const
The documentation for this class was generated from the following files: