Geant4-11
Public Member Functions | Protected Attributes | Private Attributes
G4ParticleHPFissionData Class Reference

#include <G4ParticleHPFissionData.hh>

Inheritance diagram for G4ParticleHPFissionData:
G4VCrossSectionDataSet

Public Member Functions

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

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Attributes

const G4Elementelement_cache
 
G4bool instanceOfWorker
 
G4bool isForAllAtomsAndEnergies
 
G4double ke_cache
 
const G4Materialmaterial_cache
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4CrossSectionDataSetRegistryregistry
 
G4PhysicsTabletheCrossSections
 
G4double xs_cache
 

Detailed Description

Definition at line 48 of file G4ParticleHPFissionData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFissionData()

G4ParticleHPFissionData::G4ParticleHPFissionData ( )

Definition at line 48 of file G4ParticleHPFissionData.cc.

49:G4VCrossSectionDataSet("NeutronHPFissionXS")
50{
51 SetMinKinEnergy( 0*MeV );
52 SetMaxKinEnergy( 20*MeV );
53
55 instanceOfWorker = false;
57 instanceOfWorker = true;
58 }
59 element_cache = NULL;
60 material_cache = NULL;
61 ke_cache = 0.0;
62 xs_cache = 0.0;
63 //BuildPhysicsTable(*G4Neutron::Neutron());
64}
static constexpr double MeV
Definition: G4SIunits.hh:200
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

References element_cache, instanceOfWorker, G4Threading::IsWorkerThread(), ke_cache, material_cache, MeV, G4VCrossSectionDataSet::SetMaxKinEnergy(), G4VCrossSectionDataSet::SetMinKinEnergy(), theCrossSections, and xs_cache.

◆ ~G4ParticleHPFissionData()

G4ParticleHPFissionData::~G4ParticleHPFissionData ( )

Definition at line 66 of file G4ParticleHPFissionData.cc.

67{
68 if ( theCrossSections != NULL && instanceOfWorker != true ) {
70 delete theCrossSections;
71 theCrossSections = NULL;
72 }
73}
void clearAndDestroy()

References G4PhysicsTable::clearAndDestroy(), instanceOfWorker, and theCrossSections.

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPFissionData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 114 of file G4ParticleHPFissionData.cc.

115{
116 if(&aP!=G4Neutron::Neutron())
117 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
118
121 return;
122 }
123
124 size_t numberOfElements = G4Element::GetNumberOfElements();
125 //theCrossSections = new G4PhysicsTable( numberOfElements );
126 // TKDB
127 //if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
128 if ( theCrossSections == NULL )
129 theCrossSections = new G4PhysicsTable( numberOfElements );
130 else
132
133 // make a PhysicsVector for each element
134
135 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
136 for( size_t i=0; i<numberOfElements; ++i )
137 {
139 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
140 theCrossSections->push_back(physVec);
141 }
142
144}
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
G4PhysicsTable * GetFissionCrossSections()
static G4ParticleHPManager * GetInstance()
void RegisterFissionCrossSections(G4PhysicsTable *val)
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition: tls.hh:77

References G4PhysicsTable::clearAndDestroy(), G4ThreadLocal, G4Element::GetElementTable(), G4ParticleHPManager::GetFissionCrossSections(), G4ParticleHPManager::GetInstance(), G4Element::GetNumberOfElements(), G4ParticleHPData::Instance(), G4Threading::IsWorkerThread(), G4ParticleHPData::MakePhysicsVector(), G4Neutron::Neutron(), G4PhysicsTable::push_back(), G4ParticleHPManager::RegisterFissionCrossSections(), and theCrossSections.

◆ 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
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
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 G4ParticleHPFissionData::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 286 of file G4ParticleHPFissionData.cc.

287{
288 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n" ;
289}

◆ DumpPhysicsTable()

void G4ParticleHPFissionData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 146 of file G4ParticleHPFissionData.cc.

