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

#include <G4Material.hh>

+ Inheritance diagram for G4Material:

Public Member Functions

 G4Material (const G4String &name, G4double z, G4double a, G4double density, G4State state=kStateUndefined, G4double temp=NTP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
 G4Material (const G4String &name, G4double density, G4int nComponents, G4State state=kStateUndefined, G4double temp=NTP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
 G4Material (const G4String &name, G4double density, const G4Material *baseMaterial, G4State state=kStateUndefined, G4double temp=NTP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
virtual ~G4Material ()
 
void SetChemicalFormula (const G4String &chF)
 
void SetFreeElectronDensity (G4double val)
 
void ComputeDensityEffectOnFly (G4bool val)
 
 G4Material (const G4Material &)=delete
 
const G4Materialoperator= (const G4Material &)=delete
 
void AddElementByNumberOfAtoms (const G4Element *elm, G4int nAtoms)
 
void AddElement (G4Element *elm, G4int nAtoms)
 
void AddElementByMassFraction (const G4Element *elm, G4double fraction)
 
void AddElement (G4Element *elm, G4double frac)
 
void AddMaterial (G4Material *material, G4double fraction)
 
const G4StringGetName () const
 
const G4StringGetChemicalFormula () const
 
G4double GetFreeElectronDensity () const
 
G4double GetDensity () const
 
G4State GetState () const
 
G4double GetTemperature () const
 
G4double GetPressure () const
 
std::size_t GetNumberOfElements () const
 
const G4ElementVectorGetElementVector () const
 
const G4doubleGetFractionVector () const
 
const G4intGetAtomsVector () const
 
const G4ElementGetElement (G4int iel) const
 
const G4doubleGetVecNbOfAtomsPerVolume () const
 
G4double GetTotNbOfAtomsPerVolume () const
 
G4double GetTotNbOfElectPerVolume () const
 
const G4doubleGetAtomicNumDensityVector () const
 
G4double GetElectronDensity () const
 
G4double GetRadlen () const
 
G4double GetNuclearInterLength () const
 
G4IonisParamMatGetIonisation () const
 
G4SandiaTableGetSandiaTable () const
 
const G4MaterialGetBaseMaterial () const
 
const std::map< G4Material *, G4double > & GetMatComponents () const
 
G4double GetMassOfMolecule () const
 
G4double GetZ () const
 
G4double GetA () const
 
void SetMaterialPropertiesTable (G4MaterialPropertiesTable *anMPT)
 
G4MaterialPropertiesTableGetMaterialPropertiesTable () const
 
std::size_t GetIndex () const
 
void SetName (const G4String &name)
 
virtual G4bool IsExtended () const
 
G4bool operator== (const G4Material &) const =delete
 
G4bool operator!= (const G4Material &) const =delete
 

Static Public Member Functions

static G4MaterialTableGetMaterialTable ()
 
static std::size_t GetNumberOfMaterials ()
 
static G4MaterialGetMaterial (const G4String &name, G4bool warning=true)
 
static G4MaterialGetMaterial (G4double z, G4double a, G4double dens)
 
static G4MaterialGetMaterial (std::size_t nComp, G4double dens)
 

Friends

std::ostream & operator<< (std::ostream &flux, const G4Material *material)
 
std::ostream & operator<< (std::ostream &flux, const G4Material &material)
 
std::ostream & operator<< (std::ostream &flux, const G4MaterialTable &MaterialTable)
 

Detailed Description

Definition at line 116 of file G4Material.hh.

Constructor & Destructor Documentation

◆ G4Material() [1/4]

G4Material::G4Material ( const G4String & name,
G4double z,
G4double a,
G4double density,
G4State state = kStateUndefined,
G4double temp = NTP_Temperature,
G4double pressure = CLHEP::STP_Pressure )

Definition at line 85 of file G4Material.cc.

87 : fName(name)
88{
89 InitializePointers();
90
91 if (density < universe_mean_density) {
92 G4cout << " G4Material WARNING:"
93 << " define a material with density=0 is not allowed. \n"
94 << " The material " << name << " will be constructed with the"
95 << " default minimal density: " << universe_mean_density / (g / cm3) << "g/cm3"
96 << G4endl;
97 density = universe_mean_density;
98 }
99
100 fDensity = density;
101 fState = state;
102 fTemp = temp;
103 fPressure = pressure;
104
105 // Initialize theElementVector allocating one
106 // element corresponding to this material
107 fNbComponents = fNumberOfElements = 1;
108 theElementVector = new G4ElementVector();
109
110 // take element from DB
111 G4NistManager* nist = G4NistManager::Instance();
112 G4int iz = G4lrint(z);
113 auto elm = nist->FindOrBuildElement(iz);
114 if (elm == nullptr) {
115 elm = new G4Element("ELM_" + name, name, z, a);
116 }
117 theElementVector->push_back(elm);
118
119 fMassFractionVector = new G4double[1];
120 fMassFractionVector[0] = 1.;
121 fMassOfMolecule = a / CLHEP::Avogadro;
122
123 if (fState == kStateUndefined) {
124 if (fDensity > kGasThreshold) {
125 fState = kStateSolid;
126 }
127 else {
128 fState = kStateGas;
129 }
130 }
131
132 ComputeDerivedQuantities();
133}
std::vector< const G4Element * > G4ElementVector
@ kStateSolid
@ kStateGas
@ kStateUndefined
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4NistManager * Instance()
G4Element * FindOrBuildElement(G4int Z, G4bool isotopes=true)
const char * name(G4int ptype)
int G4lrint(double ad)
Definition templates.hh:134

Referenced by AddMaterial(), G4ExtendedMaterial::G4ExtendedMaterial(), G4ExtendedMaterial::G4ExtendedMaterial(), G4ExtendedMaterial::G4ExtendedMaterial(), G4ExtendedMaterial::G4ExtendedMaterial(), G4Material(), G4Material(), GetBaseMaterial(), GetMaterial(), GetMaterial(), GetMaterial(), operator!=(), operator<<, operator<<, operator=(), and operator==().

◆ G4Material() [2/4]

G4Material::G4Material ( const G4String & name,
G4double density,
G4int nComponents,
G4State state = kStateUndefined,
G4double temp = NTP_Temperature,
G4double pressure = CLHEP::STP_Pressure )

Definition at line 140 of file G4Material.cc.

142 : fName(name)
143{
144 InitializePointers();
145
146 if (density < universe_mean_density) {
147 G4cout << "--- Warning from G4Material::G4Material()"
148 << " define a material with density=0 is not allowed. \n"
149 << " The material " << name << " will be constructed with the"
150 << " default minimal density: " << universe_mean_density / (g / cm3) << "g/cm3"
151 << G4endl;
152 density = universe_mean_density;
153 }
154
155 fDensity = density;
156 fState = state;
157 fTemp = temp;
158 fPressure = pressure;
159
160 fNbComponents = nComponents;
161 fMassFraction = true;
162
163 if (fState == kStateUndefined) {
164 if (fDensity > kGasThreshold) {
165 fState = kStateSolid;
166 }
167 else {
168 fState = kStateGas;
169 }
170 }
171}

◆ G4Material() [3/4]

G4Material::G4Material ( const G4String & name,
G4double density,
const G4Material * baseMaterial,
G4State state = kStateUndefined,
G4double temp = NTP_Temperature,
G4double pressure = CLHEP::STP_Pressure )

Definition at line 177 of file G4Material.cc.

179 : fName(name)
180{
181 InitializePointers();
182
183 if (density < universe_mean_density) {
184 G4cout << "--- Warning from G4Material::G4Material()"
185 << " define a material with density=0 is not allowed. \n"
186 << " The material " << name << " will be constructed with the"
187 << " default minimal density: " << universe_mean_density / (g / cm3) << "g/cm3"
188 << G4endl;
189 density = universe_mean_density;
190 }
191
192 fDensity = density;
193 fState = state;
194 fTemp = temp;
195 fPressure = pressure;
196
197 fBaseMaterial = bmat;
198 auto ptr = bmat;
199 if (nullptr != ptr) {
200 while (true) {
201 ptr = ptr->GetBaseMaterial();
202 if (nullptr == ptr) {
203 break;
204 }
205 fBaseMaterial = ptr;
206 }
207 }
208
209 fChemicalFormula = fBaseMaterial->GetChemicalFormula();
210 fMassOfMolecule = fBaseMaterial->GetMassOfMolecule();
211
212 fNumberOfElements = (G4int)fBaseMaterial->GetNumberOfElements();
213 fNbComponents = fNumberOfElements;
214
215 CopyPointersOfBaseMaterial();
216}

◆ ~G4Material()

G4Material::~G4Material ( )
virtual

Definition at line 220 of file G4Material.cc.

221{
222 if (fBaseMaterial == nullptr) {
223 delete theElementVector;
224 delete fSandiaTable;
225 delete[] fMassFractionVector;
226 delete[] fAtomsVector;
227 }
228 delete fIonisation;
229 delete[] fVecNbOfAtomsPerVolume;
230
231 // Remove this material from the MaterialTable.
232 //
233 (*GetMaterialTable())[fIndexInTable] = nullptr;
234}
static G4MaterialTable * GetMaterialTable()

◆ G4Material() [4/4]

G4Material::G4Material ( const G4Material & )
delete

Member Function Documentation

◆ AddElement() [1/2]

void G4Material::AddElement ( G4Element * elm,
G4double frac )
inline

Definition at line 164 of file G4Material.hh.

164{ AddElementByMassFraction(elm, frac); }
void AddElementByMassFraction(const G4Element *elm, G4double fraction)

◆ AddElement() [2/2]

void G4Material::AddElement ( G4Element * elm,
G4int nAtoms )
inline

◆ AddElementByMassFraction()

void G4Material::AddElementByMassFraction ( const G4Element * elm,
G4double fraction )

Definition at line 423 of file G4Material.cc.

424{
425 // perform checks consistency
426 if (fraction < 0.0 || fraction > 1.0) {
428 ed << "For material " << fName << " and added element " << elm->GetName()
429 << " massFraction= " << fraction << " is wrong ";
430 G4Exception("G4Material::AddElementByMassFraction()", "mat031", FatalException, ed, "");
431 }
432 if (! fMassFraction) {
434 ed << "For material " << fName << " and added element " << elm->GetName()
435 << ", massFraction= " << fraction << ", fIdxComponent=" << fIdxComponent
436 << " problem: cannot add by mass fraction after "
437 << "addition of elements by number of atoms";
438 G4Exception("G4Material::AddElementByMassFraction()", "mat031", FatalException, ed, "");
439 }
440 if (fIdxComponent >= fNbComponents) {
442 ed << "For material " << fName << " and added element " << elm->GetName()
443 << ", massFraction= " << fraction << ", fIdxComponent=" << fIdxComponent
444 << "; attempt to add more than the declared number of components " << fIdxComponent
445 << " >= " << fNbComponents;
446 G4Exception("G4Material::AddElementByMassFraction()", "mat031", FatalException, ed, "");
447 }
448 if (0 == fIdxComponent) {
449 fElmFrac = new std::vector<G4double>;
450 fElm = new std::vector<const G4Element*>;
451 }
452
453 // filling
454 G4bool isAdded = false;
455 if (! fElm->empty()) {
456 for (G4int i = 0; i < fNumberOfElements; ++i) {
457 if (elm == (*fElm)[i]) {
458 (*fElmFrac)[i] += fraction;
459 isAdded = true;
460 break;
461 }
462 }
463 }
464 if (! isAdded) {
465 fElm->push_back(elm);
466 fElmFrac->push_back(fraction);
467 ++fNumberOfElements;
468 }
469 ++fIdxComponent;
470
471 // is filled
472 if (fIdxComponent == fNbComponents) {
473 FillVectors();
474 }
475}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
bool G4bool
Definition G4Types.hh:86
const G4String & GetName() const
Definition G4Element.hh:115

Referenced by AddElement().

◆ AddElementByNumberOfAtoms()

void G4Material::AddElementByNumberOfAtoms ( const G4Element * elm,
G4int nAtoms )

Definition at line 348 of file G4Material.cc.

349{
350 // perform checks consistency
351 if (0 == fIdxComponent) {
352 fMassFraction = false;
353 fAtoms = new std::vector<G4int>;
354 fElm = new std::vector<const G4Element*>;
355 }
356 if (fIdxComponent >= fNbComponents) {
358 ed << "For material " << fName << " and added element " << elm->GetName()
359 << " with Natoms=" << nAtoms
360 << " wrong attempt to add more than the declared number of elements " << fIdxComponent
361 << " >= " << fNbComponents;
362 G4Exception("G4Material::AddElementByNumberOfAtoms()", "mat031", FatalException, ed, "");
363 }
364 if (fMassFraction) {
366 ed << "For material " << fName << " and added element " << elm->GetName()
367 << " with Natoms=" << nAtoms << " problem: cannot add by number of atoms after "
368 << "addition of elements by mass fraction";
369 G4Exception("G4Material::AddElementByNumberOfAtoms()", "mat031", FatalException, ed, "");
370 }
371 if (0 >= nAtoms) {
373 ed << "For material " << fName << " and added element " << elm->GetName()
374 << " with Natoms=" << nAtoms << " problem: number of atoms should be above zero";
375 G4Exception("G4Material::AddElementByNumberOfAtoms()", "mat031", FatalException, ed, "");
376 }
377
378 // filling
379 G4bool isAdded = false;
380 if (! fElm->empty()) {
381 for (G4int i = 0; i < fNumberOfElements; ++i) {
382 if (elm == (*fElm)[i]) {
383 (*fAtoms)[i] += nAtoms;
384 isAdded = true;
385 break;
386 }
387 }
388 }
389 if (! isAdded) {
390 fElm->push_back(elm);
391 fAtoms->push_back(nAtoms);
392 ++fNumberOfElements;
393 }
394 ++fIdxComponent;
395
396 // is filled - complete composition of atoms
397 if (fIdxComponent == fNbComponents) {
398 theElementVector = new G4ElementVector();
399 theElementVector->reserve(fNumberOfElements);
400 fAtomsVector = new G4int[fNumberOfElements];
401 fMassFractionVector = new G4double[fNumberOfElements];
402
403 G4double Amol = 0.;
404 for (G4int i = 0; i < fNumberOfElements; ++i) {
405 theElementVector->push_back((*fElm)[i]);
406 fAtomsVector[i] = (*fAtoms)[i];
407 G4double w = fAtomsVector[i] * (*fElm)[i]->GetA();
408 Amol += w;
409 fMassFractionVector[i] = w;
410 }
411 for (G4int i = 0; i < fNumberOfElements; ++i) {
412 fMassFractionVector[i] /= Amol;
413 }
414 delete fAtoms;
415 delete fElm;
416 fMassOfMolecule = Amol / CLHEP::Avogadro;
417 ComputeDerivedQuantities();
418 }
419}

Referenced by AddElement().

◆ AddMaterial()

void G4Material::AddMaterial ( G4Material * material,
G4double fraction )

Definition at line 480 of file G4Material.cc.

481{
482 if (fraction < 0.0 || fraction > 1.0) {
484 ed << "For material " << fName << " and added material " << material->GetName()
485 << ", massFraction= " << fraction << " is wrong ";
486 G4Exception("G4Material::AddMaterial()", "mat031", FatalException, ed, "");
487 }
488 if (! fMassFraction) {
490 ed << "For material " << fName << " and added material " << material->GetName()
491 << ", massFraction= " << fraction << ", fIdxComponent=" << fIdxComponent
492 << " problem: cannot add by mass fraction after "
493 << "addition of elements by number of atoms";
494 G4Exception("G4Material::AddMaterial()", "mat031", FatalException, ed, "");
495 }
496 if (fIdxComponent >= fNbComponents) {
498 ed << "For material " << fName << " and added material " << material->GetName()
499 << ", massFraction= " << fraction
500 << "; attempt to add more than the declared number of components " << fIdxComponent
501 << " >= " << fNbComponents;
502 G4Exception("G4Material::AddMaterial()", "mat031", FatalException, ed, "");
503 }
504 if (0 == fIdxComponent) {
505 fElmFrac = new std::vector<G4double>;
506 fElm = new std::vector<const G4Element*>;
507 }
508
509 // filling
510 auto nelm = (G4int)material->GetNumberOfElements();
511 for (G4int j = 0; j < nelm; ++j) {
512 auto elm = material->GetElement(j);
513 auto frac = material->GetFractionVector();
514 G4bool isAdded = false;
515 if (! fElm->empty()) {
516 for (G4int i = 0; i < fNumberOfElements; ++i) {
517 if (elm == (*fElm)[i]) {
518 (*fElmFrac)[i] += fraction * frac[j];
519 isAdded = true;
520 break;
521 }
522 }
523 }
524 if (! isAdded) {
525 fElm->push_back(elm);
526 fElmFrac->push_back(fraction * frac[j]);
527 ++fNumberOfElements;
528 }
529 }
530
531 fMatComponents[material] = fraction;
532 ++fIdxComponent;
533
534 // is filled
535 if (fIdxComponent == fNbComponents) {
536 FillVectors();
537 }
538}
const G4Element * GetElement(G4int iel) const
const G4double * GetFractionVector() const
std::size_t GetNumberOfElements() const
const G4String & GetName() const

Referenced by G4tgbMaterialMixtureByVolume::BuildG4Material(), G4tgbMaterialMixtureByWeight::BuildG4Material(), and G4GDMLReadMaterials::MixtureRead().

◆ ComputeDensityEffectOnFly()

void G4Material::ComputeDensityEffectOnFly ( G4bool val)

Definition at line 632 of file G4Material.cc.

633{
634 if (! IsLocked()) {
635 if (nullptr == fIonisation) {
636 fIonisation = new G4IonisParamMat(this);
637 }
638 fIonisation->ComputeDensityEffectOnFly(val);
639 }
640}

Referenced by G4NistManager::SetDensityEffectCalculatorFlag().

◆ GetA()

G4double G4Material::GetA ( ) const

Definition at line 723 of file G4Material.cc.

724{
725 if (fNumberOfElements > 1) {
727 ed << "For material " << fName << " ERROR in GetA() - Nelm=" << fNumberOfElements
728 << " > 1, which is not allowed";
729 G4Exception("G4Material::GetA()", "mat036", FatalException, ed, "");
730 }
731 return (*theElementVector)[0]->GetA();
732}

Referenced by G4tgbGeometryDumper::DumpMaterial(), G4DNAQuinnPlasmonExcitationModel::GetCrossSection(), GVFlashShowerParameterisation::GetEffA(), G4GDMLWriteMaterials::MaterialWrite(), and G4DNAQuinnPlasmonExcitationModel::SampleSecondaries().

◆ GetAtomicNumDensityVector()

◆ GetAtomsVector()

const G4int * G4Material::GetAtomsVector ( ) const
inline

Definition at line 189 of file G4Material.hh.

189{ return fAtomsVector; }

Referenced by G4VLEPTSModel::ReadParam().

◆ GetBaseMaterial()

◆ GetChemicalFormula()

const G4String & G4Material::GetChemicalFormula ( ) const
inline

◆ GetDensity()

G4double G4Material::GetDensity ( ) const
inline

◆ GetElectronDensity()

G4double G4Material::GetElectronDensity ( ) const
inline

Definition at line 203 of file G4Material.hh.

203{ return fTotNbOfElectPerVolume; }

Referenced by G4ForwardXrayTR::BuildXrayTRtables(), G4BetheBlochModel::ComputeDEDXPerVolume(), G4BraggIonModel::ComputeDEDXPerVolume(), G4BraggModel::ComputeDEDXPerVolume(), G4ICRU73QOModel::ComputeDEDXPerVolume(), G4MollerBhabhaModel::ComputeDEDXPerVolume(), G4MuBetheBlochModel::ComputeDEDXPerVolume(), G4LindhardSorensenIonModel::CorrectionsAlongStep(), G4AtimaEnergyLossModel::CrossSectionPerVolume(), G4BetheBlochModel::CrossSectionPerVolume(), G4BraggIonModel::CrossSectionPerVolume(), G4BraggModel::CrossSectionPerVolume(), G4eeToHadronsModel::CrossSectionPerVolume(), G4eeToHadronsMultiModel::CrossSectionPerVolume(), G4eeToTwoGammaModel::CrossSectionPerVolume(), G4eplusTo2or3GammaModel::CrossSectionPerVolume(), G4eplusTo3GammaOKVIModel::CrossSectionPerVolume(), G4ICRU73QOModel::CrossSectionPerVolume(), G4LindhardSorensenIonModel::CrossSectionPerVolume(), G4MollerBhabhaModel::CrossSectionPerVolume(), G4MuBetheBlochModel::CrossSectionPerVolume(), G4DynamicParticleFluctuation::Dispersion(), G4IonFluctuations::Dispersion(), G4mplIonisationModel::Dispersion(), G4mplIonisationWithDeltaModel::Dispersion(), G4PAIModel::Dispersion(), G4PAIPhotModel::Dispersion(), G4UniversalFluctuation::Dispersion(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4mplIonisationModel::Initialise(), G4mplIonisationWithDeltaModel::Initialise(), G4PAIxSection::Initialize(), G4PAIySection::Initialize(), G4DynamicParticleFluctuation::SampleFluctuations(), G4UniversalFluctuation::SampleFluctuations(), G4eBremParametrizedModel::SetupForMaterial(), G4eBremsstrahlungRelModel::SetupForMaterial(), and G4SeltzerBergerModel::SetupForMaterial().

◆ GetElement()

◆ GetElementVector()

const G4ElementVector * G4Material::GetElementVector ( ) const
inline

Definition at line 183 of file G4Material.hh.

183{ return theElementVector; }

Referenced by G4VCrossSectionHandler::ActiveElements(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4CrossSectionHandler::BuildCrossSectionsForMaterials(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4PenelopeBremsstrahlungFS::BuildScaledXSTable(), G4Nucleus::ChooseParameters(), G4eBremParametrizedModel::ComputeDEDXPerVolume(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4LivermoreIonisationModel::ComputeDEDXPerVolume(), G4MuBremsstrahlungModel::ComputeDEDXPerVolume(), G4MuPairProductionModel::ComputeDEDXPerVolume(), G4RiGeMuPairProductionModel::ComputeDEDXPerVolume(), G4SeltzerBergerModel::ComputeDEDXPerVolume(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4tgbGeometryDumper::DumpMaterial(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(), G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(), G4IonICRU73Data::GetDEDX(), G4HadronicProcessStore::GetElasticCrossSectionPerVolume(), G4ElNeutrinoNucleusTotXsc::GetElementCrossSection(), G4MuNeutrinoNucleusTotXsc::GetElementCrossSection(), G4TauNeutrinoNucleusTotXsc::GetElementCrossSection(), G4HadronicProcessStore::GetFissionCrossSectionPerVolume(), G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(), G4BoldyshevTripletModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAQuinnPlasmonExcitationModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4IonICRU73Data::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4PenelopeRayleighModel::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4ElasticHadrNucleusHE::InitialiseModel(), G4eCoulombScatteringModel::MinPrimaryEnergy(), G4hCoulombScatteringModel::MinPrimaryEnergy(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PixeCrossSectionHandler::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomElement(), G4ElementSelector::SelectZandA(), G4MSSteppingAction::UserSteppingAction(), G4PixeCrossSectionHandler::ValueForMaterial(), and G4VCrossSectionHandler::ValueForMaterial().

◆ GetFractionVector()

◆ GetFreeElectronDensity()

G4double G4Material::GetFreeElectronDensity ( ) const
inline

Definition at line 173 of file G4Material.hh.

173{ return fFreeElecDensity; }

◆ GetIndex()

std::size_t G4Material::GetIndex ( ) const
inline

Definition at line 239 of file G4Material.hh.

239{ return fIndexInTable; }

Referenced by G4VLEPTSModel::BuildMeanFreePathTable(), G4VLEPTSModel::BuildPhysicsTable(), G4AdjointCSManager::ComputeAdjointCS(), G4EnergyLossForExtrapolator::ComputeDEDX(), G4BetheBlochModel::ComputeDEDXPerVolume(), G4EnergyLossForExtrapolator::ComputeEnergy(), G4EnergyLossForExtrapolator::ComputeRange(), G4DNABornExcitationModel1::CrossSectionPerVolume(), G4DNABornExcitationModel2::CrossSectionPerVolume(), G4DNABornIonisationModel1::CrossSectionPerVolume(), G4DNABornIonisationModel2::CrossSectionPerVolume(), G4DNAChampionElasticModel::CrossSectionPerVolume(), G4DNACPA100ElasticModel::CrossSectionPerVolume(), G4DNACPA100ExcitationModel::CrossSectionPerVolume(), G4DNACPA100IonisationModel::CrossSectionPerVolume(), G4DNADingfelderChargeDecreaseModel::CrossSectionPerVolume(), G4DNADingfelderChargeIncreaseModel::CrossSectionPerVolume(), G4DNADoubleIonisationModel::CrossSectionPerVolume(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4DNAEmfietzoglouIonisationModel::CrossSectionPerVolume(), G4DNAIonElasticModel::CrossSectionPerVolume(), G4DNAMeltonAttachmentModel::CrossSectionPerVolume(), G4DNAMillerGreenExcitationModel::CrossSectionPerVolume(), G4DNAModelInterface::CrossSectionPerVolume(), G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), G4DNAQuadrupleIonisationModel::CrossSectionPerVolume(), G4DNARPWBAExcitationModel::CrossSectionPerVolume(), G4DNARPWBAIonisationModel::CrossSectionPerVolume(), G4DNARuddIonisationExtendedModel::CrossSectionPerVolume(), G4DNARuddIonisationModel::CrossSectionPerVolume(), G4DNASancheExcitationModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNATransformElectronModel::CrossSectionPerVolume(), G4DNATripleIonisationModel::CrossSectionPerVolume(), G4DNAUeharaScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4GoudsmitSaundersonMscModel::CrossSectionPerVolume(), G4PEEffectFluoModel::CrossSectionPerVolume(), G4eCoulombScatteringModel::DefineMaterial(), G4hCoulombScatteringModel::DefineMaterial(), G4DNABrownianTransportation::Diffusion(), G4PAIxSection::G4PAIxSection(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4EnergyLossTables::GetDEDX(), G4IonICRU73Data::GetDEDX(), G4EnergyLossTables::GetDeltaLabTime(), G4EnergyLossTables::GetDeltaProperTime(), G4EnergyLossTables::GetLabTime(), G4VLEPTSModel::GetMeanFreePath(), G4DNAMolecularMaterial::GetMolecularConfiguration(), G4EnergyLossTables::GetPreciseDEDX(), G4EnergyLossTables::GetPreciseEnergyFromRange(), G4EnergyLossTables::GetPreciseRangeFromEnergy(), G4EnergyLossTables::GetProperTime(), G4EnergyLossTables::GetRange(), G4GoudsmitSaundersonMscModel::GetTransportMeanFreePath(), G4IonICRU73Data::Initialise(), G4GoudsmitSaundersonTable::InitSCPCorrection(), G4OpWLS2::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4Scintillation::PostStepDoIt(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4HadronicProcess::PostStepGetPhysicalInteractionLength(), G4DNAMolecularMaterial::RecordMolecularMaterial(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAModelInterface::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), G4DNAPTBIonisationModel::SampleSecondaries(), G4PixeCrossSectionHandler::SelectRandomAtom(), G4DNAMolecularMaterial::SetMolecularConfiguration(), G4DNAMolecularMaterial::SetMolecularConfiguration(), G4EnergyLossForExtrapolator::TrueStepLength(), and G4EmSaturation::VisibleEnergyDeposition().

◆ GetIonisation()

G4IonisParamMat * G4Material::GetIonisation ( ) const
inline

Definition at line 212 of file G4Material.hh.

212{ return fIonisation; }

Referenced by G4hImpactIonisation::BuildPhysicsTable(), G4AtimaEnergyLossModel::ComputeDEDXPerVolume(), G4BetheBlochModel::ComputeDEDXPerVolume(), G4MollerBhabhaModel::ComputeDEDXPerVolume(), G4MuBetheBlochModel::ComputeDEDXPerVolume(), G4AtimaFluctuations::Dispersion(), G4IonFluctuations::Dispersion(), G4EmSaturation::DumpBirksCoefficients(), G4tgbGeometryDumper::DumpMaterial(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), G4ionEffectiveCharge::EffectiveCharge(), G4ElectronIonPair::FindG4MeanEnergyPerIonPair(), G4tgbMaterialMgr::FindOrBuildG4Material(), G4eDPWAElasticDCS::InitSCPCorrection(), G4GoudsmitSaundersonTable::InitSCPCorrection(), G4hBetheBlochModel::LowEnergyLimit(), G4GDMLReadMaterials::MaterialRead(), G4GDMLWriteMaterials::MaterialWrite(), G4ElectronIonPair::MeanNumberOfIonsAlongStep(), G4AtimaEnergyLossModel::MinEnergyCut(), G4BetheBlochModel::MinEnergyCut(), G4BraggModel::MinEnergyCut(), G4IonParametrisedLossModel::MinEnergyCut(), G4LindhardSorensenIonModel::MinEnergyCut(), G4mplIonisationWithDeltaModel::MinEnergyCut(), G4MuBetheBlochModel::MinEnergyCut(), G4CoulombScattering::MinPrimaryEnergy(), operator<<, G4hImpactIonisation::PrintInfoDefinition(), G4DynamicParticleFluctuation::SampleFluctuations(), G4UniversalFluctuation::SampleFluctuations(), G4UrbanFluctuation::SampleGlandz(), G4AllisonPositronAtRestModel::SampleSecondaries(), G4IonisParamMat::SetDensityEffectParameters(), G4ChannelingFastSimCrystalData::SetMaterialProperties(), G4WentzelOKandVIxSection::SetupKinematic(), G4WentzelVIRelXSection::SetupKinematic(), and G4EmSaturation::VisibleEnergyDeposition().

◆ GetMassOfMolecule()

G4double G4Material::GetMassOfMolecule ( ) const
inline

◆ GetMatComponents()

const std::map< G4Material *, G4double > & G4Material::GetMatComponents ( ) const
inline

◆ GetMaterial() [1/3]

G4Material * G4Material::GetMaterial ( const G4String & name,
G4bool warning = true )
static

Definition at line 663 of file G4Material.cc.

664{
665 // search the material by its name
666 for (auto const & j : *GetMaterialTable()) {
667 if (j->GetName() == materialName) {
668 return j;
669 }
670 }
671
672 // the material does not exist in the table
673 if (warn) {
674 G4cout << "G4Material::GetMaterial() WARNING: The material: " << materialName
675 << " does not exist in the table. Return NULL pointer." << G4endl;
676 }
677 return nullptr;
678}

Referenced by G4DNABrownianTransportation::BuildPhysicsTable(), G4LossTableBuilder::BuildTableForModel(), G4Track::CalculateVelocityForOpticalPhoton(), G4ProductionCutsTable::CheckMaterialInfo(), G4DNAMolecularDissociation::DecayIt(), G4EmCalculator::FindMaterial(), G4DNACPA100ElasticModel::G4DNACPA100ElasticModel(), G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel(), G4DNACPA100ExcitationStructure::G4DNACPA100ExcitationStructure(), G4DNACPA100IonisationModel::G4DNACPA100IonisationModel(), G4DNACPA100IonisationStructure::G4DNACPA100IonisationStructure(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4DNAPTBElasticModel::G4DNAPTBElasticModel(), G4DNAPTBExcitationModel::G4DNAPTBExcitationModel(), G4DNAPTBExcitationStructure::G4DNAPTBExcitationStructure(), G4DNAPTBIonisationModel::G4DNAPTBIonisationModel(), G4DNAPTBIonisationStructure::G4DNAPTBIonisationStructure(), G4GDMLReadMaterials::GetMaterial(), G4XrayReflection::GetMeanFreePath(), G4GoudsmitSaundersonMscModel::GetTransportMeanFreePath(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNAChampionElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNADoubleIonisationModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNAIonElasticModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNAModelInterface::Initialise(), G4DNAQuadrupleIonisationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4DNARPWBAExcitationModel::Initialise(), G4DNARPWBAIonisationModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNATransformElectronModel::Initialise(), G4DNATripleIonisationModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4DNAVacuumModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4TDNAOneStepThermalizationModel< MODEL >::Initialise(), G4VDNAModel::LoadCrossSectionData(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAMolecularMaterial::SetMolecularConfiguration(), G4EnergySplitter::SplitEnergyInVolumes(), and G4ProductionCutsTable::UpdateCoupleTable().

◆ GetMaterial() [2/3]

G4Material * G4Material::GetMaterial ( G4double z,
G4double a,
G4double dens )
static

Definition at line 682 of file G4Material.cc.

683{
684 // search the material by its name
685 for (auto const & mat : *GetMaterialTable()) {
686 if (1 == mat->GetNumberOfElements() && z == mat->GetZ() && a == mat->GetA() &&
687 dens == mat->GetDensity())
688 {
689 return mat;
690 }
691 }
692 return nullptr;
693}

◆ GetMaterial() [3/3]

G4Material * G4Material::GetMaterial ( std::size_t nComp,
G4double dens )
static

Definition at line 697 of file G4Material.cc.

698{
699 // search the material by its name
700 for (auto const & mat : *GetMaterialTable()) {
701 if (nComp == mat->GetNumberOfElements() && dens == mat->GetDensity()) {
702 return mat;
703 }
704 }
705 return nullptr;
706}

◆ GetMaterialPropertiesTable()

◆ GetMaterialTable()

G4MaterialTable * G4Material::GetMaterialTable ( )
static

Definition at line 644 of file G4Material.cc.

645{
646 struct Holder {
648 ~Holder() {
649 for(auto item : instance)
650 delete item;
651 }
652 };
653 static Holder _holder;
654 return &_holder.instance;
655}
std::vector< G4Material * > G4MaterialTable
G4TemplateRNGHelper< G4long > * G4TemplateRNGHelper< G4long >::instance

Referenced by G4VCrossSectionHandler::ActiveElements(), G4AdjointCSManager::BuildCrossSectionMatrices(), G4Cerenkov::BuildPhysicsTable(), G4CrossSectionDataStore::BuildPhysicsTable(), G4CrossSectionHP::BuildPhysicsTable(), G4GammaConversionToMuons::BuildPhysicsTable(), G4NeutronGeneralProcess::BuildPhysicsTable(), G4OpRayleigh::BuildPhysicsTable(), G4OpWLS2::BuildPhysicsTable(), G4OpWLS::BuildPhysicsTable(), G4ParticleHPThermalScatteringData::BuildPhysicsTable(), G4Scintillation::BuildPhysicsTable(), G4VLEPTSModel::BuildPhysicsTable(), G4VXTRenergyLoss::ComputeGasPhotoAbsCof(), G4PAIxSection::ComputeLowEnergyCof(), G4StrawTubeXTRadiator::ComputeMediumPhotoAbsCof(), G4VXTRenergyLoss::ComputePlatePhotoAbsCof(), G4DNACPA100ElasticModel::CrossSectionPerVolume(), G4DNACPA100ExcitationModel::CrossSectionPerVolume(), G4DNACPA100IonisationModel::CrossSectionPerVolume(), G4EmSaturation::DumpBirksCoefficients(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), G4HadXSHelper::FillPeaksStructure(), G4HadXSHelper::FindCrossSectionMax(), G4NistMaterialBuilder::FindMaterial(), G4PAIxSection::G4PAIxSection(), G4PAIxSection::G4PAIxSection(), G4PAIxSection::G4PAIxSection(), G4PAIxSection::G4PAIxSection(), G4SandiaTable::G4SandiaTable(), G4DNAMolecularMaterial::GetDensityTableFor(), G4VXTRenergyLoss::GetGasCompton(), GetMaterial(), GetMaterial(), GetMaterial(), GetNumberOfMaterials(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VXTRenergyLoss::GetPlateCompton(), G4TablesForExtrapolator::Initialisation(), G4ASTARStopping::Initialise(), G4ICRU73QOModel::Initialise(), G4ICRU90StoppingData::Initialise(), G4IonICRU73Data::Initialise(), G4LEPTSElasticModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4PEEffectFluoModel::Initialise(), G4PSTARStopping::Initialise(), G4EmSaturation::InitialiseG4Saturation(), G4DNAMolecularMaterial::Initialize(), G4DNAMolecularMaterial::InitializeDensity(), G4VDNAModel::IsMaterialDefine(), G4NeutronGeneralProcess::PreparePhysicsTable(), G4NistManager::PrintG4Material(), G4VDNAModel::RandomSelectShell(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4XrayReflection::SaveHenkeDataAsMaterialProperty(), G4NistMessenger::SetNewValue(), G4ProductionCutsTable::StoreMaterialInfo(), G4DNAModelInterface::StreamInfo(), G4GDMLRead::StripNames(), and ~G4Material().

◆ GetName()

const G4String & G4Material::GetName ( ) const
inline

Definition at line 171 of file G4Material.hh.

171{ return fName; }

Referenced by AddMaterial(), G4GMocrenFileSceneHandler::AddSolid(), G4GMocrenFileSceneHandler::AddSolid(), G4ErrorEnergyLoss::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4ParticleHPInelastic::ApplyYourself(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4GammaGeneralProcess::BuildPhysicsTable(), G4MicroElecSurface::BuildPhysicsTable(), G4NeutronGeneralProcess::BuildPhysicsTable(), G4ParticleHPThermalScatteringData::BuildPhysicsTable(), G4VLEPTSModel::BuildPhysicsTable(), G4PenelopeBremsstrahlungFS::BuildScaledXSTable(), G4PenelopeIonisationXSHandler::BuildXSTable(), G4ProductionCutsTable::CheckMaterialCutsCoupleInfo(), G4EmCalculator::ComputeCrossSectionPerVolume(), G4PenelopeIonisationModel::ComputeDEDXPerVolume(), G4VCrossSectionDataSet::ComputeIsoCrossSection(), G4EmCalculator::ComputeMeanFreePath(), G4tgbVolume::ConstructG4LogVol(), G4tgbVolume::ConstructG4PhysVol(), G4VRangeToEnergyConverter::Convert(), G4PenelopeIonisationCrossSection::CrossSection(), G4DNAModelInterface::CrossSectionPerVolume(), G4MicroElecElasticModel_new::CrossSectionPerVolume(), G4MicroElecInelasticModel_new::CrossSectionPerVolume(), G4MicroElecLOPhononModel::CrossSectionPerVolume(), G4PenelopeBremsstrahlungModel::CrossSectionPerVolume(), G4PenelopeComptonModel::CrossSectionPerVolume(), G4PenelopeIonisationModel::CrossSectionPerVolume(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4PenelopeOscillatorManager::Dump(), G4EmSaturation::DumpBirksCoefficients(), G4ProductionCutsTable::DumpCouples(), G4PenelopeRayleighModel::DumpFormFactorTable(), G4PenelopeRayleighModelMI::DumpFormFactorTable(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), G4HadronicProcess::DumpState(), G4EmCorrections::EffectiveChargeCorrection(), G4EmModelManager::FillDEDXVector(), G4EmModelManager::FillLambdaVector(), G4EmSaturation::FindG4BirksCoefficient(), G4ElectronIonPair::FindG4MeanEnergyPerIonPair(), G4IonisParamMat::FindMeanExcitationEnergy(), G4tgbMaterialMgr::FindOrBuildG4Material(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4ForwardXrayTR::G4ForwardXrayTR(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4PenelopeOscillatorManager::GetAtomsPerMolecule(), G4CrossSectionDataStore::GetCrossSection(), G4EmCalculator::GetCrossSectionPerVolume(), G4EmCalculator::GetCSDARange(), G4EmCalculator::GetDEDX(), G4PenelopeIonisationXSHandler::GetDensityCorrection(), G4PenelopeBremsstrahlungFS::GetEffectiveZSquared(), G4VCrossSectionDataSet::GetElementCrossSection(), G4ESTARStopping::GetIndex(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4EmCalculator::GetKinEnergy(), G4LatticeManager::GetLattice(), G4PenelopeOscillatorManager::GetMeanExcitationEnergy(), G4EmCalculator::GetMeanFreePath(), G4MicroElecCapture::GetMeanFreePath(), G4XrayReflection::GetMeanFreePath(), G4PenelopeOscillatorManager::GetNumberOfZAtomsPerMolecule(), G4PenelopeOscillatorManager::GetOscillatorCompton(), G4PenelopeOscillatorManager::GetOscillatorIonisation(), G4PenelopeOscillatorManager::GetOscillatorTableCompton(), G4PenelopeOscillatorManager::GetOscillatorTableIonisation(), G4PenelopeOscillatorManager::GetPlasmaEnergySquared(), G4EmCalculator::GetRangeFromRestricteDEDX(), G4PenelopeOscillatorManager::GetTotalA(), G4PenelopeOscillatorManager::GetTotalZ(), G4EmCorrections::HighOrderCorrections(), G4ASTARStopping::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4EmModelManager::Initialise(), G4ICRU90StoppingData::Initialise(), G4IonICRU73Data::Initialise(), G4LEPTSElasticModel::Initialise(), G4MicroElecCapture::Initialise(), G4MicroElecElasticModel_new::Initialise(), G4MicroElecInelasticModel_new::Initialise(), G4MicroElecSurface::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4PSTARStopping::Initialise(), G4LatticeManager::LoadLattice(), G4GDMLWriteMaterials::MaterialWrite(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4MicroElecCapture::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4NeutrinoElectronProcess::PostStepDoIt(), G4TauNeutrinoNucleusProcess::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4Decay::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4hImpactIonisation::PrintInfoDefinition(), GVFlashShowerParameterisation::PrintMaterial(), G4DNAMolecularMaterial::PrintNotAMolecularMaterial(), G4ExtendedMaterial::RegisterExtension(), G4TransportationLogger::ReportLoopingTrack(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4ExtendedMaterial::RetrieveExtension(), G4PenelopeBremsstrahlungAngular::SampleDirection(), G4VLEPTSModel::SampleEnergyLoss(), G4PenelopeBremsstrahlungFS::SampleGammaEnergy(), G4DNAPTBIonisationModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4PenelopeBremsstrahlungModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4ChannelingFastSimCrystalData::SetMaterialProperties(), G4EmCalculator::SetupMaterial(), G4EnergySplitter::SplitEnergyInVolumes(), G4ProductionCutsTable::StoreMaterialCutsCoupleInfo(), G4GDMLRead::StripNames(), G4ParallelWorldProcess::SwitchMaterial(), G4tgbMaterialMixtureByVolume::TransformToFractionsByWeight(), G4GDMLWriteStructure::TraverseVolumeTree(), and G4MSSteppingAction::UserSteppingAction().

◆ GetNuclearInterLength()

G4double G4Material::GetNuclearInterLength ( ) const
inline

Definition at line 209 of file G4Material.hh.

209{ return fNuclInterLen; }

Referenced by G4MSSteppingAction::UserSteppingAction().

◆ GetNumberOfElements()

std::size_t G4Material::GetNumberOfElements ( ) const
inline

Definition at line 180 of file G4Material.hh.

180{ return fNumberOfElements; }

Referenced by G4VCrossSectionHandler::ActiveElements(), AddMaterial(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4FissLib::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4ParticleHPCaptureURR::ApplyYourself(), G4ParticleHPElastic::ApplyYourself(), G4ParticleHPElasticURR::ApplyYourself(), G4ParticleHPFission::ApplyYourself(), G4ParticleHPFissionURR::ApplyYourself(), G4ParticleHPInelastic::ApplyYourself(), G4ParticleHPInelasticURR::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4CrossSectionHandler::BuildCrossSectionsForMaterials(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4CrossSectionDataStore::BuildPhysicsTable(), G4ParticleHPThermalScatteringData::BuildPhysicsTable(), G4PenelopeBremsstrahlungFS::BuildScaledXSTable(), G4Nucleus::ChooseParameters(), G4AdjointCSManager::ComputeAdjointCS(), G4CrossSectionDataStore::ComputeCrossSection(), G4eBremParametrizedModel::ComputeDEDXPerVolume(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4LivermoreIonisationModel::ComputeDEDXPerVolume(), G4MuBremsstrahlungModel::ComputeDEDXPerVolume(), G4MuPairProductionModel::ComputeDEDXPerVolume(), G4RiGeMuPairProductionModel::ComputeDEDXPerVolume(), G4PAIxSection::ComputeLowEnergyCof(), G4PAIySection::ComputeLowEnergyCof(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4DNADiracRMatrixExcitationModel::CrossSectionPerVolume(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4DNAQuinnPlasmonExcitationModel::CrossSectionPerVolume(), G4DNARelativisticIonisationModel::CrossSectionPerVolume(), G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4VEmModel::CrossSectionPerVolume(), G4tgbGeometryDumper::DumpMaterial(), G4DNAELSEPAElasticModel::G4DNAELSEPAElasticModel(), G4ShellVacancy::GenerateNumberOfIonisations(), G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(), G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(), G4IonICRU73Data::GetDEDX(), GVFlashShowerParameterisation::GetEffA(), GVFlashShowerParameterisation::GetEffZ(), G4HadronicProcessStore::GetElasticCrossSectionPerVolume(), G4HadronicProcessStore::GetFissionCrossSectionPerVolume(), G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(), G4hICRU49He::HasMaterial(), G4hICRU49p::HasMaterial(), G4hZiegler1985p::HasMaterial(), G4BoldyshevTripletModel::Initialise(), G4DNAELSEPAElasticModel::Initialise(), G4DNAQuinnPlasmonExcitationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4IonICRU73Data::Initialise(), G4JAEAElasticScatteringModel::Initialise(), G4JAEAPolarizedElasticScatteringModel::Initialise(), G4LivermoreNuclearGammaConversionModel::Initialise(), G4LivermorePolarizedComptonModel::Initialise(), G4LivermorePolarizedGammaConversionModel::Initialise(), G4LowEPComptonModel::Initialise(), G4LowEPPolarizedComptonModel::Initialise(), G4PenelopeGammaConversionModel::Initialise(), G4PenelopePhotoElectricModel::Initialise(), G4PenelopeRayleighModel::Initialise(), G4PenelopeRayleighModelMI::Initialise(), G4WentzelVIModel::Initialise(), G4VEmModel::InitialiseForMaterial(), G4QAOLowEnergyLoss::IsInCharge(), G4QAOLowEnergyLoss::IsInCharge(), G4GDMLWriteMaterials::MaterialWrite(), G4eCoulombScatteringModel::MinPrimaryEnergy(), G4hCoulombScatteringModel::MinPrimaryEnergy(), G4MicroElecCapture::PostStepDoIt(), G4VLEPTSModel::ReadParam(), G4EmUtility::SampleRandomElement(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4CrossSectionDataStore::SampleZandA(), G4PixeCrossSectionHandler::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomAtom(), G4VEmModel::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomElement(), G4ElementSelector::SelectZandA(), G4ChannelingFastSimCrystalData::SetMaterialProperties(), G4hICRU49He::StoppingPower(), G4hICRU49p::StoppingPower(), G4hZiegler1985p::StoppingPower(), G4PixeCrossSectionHandler::ValueForMaterial(), and G4VCrossSectionHandler::ValueForMaterial().

◆ GetNumberOfMaterials()

◆ GetPressure()

◆ GetRadlen()

◆ GetSandiaTable()

◆ GetState()

◆ GetTemperature()

G4double G4Material::GetTemperature ( ) const
inline

Definition at line 176 of file G4Material.hh.

176{ return fTemp; }

Referenced by G4FissionLibrary::ApplyYourself(), G4FissLib::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4LENDModel::ApplyYourself(), G4NeutronFissionVI::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4NeutronHPCaptureFS::ApplyYourself(), G4NeutronRadCaptureHP::ApplyYourself(), G4ParticleHPChannel::ApplyYourself(), G4ParticleHPChannelList::ApplyYourself(), G4ParticleHPChannelList::ApplyYourself(), G4ParticleHPElastic::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4ParticleHPFission::ApplyYourself(), G4ParticleHPFissionFS::ApplyYourself(), G4ParticleHPInelastic::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4tgbMaterialMixtureByWeight::BuildG4Material(), G4NistManager::BuildMaterialWithNewDensity(), G4ParticleHPInelasticCompFS::CompositeApply(), G4CrossSectionHP::ComputeIsoCrossSection(), G4DNABrownianTransportation::ComputeStep(), G4NistMaterialBuilder::ConstructNewGasMaterial(), G4tgbGeometryDumper::DumpMaterial(), G4ExtendedMaterial::G4ExtendedMaterial(), G4Nucleus::G4Nucleus(), G4ParticleHPThermalScatteringData::GetCoherentCrossSection(), G4ParticleHPThermalScatteringData::GetCrossSection(), G4ParticleHPThermalScatteringData::GetIncoherentCrossSection(), G4ParticleHPThermalScatteringData::GetInelasticCrossSection(), G4CrossSectionHP::GetIsoCrossSection(), G4LENDCrossSection::GetIsoCrossSection(), G4NeutronHPCaptureData::GetIsoCrossSection(), G4ParticleHPElasticData::GetIsoCrossSection(), G4ParticleHPFissionData::GetIsoCrossSection(), G4ParticleHPInelasticData::GetIsoCrossSection(), G4ParticleHPProbabilityTablesStore::GetIsoCrossSectionPT(), and G4GDMLWriteMaterials::MaterialWrite().

◆ GetTotNbOfAtomsPerVolume()

◆ GetTotNbOfElectPerVolume()

G4double G4Material::GetTotNbOfElectPerVolume ( ) const
inline

◆ GetVecNbOfAtomsPerVolume()

◆ GetZ()

◆ IsExtended()

G4bool G4Material::IsExtended ( ) const
virtual

Reimplemented in G4ExtendedMaterial.

Definition at line 799 of file G4Material.cc.

799{ return false; }

Referenced by G4PenelopeRayleighModelMI::CrossSectionPerVolume(), and operator<<.

◆ operator!=()

G4bool G4Material::operator!= ( const G4Material & ) const
delete

◆ operator=()

const G4Material & G4Material::operator= ( const G4Material & )
delete

◆ operator==()

G4bool G4Material::operator== ( const G4Material & ) const
delete

◆ SetChemicalFormula()

void G4Material::SetChemicalFormula ( const G4String & chF)

Definition at line 614 of file G4Material.cc.

615{
616 if (! IsLocked()) {
617 fChemicalFormula = chF;
618 }
619}

◆ SetFreeElectronDensity()

void G4Material::SetFreeElectronDensity ( G4double val)

Definition at line 623 of file G4Material.cc.

624{
625 if (val >= 0. && ! IsLocked()) {
626 fFreeElecDensity = val;
627 }
628}

◆ SetMaterialPropertiesTable()

void G4Material::SetMaterialPropertiesTable ( G4MaterialPropertiesTable * anMPT)

Definition at line 803 of file G4Material.cc.

804{
805 if (fMaterialPropertiesTable != anMPT && ! IsLocked()) {
806 delete fMaterialPropertiesTable;
807 fMaterialPropertiesTable = anMPT;
808 }
809}

Referenced by G4GDMLReadMaterials::PropertyRead().

◆ SetName()

void G4Material::SetName ( const G4String & name)
inline

Definition at line 260 of file G4Material.hh.

260{ fName = name; }

Referenced by G4GDMLRead::StripNames().

Friends And Related Symbol Documentation

◆ operator<< [1/3]

std::ostream & operator<< ( std::ostream & flux,
const G4Material & material )
friend

Definition at line 777 of file G4Material.cc.

778{
779 flux << &material;
780 return flux;
781}

◆ operator<< [2/3]

std::ostream & operator<< ( std::ostream & flux,
const G4Material * material )
friend

Definition at line 736 of file G4Material.cc.

737{
738 std::ios::fmtflags mode = flux.flags();
739 flux.setf(std::ios::fixed, std::ios::floatfield);
740 G4long prec = flux.precision(3);
741
742 flux << " Material: " << std::setw(8) << material->fName << " " << material->fChemicalFormula
743 << " "
744 << " density: " << std::setw(6) << std::setprecision(3)
745 << G4BestUnit(material->fDensity, "Volumic Mass") << " RadL: " << std::setw(7)
746 << std::setprecision(3) << G4BestUnit(material->fRadlen, "Length")
747 << " Nucl.Int.Length: " << std::setw(7) << std::setprecision(3)
748 << G4BestUnit(material->fNuclInterLen, "Length") << "\n"
749 << std::setw(30) << " Imean: " << std::setw(7) << std::setprecision(3)
750 << G4BestUnit(material->GetIonisation()->GetMeanExcitationEnergy(), "Energy")
751 << " temperature: " << std::setw(6) << std::setprecision(2)
752 << (material->fTemp) / CLHEP::kelvin << " K"
753 << " pressure: " << std::setw(6) << std::setprecision(2)
754 << (material->fPressure) / CLHEP::atmosphere << " atm"
755 << "\n";
756
757 for (G4int i = 0; i < material->fNumberOfElements; i++) {
758 flux << "\n ---> " << (*(material->theElementVector))[i]
759 << "\n ElmMassFraction: " << std::setw(6) << std::setprecision(2)
760 << (material->fMassFractionVector[i]) / perCent << " %"
761 << " ElmAbundance " << std::setw(6) << std::setprecision(2)
762 << 100 * (material->fVecNbOfAtomsPerVolume[i]) / (material->fTotNbOfAtomsPerVolume)
763 << " % \n";
764 }
765 flux.precision(prec);
766 flux.setf(mode, std::ios::floatfield);
767
768 if (material->IsExtended()) {
769 static_cast<const G4ExtendedMaterial*>(material)->Print(flux);
770 }
771
772 return flux;
773}
#define G4BestUnit(a, b)
long G4long
Definition G4Types.hh:87
G4double GetMeanExcitationEnergy() const
virtual G4bool IsExtended() const
G4IonisParamMat * GetIonisation() const

◆ operator<< [3/3]

std::ostream & operator<< ( std::ostream & flux,
const G4MaterialTable & MaterialTable )
friend

Definition at line 785 of file G4Material.cc.

786{
787 // Dump info for all known materials
788 flux << "\n***** Table : Nb of materials = " << MaterialTable.size() << " *****\n" << G4endl;
789
790 for (auto i : MaterialTable) {
791 flux << i << G4endl << G4endl;
792 }
793
794 return flux;
795}

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