Geant4-11
Public Member Functions | Static Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4NeutronElasticXS Class Referencefinal

#include <G4NeutronElasticXS.hh>

Inheritance diagram for G4NeutronElasticXS:
G4VCrossSectionDataSet

Public Member Functions

void BuildPhysicsTable (const G4ParticleDefinition &) final
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
void CrossSectionDescription (std::ostream &) const final
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
bool ForAllAtomsAndEnergies () const
 
 G4NeutronElasticXS ()
 
 G4NeutronElasticXS (const G4NeutronElasticXS &)=delete
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
 
G4double GetMaxKinEnergy () const
 
G4double GetMinKinEnergy () const
 
const G4StringGetName () const
 
virtual G4int GetVerboseLevel () const
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
 
G4NeutronElasticXSoperator= (const G4NeutronElasticXS &right)=delete
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void SetForAllAtomsAndEnergies (G4bool val)
 
void SetMaxKinEnergy (G4double value)
 
void SetMinKinEnergy (G4double value)
 
void SetName (const G4String &nam)
 
virtual void SetVerboseLevel (G4int value)
 
 ~G4NeutronElasticXS () final
 

Static Public Member Functions

static const char * Default_Name ()
 

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Member Functions

const G4StringFindDirectoryPath ()
 
G4PhysicsVectorGetPhysicsVector (G4int Z)
 
void Initialise (G4int Z)
 
void InitialiseOnFly (G4int Z)
 

Private Attributes

G4VComponentCrossSectionggXsection = nullptr
 
G4bool isForAllAtomsAndEnergies
 
G4bool isMaster = false
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
const G4ParticleDefinitionneutron
 
G4CrossSectionDataSetRegistryregistry
 

Static Private Attributes

static G4double coeff [MAXZEL] = {0.0}
 
static G4PhysicsVectordata [MAXZEL] = {nullptr}
 
static G4String gDataDirectory = ""
 
static const G4int MAXZEL = 93
 

Detailed Description

Definition at line 56 of file G4NeutronElasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronElasticXS() [1/2]

G4NeutronElasticXS::G4NeutronElasticXS ( )
explicit

Definition at line 62 of file G4NeutronElasticXS.cc.

65{
66 // verboseLevel = 0;
67 if(verboseLevel > 0){
68 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
69 << MAXZEL << G4endl;
70 }
74}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
const G4ParticleDefinition * neutron
static const G4int MAXZEL
G4VComponentCrossSection * ggXsection
static const char * Default_Name()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

References G4cout, G4endl, G4CrossSectionDataSetRegistry::GetComponentCrossSection(), ggXsection, G4CrossSectionDataSetRegistry::Instance(), MAXZEL, G4VCrossSectionDataSet::SetForAllAtomsAndEnergies(), and G4VCrossSectionDataSet::verboseLevel.

◆ ~G4NeutronElasticXS()

G4NeutronElasticXS::~G4NeutronElasticXS ( )
final

Definition at line 76 of file G4NeutronElasticXS.cc.

77{
78 if(isMaster) {
79 for(G4int i=0; i<MAXZEL; ++i) {
80 delete data[i];
81 data[i] = nullptr;
82 }
83 }
84}
int G4int
Definition: G4Types.hh:85
static G4PhysicsVector * data[MAXZEL]

References data, isMaster, and MAXZEL.

◆ G4NeutronElasticXS() [2/2]

G4NeutronElasticXS::G4NeutronElasticXS ( const G4NeutronElasticXS )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 177 of file G4NeutronElasticXS.cc.

