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

#include <G4Element.hh>

Public Member Functions

 G4Element (const G4String &name, const G4String &symbol, G4double Zeff, G4double Aeff)
 
 G4Element (const G4String &name, const G4String &symbol, G4int nbIsotopes)
 
void AddIsotope (G4Isotope *isotope, G4double RelativeAbundance)
 
virtual ~G4Element ()
 
const G4StringGetName () const
 
const G4StringGetSymbol () const
 
G4double GetZ () const
 
G4double GetN () const
 
G4double GetAtomicMassAmu () const
 
G4double GetA () const
 
G4bool GetNaturalAbandancesFlag ()
 
void SetNaturalAbandancesFlag (G4bool)
 
G4int GetNbOfAtomicShells () const
 
G4double GetAtomicShell (G4int index) const
 
G4int GetNbOfShellElectrons (G4int index) const
 
size_t GetNumberOfIsotopes () const
 
G4IsotopeVectorGetIsotopeVector () const
 
G4doubleGetRelativeAbundanceVector () const
 
const G4IsotopeGetIsotope (G4int iso) const
 
size_t GetIndex () const
 
G4int GetCountUse () const
 
void increaseCountUse ()
 
void decreaseCountUse ()
 
G4int GetIndexZ () const
 
G4double GetfCoulomb () const
 
G4double GetfRadTsai () const
 
G4IonisParamElmGetIonisation () const
 
G4int operator== (const G4Element &) const
 
G4int operator!= (const G4Element &) const
 
 G4Element (__void__ &)
 
void SetName (const G4String &name)
 

Static Public Member Functions

static const G4ElementTableGetElementTable ()
 
static size_t GetNumberOfElements ()
 
static G4ElementGetElement (G4String name, G4bool warning=true)
 

Friends

std::ostream & operator<< (std::ostream &, G4Element *)
 
std::ostream & operator<< (std::ostream &, G4Element &)
 
std::ostream & operator<< (std::ostream &, G4ElementTable)
 

Detailed Description

Definition at line 97 of file G4Element.hh.

Constructor & Destructor Documentation

◆ G4Element() [1/3]

G4Element::G4Element ( const G4String name,
const G4String symbol,
G4double  Zeff,
G4double  Aeff 
)

Definition at line 74 of file G4Element.cc.

76 : fName(name), fSymbol(symbol)
77{
78 G4int iz = (G4int)zeff;
79 if (zeff<1.) {
81 ed << "Fail to create G4Element " << name
82 << " Z= " << zeff << " < 1 !" << G4endl;
83 G4Exception ("G4Element::G4Element()", "mat011", FatalException, ed);
84 }
85 if (std::abs(zeff - iz) > perMillion) {
87 ed << "G4Element Warning: " << name << " Z= " << zeff
88 << " A= " << aeff/(g/mole) << G4endl;
89 G4Exception("G4Element::G4Element()", "mat017", JustWarning, ed);
90 }
91
92 InitializePointers();
93
94 fZeff = zeff;
95 fNeff = aeff/(g/mole);
96 fAeff = aeff;
97
98 if(fNeff < 1.0) fNeff = 1.0;
99
100 if (fNeff < zeff) {
102 ed << "Fail to create G4Element " << name
103 << " with Z= " << zeff << " N= " << fNeff
104 << " N < Z is not allowed" << G4endl;
105 G4Exception ("G4Element::G4Element()", "mat012", FatalException, ed);
106 }
107
108 fNbOfAtomicShells = G4AtomicShells::GetNumberOfShells(iz);
109 fAtomicShells = new G4double[fNbOfAtomicShells];
110 fNbOfShellElectrons = new G4int[fNbOfAtomicShells];
111
112 AddNaturalIsotopes();
113
114 for (G4int i=0;i<fNbOfAtomicShells;i++)
115 {
116 fAtomicShells[i] = G4AtomicShells::GetBindingEnergy(iz, i);
117 fNbOfShellElectrons[i] = G4AtomicShells::GetNumberOfElectrons(iz, i);
118 }
119 ComputeDerivedQuantities();
120}
@ JustWarning
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
static G4int GetNumberOfShells(G4int Z)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ G4Element() [2/3]

G4Element::G4Element ( const G4String name,
const G4String symbol,
G4int  nbIsotopes 
)

Definition at line 127 of file G4Element.cc.