147{
148 if(&aP!=G4Neutron::Neutron())
149 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
150
151 #ifdef G4VERBOSE
152 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
153
154//
155// Dump element based cross section
156// range 10e-5 eV to 20 MeV
157// 10 point per decade
158// in barn
159//
160
161 G4cout << G4endl;
162 G4cout << G4endl;
163 G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
164 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
165 G4cout << G4endl;
166 G4cout << "Name of Element" << G4endl;
167 G4cout << "Energy[eV] XS[barn]" << G4endl;
168 G4cout << G4endl;
169
170 size_t numberOfElements = G4Element::GetNumberOfElements();
171 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
172
173 for ( size_t i = 0 ; i < numberOfElements ; ++i )
174 {
175
176 G4cout << (*theElementTable)[i]->GetName() << G4endl;
177
178 if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
179 {
180 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
181 G4cout << G4endl;
182 continue;
183 }
184
185 G4int ie = 0;
186
187 for ( ie = 0 ; ie < 130 ; ie++ )
188 {
189 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
190 G4bool outOfRange = false;
191
192 if ( eKinetic < 20*MeV )
193 {
194 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
195 }
196
197 }
198
199 G4cout << G4endl;
200 }
201
202 //G4cout << "G4ParticleHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
203 #endif
204}
static constexpr double barn
Definition: G4SIunits.hh:85
static constexpr double eV
Definition: G4SIunits.hh:201
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4HadronicParameters * Instance()
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230

References barn, eV, G4cout, G4endl, G4ThreadLocal, G4Element::GetElementTable(), G4Pow::GetInstance(), G4Element::GetNumberOfElements(), GetVerboseLevel(), G4HadronicParameters::Instance(), MeV, G4Neutron::Neutron(), G4Pow::powA(), and theCrossSections.

◆ ForAllAtomsAndEnergies()

bool G4VCrossSectionDataSet::ForAllAtomsAndEnergies ( ) const
inlineinherited

◆ GetCrossSection() [1/2]

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().

◆ GetCrossSection() [2/2]

G4double G4ParticleHPFissionData::GetCrossSection ( const G4DynamicParticle aP,
const G4Element anE,
G4double  aT 
)

Definition at line 208 of file G4ParticleHPFissionData.cc.

210{
211 G4double result = 0;
212 if(anE->GetZ()<88) return result;
213 G4bool outOfRange;
214 G4int index = anE->GetIndex();
215
216// 100729 TK add safety
217if ( ( ( *theCrossSections )( index ) )->GetVectorLength() == 0 ) return result;
218
219 // prepare neutron
220 G4double eKinetic = aP->GetKineticEnergy();
221 G4ReactionProduct theNeutronRP( aP->GetDefinition() );
222 theNeutronRP.SetMomentum( aP->GetMomentum() );
223 theNeutronRP.SetKineticEnergy( eKinetic );
224
225 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
226 {
227 //NEGLECT_DOPPLER
228 G4double factor = 1.0;
229 if ( eKinetic < aT * k_Boltzmann ) {
230 // below 0.1 eV neutrons
231 // Have to do some, but now just igonre.
232 // Will take care after performance check.
233 // factor = factor * targetV;
234 }
235 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
236 }
237
238 // prepare thermal nucleus
239 G4Nucleus aNuc;
240 G4double eps = 0.0001;
241 G4double theA = anE->GetN();
242 G4double theZ = anE->GetZ();
243 G4double eleMass;
244 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
246
247 G4ReactionProduct boosted;
248 G4double aXsection;
249
250 // MC integration loop
251 G4int counter = 0;
252 G4double buffer = 0;
253 G4int size = G4int(std::max(10., aT/60*kelvin));
254 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutronRP.GetMomentum();
255 G4double neutronVMag = neutronVelocity.mag();
256
257 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
258 {
259 if(counter) buffer = result/counter;
260 while (counter<size) // Loop checking, 11.05.2015, T. Koi
261 {
262 counter ++;
263 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
264 boosted.Lorentz(theNeutronRP, aThermalNuc);
265 G4double theEkin = boosted.GetKineticEnergy();
266 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
267 // velocity correction.
268 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
269 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
270 result += aXsection;
271 }
272 size += size;
273 }
274 result /= counter;
275 return result;
276}
static const G4double eps
static constexpr double kelvin
Definition: G4SIunits.hh:274
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
size_t GetIndex() const
Definition: G4Element.hh:182
G4double GetN() const
Definition: G4Element.hh:135
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:236
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
float k_Boltzmann
Definition: hepunit.py:298
#define buffer
Definition: xmlparse.cc:628