178{
179 if(verboseLevel > 0){
180 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
181 << p.GetParticleName() << G4endl;
182 }
183 if(p.GetParticleName() != "neutron") {
185 ed << p.GetParticleName() << " is a wrong particle type -"
186 << " only neutron is allowed";
187 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
188 FatalException, ed, "");
189 return;
190 }
191 if(0. == coeff[0]) {
192#ifdef G4MULTITHREADED
193 G4MUTEXLOCK(&neutronElasticXSMutex);
194 if(0. == coeff[0]) {
195#endif
196 coeff[0] = 1.0;
197 isMaster = true;
199#ifdef G4MULTITHREADED
200 }
201 G4MUTEXUNLOCK(&neutronElasticXSMutex);
202#endif
203 }
204
205 // it is possible re-initialisation for the second run
206 if(isMaster) {
207
208 // Access to elements
210 for ( auto & elm : *table ) {
211 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZEL-1) );
212 if ( nullptr == data[Z] ) { Initialise(Z); }
213 }
214 }
215}
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
const G4int Z[17]
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
const G4String & FindDirectoryPath()
static G4double coeff[MAXZEL]
const G4String & GetParticleName() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References coeff, data, FatalException, FindDirectoryPath(), G4cout, G4endl, G4Exception(), G4MUTEXLOCK, G4MUTEXUNLOCK, G4Element::GetElementTable(), G4ParticleDefinition::GetParticleName(), Initialise(), isMaster, G4INCL::Math::max(), MAXZEL, G4INCL::Math::min(), G4VCrossSectionDataSet::verboseLevel, and Z.

◆ ComputeCrossSection()

G4double G4VCrossSectionDataSet::ComputeCrossSection ( const G4DynamicParticle part,
const G4Element elm,
const G4Material mat = nullptr 
)
inherited

Definition at line 81 of file G4VCrossSectionDataSet.cc.

84{
85 G4int Z = elm->GetZasInt();
86
87 if (IsElementApplicable(part, Z, mat)) {
88 return GetElementCrossSection(part, Z, mat);
89 }
90
91 // isotope-wise cross section making sum over available
92 // isotope cross sections, which may be incomplete, so
93 // the result is corrected
94 size_t nIso = elm->GetNumberOfIsotopes();
95 G4double fact = 0.0;
96 G4double xsec = 0.0;
97
98 // user-defined isotope abundances
99 const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
100 const G4double* abundVector = elm->GetRelativeAbundanceVector();
101
102 for (size_t j=0; j<nIso; ++j) {
103 const G4Isotope* iso = (*isoVector)[j];
104 G4int A = iso->GetN();
105 if(abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
106 fact += abundVector[j];
107 xsec += abundVector[j]*GetIsoCrossSection(part, Z, A, iso, elm, mat);
108 }
109 }
110 return (fact > 0.0) ? xsec/fact : 0.0;
111}
std::vector< G4Isotope * > G4IsotopeVector
double G4double
Definition: G4Types.hh:83
const G4double A[17]
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZasInt() const
Definition: G4Element.hh:132
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
G4int GetN() const
Definition: G4Isotope.hh:93
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)