129 : fName(name),fSymbol(symbol)
130{
131 InitializePointers();
132
133 size_t n = size_t(nIsotopes);
134
135 if(0 == nIsotopes) {
136 AddNaturalIsotopes();
137 } else {
138 theIsotopeVector = new G4IsotopeVector(n,0);
139 fRelativeAbundanceVector = new G4double[nIsotopes];
140 }
141}
std::vector< G4Isotope * > G4IsotopeVector

◆ ~G4Element()

G4Element::~G4Element ( )
virtual

Definition at line 255 of file G4Element.cc.

256{
257 // G4cout << "### Destruction of element " << fName << " started" <<G4endl;
258
259 if (theIsotopeVector) { delete theIsotopeVector; }
260 if (fRelativeAbundanceVector) { delete [] fRelativeAbundanceVector; }
261 if (fAtomicShells) { delete [] fAtomicShells; }
262 if (fNbOfShellElectrons) { delete [] fNbOfShellElectrons; }
263 if (fIonisation) { delete fIonisation; }
264
265 //remove this element from theElementTable
266 theElementTable[fIndexInTable] = 0;
267}

◆ G4Element() [3/3]

G4Element::G4Element ( __void__ &  )

Definition at line 247 of file G4Element.cc.

248 : fZeff(0), fNeff(0), fAeff(0)
249{
250 InitializePointers();
251}

Member Function Documentation

◆ AddIsotope()

void G4Element::AddIsotope ( G4Isotope isotope,
G4double  RelativeAbundance 
)

Definition at line 147 of file G4Element.cc.

148{
149 if (theIsotopeVector == 0) {
151 ed << "Fail to add Isotope to G4Element " << fName
152 << " with Z= " << fZeff << " N= " << fNeff << G4endl;
153 G4Exception ("G4Element::AddIsotope()", "mat013", FatalException, ed);
154 return;
155 }
156 G4int iz = isotope->GetZ();
157
158 // filling ...
159 if ( fNumberOfIsotopes < theIsotopeVector->size() ) {
160 // check same Z
161 if (fNumberOfIsotopes==0) { fZeff = G4double(iz); }
162 else if (G4double(iz) != fZeff) {
164 ed << "Fail to add Isotope Z= " << iz << " to G4Element " << fName
165 << " with different Z= " << fZeff << fNeff
166 << G4endl;
167 G4Exception ("G4Element::AddIsotope()", "mat014", FatalException, ed);
168 return;
169 }
170 //Z ok
171 fRelativeAbundanceVector[fNumberOfIsotopes] = abundance;
172 (*theIsotopeVector)[fNumberOfIsotopes] = isotope;
173 ++fNumberOfIsotopes;
174 isotope->increaseCountUse();
175 } else {
177 ed << "Fail to add Isotope Z= " << iz << " to G4Element " << fName
178 << " - more isotopes than declaired " << G4endl;
179 G4Exception ("G4Element::AddIsotope()", "mat015", FatalException, ed);
180 return;
181 }
182
183 // filled.
184 if ( fNumberOfIsotopes == theIsotopeVector->size() ) {
185 // Compute Neff, Aeff
186 G4double wtSum=0.0;
187
188 fNeff = fAeff = 0.0;
189 for (size_t i=0;i<fNumberOfIsotopes;i++) {
190 fAeff += fRelativeAbundanceVector[i]*(*theIsotopeVector)[i]->GetA();
191 fNeff += fRelativeAbundanceVector[i]*(*theIsotopeVector)[i]->GetN();
192 wtSum += fRelativeAbundanceVector[i];
193 }
194 fAeff /= wtSum;
195 fNeff /= wtSum;
196
197 if(wtSum != 1.0) {
198 for(size_t i=0; i<fNumberOfIsotopes; ++i) {
199 fRelativeAbundanceVector[i] /= wtSum;
200 }
201 }
202
203 fNbOfAtomicShells = G4AtomicShells::GetNumberOfShells(iz);
204 fAtomicShells = new G4double[fNbOfAtomicShells];
205 fNbOfShellElectrons = new G4int[fNbOfAtomicShells];
206
207 for ( G4int j = 0; j < fNbOfAtomicShells; j++ )
208 {
209 fAtomicShells[j] = G4AtomicShells::GetBindingEnergy(iz, j);
210 fNbOfShellElectrons[j] = G4AtomicShells::GetNumberOfElectrons(iz, j);
211 }
212 ComputeDerivedQuantities();
213
214 }
215}
void increaseCountUse()
Definition: G4Isotope.hh:135
G4int GetZ() const
Definition: G4Isotope.hh:91