References buffer, eps, G4DynamicParticle::GetDefinition(), G4Element::GetIndex(), G4ParticleHPManager::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4ReactionProduct::GetKineticEnergy(), G4ReactionProduct::GetMass(), G4DynamicParticle::GetMomentum(), G4ReactionProduct::GetMomentum(), G4Element::GetN(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGMass(), G4Nucleus::GetThermalNucleus(), G4Element::GetZ(), source.hepunit::k_Boltzmann, kelvin, G4ReactionProduct::Lorentz(), CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), G4Neutron::Neutron(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMomentum(), and theCrossSections.

Referenced by GetIsoCrossSection().

◆ GetElementCrossSection()

G4double G4VCrossSectionDataSet::GetElementCrossSection ( const G4DynamicParticle dynPart,
G4int  Z,
const G4Material mat = nullptr 
)
virtualinherited

Reimplemented in G4EMDissociationCrossSection, G4IonsShenCrossSection, G4NeutrinoElectronCcXsc, G4NeutrinoElectronNcXsc, G4NeutrinoElectronTotXsc, G4NeutronElectronElXsc, G4PhotoNuclearCrossSection, G4NeutronCaptureXS, G4NeutronElasticXS, G4NeutronInelasticXS, G4ElectroNuclearCrossSection, G4BGGNucleonElasticXS, G4BGGPionElasticXS, G4BGGPionInelasticXS, G4BGGNucleonInelasticXS, G4CrossSectionElastic, G4CrossSectionInelastic, G4GammaNuclearXS, G4ParticleInelasticXS, G4ZeroXS, G4NucleonNuclearCrossSection, G4MuNeutrinoNucleusTotXsc, and G4KokoulinMuonNuclearXS.

Definition at line 114 of file G4VCrossSectionDataSet.cc.

117{
119 ed << "GetElementCrossSection is not implemented in <" << name << ">\n"
120 << "Particle: " << dynPart->GetDefinition()->GetParticleName()
121 << " Ekin(MeV)= " << dynPart->GetKineticEnergy()/MeV;
122 if(mat) { ed << " material: " << mat->GetName(); }
123 ed << " target Z= " << Z << G4endl;
124 G4Exception("G4VCrossSectionDataSet::GetElementCrossSection", "had001",
125 FatalException, ed);
126 return 0.0;
127}
@ 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
const G4String & GetName() const
Definition: G4Material.hh:173
const G4String & GetParticleName() const

References FatalException, G4endl, G4Exception(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4Material::GetName(), G4ParticleDefinition::GetParticleName(), MeV, G4VCrossSectionDataSet::name, and Z.

Referenced by G4QMDReaction::ApplyYourself(), G4VCrossSectionDataSet::ComputeCrossSection(), G4GammaNuclearXS::GetElementCrossSection(), G4GammaNuclearXS::GetIsoCrossSection(), and G4GammaNuclearXS::Initialise().

◆ GetIsoCrossSection()

G4double G4ParticleHPFissionData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 88 of file G4ParticleHPFissionData.cc.

93{
94 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
95
97 element_cache = element;
99 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
100 xs_cache = xs;
101 return xs;
102}
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
string material
Definition: eplot.py:19

References element_cache, GetCrossSection(), G4DynamicParticle::GetKineticEnergy(), ke_cache, eplot::material, material_cache, and xs_cache.

◆ GetMaxKinEnergy()

G4double G4VCrossSectionDataSet::GetMaxKinEnergy ( ) const
inlineinherited

◆ GetMinKinEnergy()

G4double G4VCrossSectionDataSet::GetMinKinEnergy ( ) const
inlineinherited

◆ GetName()

const G4String & G4VCrossSectionDataSet::GetName ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4ParticleHPFissionData::GetVerboseLevel ( ) const
virtual

◆ IsElementApplicable()

G4bool G4VCrossSectionDataSet::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = nullptr 
)
virtualinherited

◆ IsIsoApplicable()

G4bool G4ParticleHPFissionData::IsIsoApplicable ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Element ,
const G4Material  
)
virtual