References A, G4VCrossSectionDataSet::GetElementCrossSection(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4Element::GetIsotopeVector(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZasInt(), G4VCrossSectionDataSet::IsElementApplicable(), G4VCrossSectionDataSet::IsIsoApplicable(), and Z.

Referenced by G4VCrossSectionDataSet::GetCrossSection().

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 86 of file G4NeutronElasticXS.cc.

87{
88 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
89 << "cross section on nuclei using data from the high precision\n"
90 << "neutron database. These data are simplified and smoothed over\n"
91 << "the resonance region in order to reduce CPU time.\n"
92 << "For high energies Glauber-Gribiv cross section is used.\n";
93}

◆ Default_Name()

static const char * G4NeutronElasticXS::Default_Name ( )
inlinestatic

Definition at line 64 of file G4NeutronElasticXS.hh.

64{return "G4NeutronElasticXS";}

◆ DumpPhysicsTable()

void G4VCrossSectionDataSet::DumpPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ FindDirectoryPath()

const G4String & G4NeutronElasticXS::FindDirectoryPath ( )
private

Definition at line 217 of file G4NeutronElasticXS.cc.

218{
219 // check environment variable
220 // build the complete string identifying the file with the data set
221 if(gDataDirectory.empty()) {
222 char* path = std::getenv("G4PARTICLEXSDATA");
223 if (nullptr != path) {
224 std::ostringstream ost;
225 ost << path << "/neutron/el";
226 gDataDirectory = ost.str();
227 } else {
228 G4Exception("G4NeutronElasticXS::Initialise(..)","had013",
230 "Environment variable G4PARTICLEXSDATA is not defined");
231 }
232 }
233 return gDataDirectory;
234}
static G4String gDataDirectory

References FatalException, G4Exception(), and gDataDirectory.

Referenced by BuildPhysicsTable(), and Initialise().

◆ ForAllAtomsAndEnergies()

bool G4VCrossSectionDataSet::ForAllAtomsAndEnergies ( ) const
inlineinherited

◆ GetCrossSection()

G4double G4VCrossSectionDataSet::GetCrossSection ( const G4DynamicParticle dp,
const G4Element elm,
const G4Material mat = nullptr 
)
inlineinherited

Definition at line 187 of file G4VCrossSectionDataSet.hh.

190{
191 return ComputeCrossSection(dp, elm, mat);
192}
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)

References G4VCrossSectionDataSet::ComputeCrossSection().

◆ GetElementCrossSection()

G4double G4NeutronElasticXS::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material  
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 110 of file G4NeutronElasticXS.cc.

112{
113 G4double xs = 0.0;
114 G4double ekin = aParticle->GetKineticEnergy();
115
116 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
117
118 auto pv = GetPhysicsVector(Z);
119 if(pv == nullptr) { return xs; }
120 // G4cout << "G4NeutronElasticXS::GetCrossSection e= " << ekin
121 // << " Z= " << Z << G4endl;
122
123 if(ekin <= pv->Energy(1)) {
124 xs = (*pv)[1];
125 } else if(ekin <= pv->GetMaxEnergy()) {
126 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
127 } else {
129 ekin, Z, aeff[Z]);
130 }
131
132#ifdef G4VERBOSE
133 if(verboseLevel > 1) {
134 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
135 << ", nElmXSel(b)= " << xs/CLHEP::barn
136 << G4endl;
137 }
138#endif
139 return xs;
140}
static const G4double aeff[]
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4PhysicsVector * GetPhysicsVector(G4int Z)
G4double GetElasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double MeV

References aeff, CLHEP::barn, coeff, G4cout, G4endl, G4VComponentCrossSection::GetElasticElementCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetLogKineticEnergy(), GetPhysicsVector(), ggXsection, MAXZEL, CLHEP::MeV, neutron, G4VCrossSectionDataSet::verboseLevel, and Z.

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

G4double G4NeutronElasticXS::GetIsoCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
G4int  A,
const G4Isotope iso,
const G4Element elm,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 142 of file G4NeutronElasticXS.cc.

147{
148 return GetElementCrossSection(aParticle, Z, mat) * A/aeff[Z];
149}
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) final

References A, aeff, GetElementCrossSection(), and Z.

◆ GetMaxKinEnergy()

G4double G4VCrossSectionDataSet::GetMaxKinEnergy ( ) const
inlineinherited

◆ GetMinKinEnergy()

G4double G4VCrossSectionDataSet::GetMinKinEnergy ( ) const
inlineinherited

◆ GetName()

const G4String & G4VCrossSectionDataSet::GetName ( ) const
inlineinherited

◆ GetPhysicsVector()

G4PhysicsVector * G4NeutronElasticXS::GetPhysicsVector ( G4int  Z)
inlineprivate

Definition at line 116 of file G4NeutronElasticXS.hh.

117{
118 if(nullptr == data[Z]) { InitialiseOnFly(Z); }
119 return data[Z];
120}
void InitialiseOnFly(G4int Z)

References data, InitialiseOnFly(), and Z.

Referenced by GetElementCrossSection().

◆ GetVerboseLevel()

G4int G4VCrossSectionDataSet::GetVerboseLevel ( ) const
inlinevirtualinherited

◆ Initialise()

void G4NeutronElasticXS::Initialise ( G4int  Z)
private

Definition at line 249 of file G4NeutronElasticXS.cc.