Referenced by G4tgbElement::BuildG4ElementFromIsotopes(), and G4GDMLReadMaterials::MixtureRead().

◆ decreaseCountUse()

void G4Element::decreaseCountUse ( )
inline

Definition at line 193 of file G4Element.hh.

193{fCountUse--;}

◆ GetA()

◆ GetAtomicMassAmu()

G4double G4Element::GetAtomicMassAmu ( ) const
inline

Definition at line 135 of file G4Element.hh.

135{return fNeff;}

◆ GetAtomicShell()

G4double G4Element::GetAtomicShell ( G4int  index) const

Definition at line 367 of file G4Element.cc.

368{
369 if (i<0 || i>=fNbOfAtomicShells) {
371 ed << "Invalid argument " << i << " in for G4Element " << fName
372 << " with Z= " << fZeff
373 << " and Nshells= " << fNbOfAtomicShells
374 << G4endl;
375 G4Exception("G4Element::GetAtomicShell()", "mat016", FatalException, ed);
376 return 0.0;
377 }
378 return fAtomicShells[i];
379}

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), and G4PEEffectModel::SampleSecondaries().

◆ GetCountUse()

G4int G4Element::GetCountUse ( ) const
inline

Definition at line 191 of file G4Element.hh.

191{return fCountUse;}

◆ GetElement()

G4Element * G4Element::GetElement ( G4String  name,
G4bool  warning = true 
)
static

Definition at line 413 of file G4Element.cc.

414{
415 // search the element by its name
416 for (size_t J=0 ; J<theElementTable.size() ; J++)
417 {
418 if (theElementTable[J]->GetName() == elementName)
419 return theElementTable[J];
420 }
421
422 // the element does not exist in the table
423 if (warning) {
424 G4cout << "\n---> warning from G4Element::GetElement(). The element: "
425 << elementName << " does not exist in the table. Return NULL pointer."
426 << G4endl;
427 }
428 return 0;
429}
G4DLLIMPORT std::ostream G4cout
const G4String & GetName() const
Definition: G4Element.hh:127

Referenced by G4GDMLReadMaterials::GetElement().

◆ GetElementTable()

const G4ElementTable * G4Element::GetElementTable ( )
static

Definition at line 399 of file G4Element.cc.

400{
401 return &theElementTable;
402}

Referenced by G4ElectroNuclearReaction::ApplyYourself(), G4AdjointCSManager::BuildCrossSectionMatrices(), G4KokoulinMuonNuclearXS::BuildCrossSectionTable(), G4VRangeToEnergyConverter::BuildLossTable(), G4NeutronCaptureXS::BuildPhysicsTable(), G4NeutronElasticXS::BuildPhysicsTable(), G4NeutronInelasticXS::BuildPhysicsTable(), G4NeutronHPCaptureData::BuildPhysicsTable(), G4NeutronHPElasticData::BuildPhysicsTable(), G4NeutronHPFissionData::BuildPhysicsTable(), G4NeutronHPInelasticData::BuildPhysicsTable(), G4NeutronHPJENDLHEData::BuildPhysicsTable(), G4NeutronHPorLEInelasticData::BuildPhysicsTable(), G4NeutronHPThermalScatteringData::BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4LENDModel::create_used_target_map(), G4NeutronHPCaptureData::DumpPhysicsTable(), G4NeutronHPElasticData::DumpPhysicsTable(), G4NeutronHPFissionData::DumpPhysicsTable(), G4NeutronHPInelasticData::DumpPhysicsTable(), G4NistElementBuilder::FindOrBuildElement(), G4FissLib::G4FissLib(), G4NeutronHPCapture::G4NeutronHPCapture(), G4NeutronHPData::G4NeutronHPData(), G4NeutronHPElastic::G4NeutronHPElastic(), G4NeutronHPFission::G4NeutronHPFission(), G4NeutronHPInelastic::G4NeutronHPInelastic(), G4NeutronHPorLCapture::G4NeutronHPorLCapture(), G4NeutronHPorLEInelastic::G4NeutronHPorLEInelastic(), G4NeutronHPorLElastic::G4NeutronHPorLElastic(), G4NeutronHPorLFission::G4NeutronHPorLFission(), G4NeutronIsotopeProduction::G4NeutronIsotopeProduction(), G4NistManager::GetElement(), G4DiffuseElastic::Initialise(), G4NuclNuclDiffuseElastic::Initialise(), G4SeltzerBergerModel::Initialise(), G4NistManager::PrintG4Element(), G4GDMLRead::StripNames(), and G4NistManager::~G4NistManager().