◆ SelectIsotope()

const G4Isotope * G4VCrossSectionDataSet::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy,
G4double  logE 
)
virtualinherited

Reimplemented in G4GammaNuclearXS, G4NeutronCaptureXS, G4NeutronElasticXS, G4NeutronInelasticXS, and G4ParticleInelasticXS.

Definition at line 149 of file G4VCrossSectionDataSet.cc.

151{
152 size_t nIso = anElement->GetNumberOfIsotopes();
153 const G4Isotope* iso = anElement->GetIsotope(0);
154
155 // more than 1 isotope
156 if(1 < nIso) {
157 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
158 G4double sum = 0.0;
160 for (size_t j=0; j<nIso; ++j) {
161 sum += abundVector[j];
162 if(q <= sum) {
163 iso = anElement->GetIsotope(j);
164 break;
165 }
166 }
167 }
168 return iso;
169}
#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 G4ParticleHPFissionData::SetVerboseLevel ( G4int  newValue)
virtual

Field Documentation

◆ element_cache

const G4Element* G4ParticleHPFissionData::element_cache
private

Definition at line 92 of file G4ParticleHPFissionData.hh.

Referenced by G4ParticleHPFissionData(), and GetIsoCrossSection().

◆ instanceOfWorker

G4bool G4ParticleHPFissionData::instanceOfWorker
private

Definition at line 88 of file G4ParticleHPFissionData.hh.

Referenced by G4ParticleHPFissionData(), and ~G4ParticleHPFissionData().

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ ke_cache

G4double G4ParticleHPFissionData::ke_cache
private

Definition at line 90 of file G4ParticleHPFissionData.hh.

Referenced by G4ParticleHPFissionData(), and GetIsoCrossSection().

◆ material_cache

const G4Material* G4ParticleHPFissionData::material_cache
private

Definition at line 93 of file G4ParticleHPFissionData.hh.

Referenced by G4ParticleHPFissionData(), and GetIsoCrossSection().

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ theCrossSections

G4PhysicsTable* G4ParticleHPFissionData::theCrossSections
private

◆ 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(), G4NeutronElasticXS::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::G4NeutronElasticXS(), G4NeutronInelasticXS::G4NeutronInelasticXS(), G4ParticleInelasticXS::G4ParticleInelasticXS(), G4NeutronCaptureXS::GetElementCrossSection(), G4NeutronElasticXS::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(), G4NeutronElasticXS::Initialise(), G4ParticleInelasticXS::IsoCrossSection(), G4NeutronCaptureXS::IsoCrossSection(), G4NeutronInelasticXS::IsoCrossSection(), G4GammaNuclearXS::RetrieveVector(), G4NeutronCaptureXS::RetrieveVector(), G4NeutronInelasticXS::RetrieveVector(), G4ParticleInelasticXS::RetrieveVector(), and G4VCrossSectionDataSet::SetVerboseLevel().

◆ xs_cache

G4double G4ParticleHPFissionData::xs_cache
private

Definition at line 91 of file G4ParticleHPFissionData.hh.

Referenced by G4ParticleHPFissionData(), and GetIsoCrossSection().


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