250{
251 if(data[Z] != nullptr) { return; }
252
253 // upload data from file
254 data[Z] = new G4PhysicsLogVector();
255
256 std::ostringstream ost;
257 ost << FindDirectoryPath() << Z ;
258 std::ifstream filein(ost.str().c_str());
259 if (!filein.is_open()) {
261 ed << "Data file <" << ost.str().c_str()
262 << "> is not opened!";
263 G4Exception("G4NeutronElasticXS::Initialise(..)","had014",
264 FatalException, ed, "Check G4PARTICLEXSDATA");
265 return;
266 }
267 if(verboseLevel > 1) {
268 G4cout << "file " << ost.str()
269 << " is opened by G4NeutronElasticXS" << G4endl;
270 }
271
272 // retrieve data from DB
273 if(!data[Z]->Retrieve(filein, true)) {
275 ed << "Data file <" << ost.str().c_str()
276 << "> is not retrieved!";
277 G4Exception("G4NeutronElasticXS::Initialise(..)","had015",
278 FatalException, ed, "Check G4PARTICLEXSDATA");
279 return;
280 }
281 // smooth transition
282 G4double sig1 = (*(data[Z]))[data[Z]->GetVectorLength()-1];
283 G4double ehigh = data[Z]->GetMaxEnergy();
285 ehigh, Z, aeff[Z]);
286 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
287}
G4double GetMaxEnergy() const
std::size_t GetVectorLength() const

References aeff, coeff, data, FatalException, FindDirectoryPath(), G4cout, G4endl, G4Exception(), G4VComponentCrossSection::GetElasticElementCrossSection(), G4PhysicsVector::GetMaxEnergy(), G4PhysicsVector::GetVectorLength(), ggXsection, neutron, G4VCrossSectionDataSet::verboseLevel, and Z.

Referenced by BuildPhysicsTable(), and InitialiseOnFly().

◆ InitialiseOnFly()

void G4NeutronElasticXS::InitialiseOnFly ( G4int  Z)
private

Definition at line 236 of file G4NeutronElasticXS.cc.

237{
238#ifdef G4MULTITHREADED
239 G4MUTEXLOCK(&neutronElasticXSMutex);
240 if(data[Z] == nullptr) {
241#endif
242 Initialise(Z);
243#ifdef G4MULTITHREADED
244 }
245 G4MUTEXUNLOCK(&neutronElasticXSMutex);
246#endif
247}

References data, G4MUTEXLOCK, G4MUTEXUNLOCK, Initialise(), and Z.

Referenced by GetPhysicsVector().

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 96 of file G4NeutronElasticXS.cc.

98{
99 return true;
100}

◆ IsIsoApplicable()

G4bool G4NeutronElasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element ,
const G4Material  
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 102 of file G4NeutronElasticXS.cc.

105{
106 return false;
107}

◆ operator=()

G4NeutronElasticXS & G4NeutronElasticXS::operator= ( const G4NeutronElasticXS right)
delete

◆ SelectIsotope()

const G4Isotope * G4NeutronElasticXS::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy,
G4double  logE 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 151 of file G4NeutronElasticXS.cc.

153{
154 size_t nIso = anElement->GetNumberOfIsotopes();
155 const G4Isotope* iso = anElement->GetIsotope(0);
156
157 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
158 if(1 == nIso) { return iso; }
159
160 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
162 G4double sum = 0.0;
163 size_t j;
164
165 // isotope wise cross section not used
166 for (j=0; j<nIso; ++j) {
167 sum += abundVector[j];
168 if(q <= sum) {
169 iso = anElement->GetIsotope(j);
170 break;
171 }
172 }
173 return iso;
174}
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170

References G4UniformRand, G4Element::GetIsotope(), G4Element::GetNumberOfIsotopes(), and G4Element::GetRelativeAbundanceVector().

◆ SetForAllAtomsAndEnergies()

void G4VCrossSectionDataSet::SetForAllAtomsAndEnergies ( G4bool  val)
inlineinherited

◆ SetMaxKinEnergy()

void G4VCrossSectionDataSet::SetMaxKinEnergy ( G4double  value)
inlineinherited

◆ SetMinKinEnergy()

void G4VCrossSectionDataSet::SetMinKinEnergy ( G4double  value)
inlineinherited

◆ SetName()

void G4VCrossSectionDataSet::SetName ( const G4String nam)
inlineinherited

Definition at line 240 of file G4VCrossSectionDataSet.hh.