◆ GetfCoulomb()

◆ GetfRadTsai()

G4double G4Element::GetfRadTsai ( ) const
inline

Definition at line 205 of file G4Element.hh.

205{return fRadTsai;}

◆ GetIndex()

◆ GetIndexZ()

G4int G4Element::GetIndexZ ( ) const
inline

Definition at line 197 of file G4Element.hh.

197{return fIndexZ;}

◆ GetIonisation()

◆ GetIsotope()

◆ GetIsotopeVector()

◆ GetN()

G4double G4Element::GetN ( ) const
inline

Definition at line 134 of file G4Element.hh.

134{return fNeff;}

Referenced by G4MuonMinusCaptureAtRest::AtRestDoIt(), G4Nucleus::ChooseParameters(), G4VEmModel::ComputeCrossSectionPerAtom(), G4EmCalculator::ComputeCrossSectionPerAtom(), G4NeutronHPCaptureData::GetCrossSection(), G4NeutronHPElasticData::GetCrossSection(), G4NeutronHPFissionData::GetCrossSection(), G4NeutronHPInelasticData::GetCrossSection(), G4NeutronHPorLCaptureData::GetCrossSection(), G4NeutronHPorLEInelasticData::GetCrossSection(), G4NeutronHPorLElasticData::GetCrossSection(), G4NeutronHPorLFissionData::GetCrossSection(), G4VComponentCrossSection::GetElasticElementCrossSection(), G4ChargeExchangeProcess::GetElementCrossSection(), G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc(), G4ComponentGGNuclNuclXsc::GetHadronNucleonXsc(), G4GGNuclNuclCrossSection::GetHadronNucleonXsc(), G4GlauberGribovCrossSection::GetHadronNucleonXsc(), G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS(), G4GlauberGribovCrossSection::GetHadronNucleonXscNS(), G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG(), G4GlauberGribovCrossSection::GetHadronNucleonXscPDG(), G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc(), G4GlauberGribovCrossSection::GetHNinelasticXsc(), G4VComponentCrossSection::GetInelasticElementCrossSection(), G4ComponentGGHadronNucleusXsc::GetNucleusRadius(), G4ComponentGGNuclNuclXsc::GetNucleusRadius(), G4GGNuclNuclCrossSection::GetNucleusRadius(), G4GlauberGribovCrossSection::GetNucleusRadius(), G4NeutronHPThermalBoost::GetThermalEnergy(), G4VComponentCrossSection::GetTotalElementCrossSection(), G4HadronNucleonXsc::IsApplicable(), and G4VEmModel::SelectIsotopeNumber().

◆ GetName()

◆ GetNaturalAbandancesFlag()

G4bool G4Element::GetNaturalAbandancesFlag ( )
inline

Definition at line 276 of file G4Element.hh.

277{
278 return fNaturalAbandances;
279}

◆ GetNbOfAtomicShells()

◆ GetNbOfShellElectrons()

G4int G4Element::GetNbOfShellElectrons ( G4int  index) const

Definition at line 383 of file G4Element.cc.

384{
385 if (i<0 || i>=fNbOfAtomicShells) {
387 ed << "Invalid argument " << i << " for G4Element " << fName
388 << " with Z= " << fZeff
389 << " and Nshells= " << fNbOfAtomicShells
390 << G4endl;
391 G4Exception("G4Element::GetNbOfShellElectrons()", "mat016", FatalException, ed);
392 return 0;
393 }
394 return fNbOfShellElectrons[i];
395}

Referenced by G4KleinNishinaModel::SampleSecondaries().

◆ GetNumberOfElements()

size_t G4Element::GetNumberOfElements ( )
static

Definition at line 406 of file G4Element.cc.

407{
408 return theElementTable.size();
409}

Referenced by G4NeutronHPCapture::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4KokoulinMuonNuclearXS::BuildCrossSectionTable(), G4VRangeToEnergyConverter::BuildLossTable(), G4NeutronCaptureXS::BuildPhysicsTable(), G4NeutronElasticXS::BuildPhysicsTable(), G4NeutronInelasticXS::BuildPhysicsTable(), G4NeutronHPCaptureData::BuildPhysicsTable(), G4NeutronHPElasticData::BuildPhysicsTable(), G4NeutronHPFissionData::BuildPhysicsTable(), G4NeutronHPInelasticData::BuildPhysicsTable(), G4NeutronHPJENDLHEData::BuildPhysicsTable(), G4NeutronHPorLEInelasticData::BuildPhysicsTable(), G4NeutronHPThermalScatteringData::BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4LENDModel::create_used_target_map(), G4NeutronHPCaptureData::DumpPhysicsTable(), G4NeutronHPElasticData::DumpPhysicsTable(), G4NeutronHPFissionData::DumpPhysicsTable(), G4NeutronHPInelasticData::DumpPhysicsTable(), G4FissLib::G4FissLib(), G4NeutronHPCapture::G4NeutronHPCapture(), G4NeutronHPData::G4NeutronHPData(), G4NeutronHPElastic::G4NeutronHPElastic(), G4NeutronHPFission::G4NeutronHPFission(), G4NeutronHPInelastic::G4NeutronHPInelastic(), G4NeutronHPorLCapture::G4NeutronHPorLCapture(), G4NeutronHPorLEInelastic::G4NeutronHPorLEInelastic(), G4NeutronHPorLElastic::G4NeutronHPorLElastic(), G4NeutronHPorLFission::G4NeutronHPorLFission(), G4NeutronIsotopeProduction::G4NeutronIsotopeProduction(), G4DiffuseElastic::Initialise(), G4NuclNuclDiffuseElastic::Initialise(), G4SeltzerBergerModel::Initialise(), and G4VRangeToEnergyConverter::operator=().

◆ GetNumberOfIsotopes()

◆ GetRelativeAbundanceVector()

◆ GetSymbol()

const G4String & G4Element::GetSymbol ( ) const
inline

Definition at line 128 of file G4Element.hh.

128{return fSymbol;}

Referenced by G4tgbGeometryDumper::DumpElement().

◆ GetZ()

G4double G4Element::GetZ ( ) const
inline

Definition at line 131 of file G4Element.hh.

131{return fZeff;}