241{
242 name = nam;
243}

References G4VCrossSectionDataSet::name.

Referenced by G4ParticleHPInelasticData::G4ParticleHPInelasticData().

◆ SetVerboseLevel()

void G4VCrossSectionDataSet::SetVerboseLevel ( G4int  value)
inlinevirtualinherited

Field Documentation

◆ coeff

G4double G4NeutronElasticXS::coeff = {0.0}
staticprivate

Definition at line 107 of file G4NeutronElasticXS.hh.

Referenced by BuildPhysicsTable(), GetElementCrossSection(), and Initialise().

◆ data

G4PhysicsVector * G4NeutronElasticXS::data = {nullptr}
staticprivate

◆ gDataDirectory

G4String G4NeutronElasticXS::gDataDirectory = ""
staticprivate

Definition at line 108 of file G4NeutronElasticXS.hh.

Referenced by FindDirectoryPath().

◆ ggXsection

G4VComponentCrossSection* G4NeutronElasticXS::ggXsection = nullptr
private

Definition at line 100 of file G4NeutronElasticXS.hh.

Referenced by G4NeutronElasticXS(), GetElementCrossSection(), and Initialise().

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ isMaster

G4bool G4NeutronElasticXS::isMaster = false
private

Definition at line 103 of file G4NeutronElasticXS.hh.

Referenced by BuildPhysicsTable(), and ~G4NeutronElasticXS().

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ MAXZEL

const G4int G4NeutronElasticXS::MAXZEL = 93
staticprivate

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ neutron

const G4ParticleDefinition* G4NeutronElasticXS::neutron
private

Definition at line 101 of file G4NeutronElasticXS.hh.

Referenced by GetElementCrossSection(), and Initialise().

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ verboseLevel

G4int G4VCrossSectionDataSet::verboseLevel
protectedinherited

Definition at line 168 of file G4VCrossSectionDataSet.hh.

Referenced by G4BGGNucleonElasticXS::BuildPhysicsTable(), G4BGGPionElasticXS::BuildPhysicsTable(), G4BGGPionInelasticXS::BuildPhysicsTable(), G4GammaNuclearXS::BuildPhysicsTable(), G4NeutronCaptureXS::BuildPhysicsTable(), BuildPhysicsTable(), G4NeutronInelasticXS::BuildPhysicsTable(), G4ParticleInelasticXS::BuildPhysicsTable(), G4BGGNucleonInelasticXS::BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(), G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS(), G4BGGPionElasticXS::G4BGGPionElasticXS(), G4BGGPionInelasticXS::G4BGGPionInelasticXS(), G4GammaNuclearXS::G4GammaNuclearXS(), G4NeutronCaptureXS::G4NeutronCaptureXS(), G4NeutronElasticXS(), G4NeutronInelasticXS::G4NeutronInelasticXS(), G4ParticleInelasticXS::G4ParticleInelasticXS(), G4NeutronCaptureXS::GetElementCrossSection(), GetElementCrossSection(), G4NeutronInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetElementCrossSection(), G4BGGPionElasticXS::GetElementCrossSection(), G4BGGPionInelasticXS::GetElementCrossSection(), G4BGGNucleonInelasticXS::GetElementCrossSection(), G4GammaNuclearXS::GetElementCrossSection(), G4ParticleInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetIsoCrossSection(), G4BGGPionElasticXS::GetIsoCrossSection(), G4BGGPionInelasticXS::GetIsoCrossSection(), G4GammaNuclearXS::GetIsoCrossSection(), G4BGGNucleonInelasticXS::GetIsoCrossSection(), G4VCrossSectionDataSet::GetVerboseLevel(), Initialise(), G4ParticleInelasticXS::IsoCrossSection(), G4NeutronCaptureXS::IsoCrossSection(), G4NeutronInelasticXS::IsoCrossSection(), G4GammaNuclearXS::RetrieveVector(), G4NeutronCaptureXS::RetrieveVector(), G4NeutronInelasticXS::RetrieveVector(), G4ParticleInelasticXS::RetrieveVector(), and G4VCrossSectionDataSet::SetVerboseLevel().


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