Referenced by G4VCrossSectionHandler::ActiveElements(), G4AdjointPhotoElectricModel::AdjointCrossSectionPerAtom(), G4NeutronHPThermalScattering::ApplyYourself(), G4MuonMinusCaptureAtRest::AtRestDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4AugerData::BuildAugerTransitionTable(), G4AdjointCSManager::BuildCrossSectionMatrices(), G4NeutronHPJENDLHEData::BuildPhysicsTable(), G4Nucleus::ChooseParameters(), G4AdjointCSManager::ComputeAdjointCS(), G4VCrossSectionDataSet::ComputeCrossSection(), G4VEmModel::ComputeCrossSectionPerAtom(), G4EmCalculator::ComputeCrossSectionPerAtom(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4PAIySection::ComputeLowEnergyCof(), G4LENDCrossSection::create_used_target_map(), G4LENDModel::create_used_target_map(), G4tgbGeometryDumper::DumpElement(), G4GDMLWriteMaterials::ElementWrite(), G4CrossSectionDataStore::GetCrossSection(), G4NeutronHPCaptureData::GetCrossSection(), G4NeutronHPElasticData::GetCrossSection(), G4NeutronHPFissionData::GetCrossSection(), G4NeutronHPInelasticData::GetCrossSection(), G4NeutronHPJENDLHEData::GetCrossSection(), G4NeutronHPorLCaptureData::GetCrossSection(), G4NeutronHPorLEInelasticData::GetCrossSection(), G4NeutronHPorLElasticData::GetCrossSection(), G4NeutronHPorLFissionData::GetCrossSection(), G4GammaConversionToMuons::GetCrossSectionPerAtom(), GVFlashShowerParameterisation::GetEffZ(), G4VComponentCrossSection::GetElasticElementCrossSection(), G4ChargeExchangeProcess::GetElementCrossSection(), G4ComponentGGHadronNucleusXsc::GetHadronNucleonXsc(), G4ComponentGGNuclNuclXsc::GetHadronNucleonXsc(), G4GGNuclNuclCrossSection::GetHadronNucleonXsc(), G4GlauberGribovCrossSection::GetHadronNucleonXsc(), G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscNS(), G4GlauberGribovCrossSection::GetHadronNucleonXscNS(), G4ComponentGGHadronNucleusXsc::GetHadronNucleonXscPDG(), G4GlauberGribovCrossSection::GetHadronNucleonXscPDG(), G4ComponentGGHadronNucleusXsc::GetHNinelasticXsc(), G4GlauberGribovCrossSection::GetHNinelasticXsc(), G4VComponentCrossSection::GetInelasticElementCrossSection(), G4NeutronIsotopeProduction::GetIsotope(), G4QCoherentChargeExchange::GetMeanFreePath(), G4QDiffraction::GetMeanFreePath(), G4QElastic::GetMeanFreePath(), G4QInelastic::GetMeanFreePath(), G4QIonIonElastic::GetMeanFreePath(), G4QLowEnergy::GetMeanFreePath(), G4QNGamma::GetMeanFreePath(), G4NeutronHPThermalBoost::GetThermalEnergy(), G4VComponentCrossSection::GetTotalElementCrossSection(), G4NeutronHPElementData::Init(), G4HadronNucleonXsc::IsApplicable(), G4EmCorrections::NuclearDEDX(), G4GammaConversionToMuons::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4NeutronHPChannel::Register(), G4eBremParametrizedModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4MuPairProductionModel::SampleSecondaries(), G4eBremsstrahlungModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4XrayRayleighModel::SampleSecondaries(), G4CrossSectionDataStore::SampleZandA(), G4NeutronCaptureXS::SelectIsotope(), G4ElementSelector::SelectZandA(), and G4hQAOModel::StoppingPower().

◆ increaseCountUse()

void G4Element::increaseCountUse ( )
inline

Definition at line 192 of file G4Element.hh.

192{fCountUse++;}

Referenced by G4Material::AddElement(), and G4Material::AddMaterial().

◆ operator!=()

G4int G4Element::operator!= ( const G4Element right) const

Definition at line 496 of file G4Element.cc.

497{
498 return (this != (G4Element*) &right);
499}

◆ operator==()

G4int G4Element::operator== ( const G4Element right) const

Definition at line 489 of file G4Element.cc.

490{
491 return (this == (G4Element*) &right);
492}

◆ SetName()

void G4Element::SetName ( const G4String name)
inline

Definition at line 227 of file G4Element.hh.

227{fName=name;}

Referenced by G4GDMLRead::StripNames().

◆ SetNaturalAbandancesFlag()

void G4Element::SetNaturalAbandancesFlag ( G4bool  val)
inline

Definition at line 281 of file G4Element.hh.

282{
283 fNaturalAbandances = val;
284}

Friends And Related Function Documentation

◆ operator<< [1/3]

std::ostream & operator<< ( std::ostream &  flux,
G4Element element 
)
friend

Definition at line 529 of file G4Element.cc.

530{
531 flux << &element;
532 return flux;
533}

◆ operator<< [2/3]

std::ostream & operator<< ( std::ostream &  flux,
G4Element element 
)
friend

Definition at line 503 of file G4Element.cc.

504{
505 std::ios::fmtflags mode = flux.flags();
506 flux.setf(std::ios::fixed,std::ios::floatfield);
507 G4long prec = flux.precision(3);
508
509 flux
510 << " Element: " << element->fName << " (" << element->fSymbol << ")"
511 << " Z = " << std::setw(4) << std::setprecision(1) << element->fZeff
512 << " N = " << std::setw(5) << std::setprecision(1) << element->fNeff
513 << " A = " << std::setw(6) << std::setprecision(2)
514 << (element->fAeff)/(g/mole) << " g/mole";
515
516 for (size_t i=0; i<element->fNumberOfIsotopes; i++)
517 flux
518 << "\n ---> " << (*(element->theIsotopeVector))[i]
519 << " abundance: " << std::setw(6) << std::setprecision(2)
520 << (element->fRelativeAbundanceVector[i])/perCent << " %";
521
522 flux.precision(prec);
523 flux.setf(mode,std::ios::floatfield);
524 return flux;
525}
long G4long
Definition: G4Types.hh:68

◆ operator<< [3/3]

std::ostream & operator<< ( std::ostream &  flux,
G4ElementTable  ElementTable 
)
friend

Definition at line 537 of file G4Element.cc.

538{
539 //Dump info for all known elements
540 flux << "\n***** Table : Nb of elements = " << ElementTable.size()
541 << " *****\n" << G4endl;
542
543 for (size_t i=0; i<ElementTable.size(); i++) flux << ElementTable[i]
544 << G4endl << G4endl;
545
546 return flux;
547}

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