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

#include <G4Material.hh>

Public Member Functions

 G4Material (const G4String &name, G4double z, G4double a, G4double density, G4State state=kStateUndefined, G4double temp=CLHEP::STP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
 G4Material (const G4String &name, G4double density, G4int nComponents, G4State state=kStateUndefined, G4double temp=CLHEP::STP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
 G4Material (const G4String &name, G4double density, const G4Material *baseMaterial, G4State state=kStateUndefined, G4double temp=CLHEP::STP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
void AddElement (G4Element *element, G4int nAtoms)
 
void AddElement (G4Element *element, G4double fraction)
 
void AddMaterial (G4Material *material, G4double fraction)
 
virtual ~G4Material ()
 
void SetChemicalFormula (const G4String &chF)
 
const G4StringGetName () const
 
const G4StringGetChemicalFormula () const
 
G4double GetDensity () const
 
G4State GetState () const
 
G4double GetTemperature () const
 
G4double GetPressure () const
 
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
 
std::map< G4Material *, G4doubleGetMatComponents () const
 
G4double GetMassOfMolecule () const
 
G4double GetZ () const
 
G4double GetA () const
 
void SetMaterialPropertiesTable (G4MaterialPropertiesTable *anMPT)
 
G4MaterialPropertiesTableGetMaterialPropertiesTable () const
 
size_t GetIndex () const
 
G4int operator== (const G4Material &) const
 
G4int operator!= (const G4Material &) const
 
 G4Material (__void__ &)
 
void SetName (const G4String &name)
 

Static Public Member Functions

static const G4MaterialTableGetMaterialTable ()
 
static size_t GetNumberOfMaterials ()
 
static G4MaterialGetMaterial (const G4String &name, G4bool warning=true)
 

Friends

std::ostream & operator<< (std::ostream &, G4Material *)
 
std::ostream & operator<< (std::ostream &, G4Material &)
 
std::ostream & operator<< (std::ostream &, G4MaterialTable)
 

Detailed Description

Definition at line 118 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 = CLHEP::STP_Temperature,
G4double  pressure = CLHEP::STP_Pressure 
)

Definition at line 88 of file G4Material.cc.

91 : fName(name)
92{
93 InitializePointers();
94
95 if (density < universe_mean_density)
96 {
97 G4cout << "--- Warning from G4Material::G4Material()"
98 << " define a material with density=0 is not allowed. \n"
99 << " The material " << name << " will be constructed with the"
100 << " default minimal density: " << universe_mean_density/(g/cm3)
101 << "g/cm3" << G4endl;
102 density = universe_mean_density;
103 }
104
105 fDensity = density;
106 fState = state;
107 fTemp = temp;
108 fPressure = pressure;
109
110 // Initialize theElementVector allocating one
111 // element corresponding to this material
112 maxNbComponents = fNumberOfComponents = fNumberOfElements = 1;
113 fArrayLength = maxNbComponents;
114 fImplicitElement = true;
115 theElementVector = new G4ElementVector();
116 theElementVector->push_back( new G4Element(name, " ", z, a));
117 fMassFractionVector = new G4double[1];
118 fMassFractionVector[0] = 1. ;
119 fMassOfMolecule = a/Avogadro;
120
121 (*theElementVector)[0] -> increaseCountUse();
122
123 if (fState == kStateUndefined)
124 {
125 if (fDensity > kGasThreshold) { fState = kStateSolid; }
126 else { fState = kStateGas; }
127 }
128
129 ComputeDerivedQuantities();
130}
std::vector< G4Element * > G4ElementVector
@ kStateSolid
Definition: G4Material.hh:114
@ kStateGas
Definition: G4Material.hh:114
@ kStateUndefined
Definition: G4Material.hh:114
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ G4Material() [2/4]

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

Definition at line 137 of file G4Material.cc.

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

◆ G4Material() [3/4]

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

Definition at line 176 of file G4Material.cc.

179 : fName(name)
180{
181 InitializePointers();
182
183 if (density < universe_mean_density)
184 {
185 G4cout << "--- Warning from G4Material::G4Material()"
186 << " define a material with density=0 is not allowed. \n"
187 << " The material " << name << " will be constructed with the"
188 << " default minimal density: " << universe_mean_density/(g/cm3)
189 << "g/cm3" << G4endl;
190 density = universe_mean_density;
191 }
192
193 fDensity = density;
194 fState = state;
195 fTemp = temp;
196 fPressure = pressure;
197
198 fBaseMaterial = bmat;
199 fChemicalFormula = fBaseMaterial->GetChemicalFormula();
200 fMassOfMolecule = fBaseMaterial->GetMassOfMolecule();
201
202 fNumberOfElements = fBaseMaterial->GetNumberOfElements();
203 maxNbComponents = fNumberOfElements;
204 fNumberOfComponents = fNumberOfElements;
205
206 fMaterialPropertiesTable = fBaseMaterial->GetMaterialPropertiesTable();
207
208 CopyPointersOfBaseMaterial();
209}
const G4String & GetChemicalFormula() const
Definition: G4Material.hh:178
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
G4double GetMassOfMolecule() const
Definition: G4Material.hh:241

◆ ~G4Material()

G4Material::~G4Material ( )
virtual

Definition at line 227 of file G4Material.cc.

228{
229 // G4cout << "### Destruction of material " << fName << " started" <<G4endl;
230 if(!fBaseMaterial) {
231 if (theElementVector) { delete theElementVector; }
232 if (fMassFractionVector) { delete [] fMassFractionVector; }
233 if (fAtomsVector) { delete [] fAtomsVector; }
234 if (fSandiaTable) { delete fSandiaTable; }
235 }
236 if (fIonisation) { delete fIonisation; }
237 if (VecNbOfAtomsPerVolume) { delete [] VecNbOfAtomsPerVolume; }
238
239 // Remove this material from theMaterialTable.
240 //
241 theMaterialTable[fIndexInTable] = 0;
242}

◆ G4Material() [4/4]

G4Material::G4Material ( __void__ &  )

Definition at line 216 of file G4Material.cc.

217 : fNumberOfComponents(0), fNumberOfElements(0), theElementVector(0),
218 fImplicitElement(false), fMassFractionVector(0), fAtomsVector(0),
219 fMaterialPropertiesTable(0), fIndexInTable(0),
220 VecNbOfAtomsPerVolume(0), fIonisation(0), fSandiaTable(0)
221{
222 InitializePointers();
223}

Member Function Documentation

◆ AddElement() [1/2]

void G4Material::AddElement ( G4Element element,
G4double  fraction 
)

Definition at line 384 of file G4Material.cc.

385{
386 if(fraction < 0.0 || fraction > 1.0) {
387 G4cout << "G4Material::AddElement ERROR for " << fName << " and "
388 << element->GetName() << " mass fraction= " << fraction
389 << " is wrong " << G4endl;
390 G4Exception ("G4Material::AddElement()", "mat032", FatalException,
391 "Attempt to add element with wrong mass fraction");
392 }
393 // initialization
394 if (fNumberOfComponents == 0) {
395 fMassFractionVector = new G4double[fArrayLength];
396 fAtomsVector = new G4int [fArrayLength];
397 }
398 // filling ...
399 if (G4int(fNumberOfComponents) < maxNbComponents) {
400 size_t el = 0;
401 while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) { ++el; }
402 if (el<fNumberOfElements) fMassFractionVector[el] += fraction;
403 else {
404 theElementVector->push_back(element);
405 fMassFractionVector[el] = fraction;
406 ++fNumberOfElements;
407 element->increaseCountUse();
408 }
409 ++fNumberOfComponents;
410 } else {
411 G4cout << "G4Material::AddElement ERROR for " << fName << " nElement= "
412 << fNumberOfElements << G4endl;
413 G4Exception ("G4Material::AddElement()", "mat033", FatalException,
414 "Attempt to add more than the declared number of elements.");
415 }
416
417 // filled.
418 if (G4int(fNumberOfComponents) == maxNbComponents) {
419
420 size_t i=0;
421 G4double Zmol(0.), Amol(0.);
422 // check sum of weights -- OK?
423 G4double wtSum(0.0);
424 for (i=0; i<fNumberOfElements; ++i) {
425 wtSum += fMassFractionVector[i];
426 Zmol += fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
427 Amol += fMassFractionVector[i]*(*theElementVector)[i]->GetA();
428 }
429 if (std::fabs(1.-wtSum) > perThousand) {
430 G4cerr << "WARNING !! for " << fName << " sum of fractional masses "
431 << wtSum << " is not 1 - results may be wrong"
432 << G4endl;
433 }
434 for (i=0; i<fNumberOfElements; ++i) {
435 fAtomsVector[i] =
436 G4int(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA()+0.5);
437 }
438
439 ComputeDerivedQuantities();
440 }
441}
@ FatalException
int G4int
Definition: G4Types.hh:66
G4DLLIMPORT std::ostream G4cerr
const G4String & GetName() const
Definition: G4Element.hh:127
void increaseCountUse()
Definition: G4Element.hh:192
G4double GetA() const
Definition: G4Material.cc:617
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ AddElement() [2/2]

void G4Material::AddElement ( G4Element element,
G4int  nAtoms 
)

Definition at line 341 of file G4Material.cc.

342{
343 // initialization
344 if ( fNumberOfElements == 0 ) {
345 fAtomsVector = new G4int [fArrayLength];
346 fMassFractionVector = new G4double[fArrayLength];
347 }
348
349 // filling ...
350 if ( G4int(fNumberOfElements) < maxNbComponents ) {
351 theElementVector->push_back(element);
352 fAtomsVector[fNumberOfElements] = nAtoms;
353 fNumberOfComponents = ++fNumberOfElements;
354 element->increaseCountUse();
355 } else {
356 G4cout << "G4Material::AddElement ERROR for " << fName << " nElement= "
357 << fNumberOfElements << G4endl;
358 G4Exception ("G4Material::AddElement()", "mat031", FatalException,
359 "Attempt to add more than the declared number of elements.");
360 }
361 // filled.
362 if ( G4int(fNumberOfElements) == maxNbComponents ) {
363 // compute proportion by mass
364 size_t i=0;
365 G4double Amol = 0.;
366 for (i=0; i<fNumberOfElements; ++i) {
367 G4double w = fAtomsVector[i]*(*theElementVector)[i]->GetA();
368 Amol += w;
369 fMassFractionVector[i] = w;
370 }
371 for (i=0; i<fNumberOfElements; ++i) {
372 fMassFractionVector[i] /= Amol;
373 }
374
375 fMassOfMolecule = Amol/Avogadro;
376 ComputeDerivedQuantities();
377 }
378}

Referenced by G4tgbMaterialMixtureByNoAtoms::BuildG4Material(), G4tgbMaterialMixtureByWeight::BuildG4Material(), G4gsmate(), G4gsmixt(), and G4GDMLReadMaterials::MixtureRead().

◆ AddMaterial()

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

store massFraction of material component

Definition at line 447 of file G4Material.cc.

448{
449 if(fraction < 0.0 || fraction > 1.0) {
450 G4cout << "G4Material::AddMaterial ERROR for " << fName << " and "
451 << material->GetName() << " mass fraction= " << fraction
452 << " is wrong " << G4endl;
453 G4Exception ("G4Material::AddMaterial()", "mat034", FatalException,
454 "Attempt to add material with wrong mass fraction");
455 }
456 // initialization
457 if (fNumberOfComponents == 0) {
458 fMassFractionVector = new G4double[fArrayLength];
459 fAtomsVector = new G4int [fArrayLength];
460 }
461
462 size_t nelm = material->GetNumberOfElements();
463
464 // arrays should be extended
465 if(nelm > 1) {
466 G4int nold = fArrayLength;
467 fArrayLength += nelm - 1;
468 G4double* v1 = new G4double[fArrayLength];
469 G4int* i1 = new G4int[fArrayLength];
470 for(G4int i=0; i<nold; ++i) {
471 v1[i] = fMassFractionVector[i];
472 i1[i] = fAtomsVector[i];
473 }
474 delete [] fAtomsVector;
475 delete [] fMassFractionVector;
476 fMassFractionVector = v1;
477 fAtomsVector = i1;
478 }
479
480 // filling ...
481 if (G4int(fNumberOfComponents) < maxNbComponents) {
482 for (size_t elm=0; elm<nelm; ++elm)
483 {
484 G4Element* element = (*(material->GetElementVector()))[elm];
485 size_t el = 0;
486 while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) el++;
487 if (el < fNumberOfElements) fMassFractionVector[el] += fraction
488 *(material->GetFractionVector())[elm];
489 else {
490 theElementVector->push_back(element);
491 fMassFractionVector[el] = fraction
492 *(material->GetFractionVector())[elm];
493 ++fNumberOfElements;
494 element->increaseCountUse();
495 }
496 }
497 ++fNumberOfComponents;
498 ///store massFraction of material component
499 fMatComponents[material] = fraction;
500
501 } else {
502 G4cout << "G4Material::AddMaterial ERROR for " << fName << " nElement= "
503 << fNumberOfElements << G4endl;
504 G4Exception ("G4Material::AddMaterial()", "mat035", FatalException,
505 "Attempt to add more than the declared number of components.");
506 }
507
508 // filled.
509 if (G4int(fNumberOfComponents) == maxNbComponents) {
510 size_t i=0;
511 G4double Zmol(0.), Amol(0.);
512 // check sum of weights -- OK?
513 G4double wtSum(0.0);
514 for (i=0; i<fNumberOfElements; ++i) {
515 wtSum += fMassFractionVector[i];
516 Zmol += fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
517 Amol += fMassFractionVector[i]*(*theElementVector)[i]->GetA();
518 }
519 if (std::fabs(1.-wtSum) > perThousand) {
520 G4cout << "G4Material::AddMaterial WARNING !! for " << fName
521 << " sum of fractional masses "
522 << wtSum << " is not 1 - results may be wrong"
523 << G4endl;
524 }
525 for (i=0;i<fNumberOfElements;i++) {
526 fAtomsVector[i] =
527 G4int(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA()+0.5);
528 }
529
530 ComputeDerivedQuantities();
531 }
532}
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
const G4double * GetFractionVector() const
Definition: G4Material.hh:193
const G4String & GetName() const
Definition: G4Material.hh:177

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

◆ GetA()

G4double G4Material::GetA ( ) const

Definition at line 617 of file G4Material.cc.

618{
619 if (fNumberOfElements > 1) {
620 G4cout << "G4Material ERROR in GetA. The material: " << fName << " is a mixture."
621 << G4endl;
622 G4Exception ("G4Material::GetA()", "mat037", FatalException,
623 "the Atomic mass is not well defined." );
624 }
625 return (*theElementVector)[0]->GetA();
626}

Referenced by AddElement(), AddMaterial(), G4tgbGeometryDumper::DumpMaterial(), GVFlashShowerParameterisation::GetEffA(), and G4GDMLWriteMaterials::MaterialWrite().

◆ GetAtomicNumDensityVector()

◆ GetAtomsVector()

const G4int * G4Material::GetAtomsVector ( ) const
inline

Definition at line 197 of file G4Material.hh.

197{return fAtomsVector;}

◆ GetBaseMaterial()

◆ GetChemicalFormula()

const G4String & G4Material::GetChemicalFormula ( ) const
inline

Definition at line 178 of file G4Material.hh.

178{return fChemicalFormula;}

Referenced by G4Material(), G4hICRU49He::HasMaterial(), and G4hICRU49p::HasMaterial().

◆ GetDensity()

◆ GetElectronDensity()

G4double G4Material::GetElectronDensity ( ) const
inline

Definition at line 216 of file G4Material.hh.

216{return TotNbOfElectPerVolume;}

Referenced by G4AdjointComptonModel::AdjointCrossSection(), G4AdjointhIonisationModel::AdjointCrossSection(), G4ForwardXrayTR::BuildXrayTRtables(), G4MuBetheBlochModel::ComputeDEDXPerVolume(), G4BetheBlochModel::ComputeDEDXPerVolume(), G4BraggIonModel::ComputeDEDXPerVolume(), G4BraggModel::ComputeDEDXPerVolume(), G4eBremsstrahlungModel::ComputeDEDXPerVolume(), G4MollerBhabhaModel::ComputeDEDXPerVolume(), G4ICRU73QOModel::ComputeDEDXPerVolume(), G4EmCorrections::ComputeIonCorrections(), G4eeToHadronsModel::CrossSectionPerVolume(), G4eeToHadronsMultiModel::CrossSectionPerVolume(), G4MuBetheBlochModel::CrossSectionPerVolume(), G4BetheBlochModel::CrossSectionPerVolume(), G4BraggIonModel::CrossSectionPerVolume(), G4BraggModel::CrossSectionPerVolume(), G4eBremsstrahlungModel::CrossSectionPerVolume(), G4eeToTwoGammaModel::CrossSectionPerVolume(), G4ICRU73QOModel::CrossSectionPerVolume(), G4MollerBhabhaModel::CrossSectionPerVolume(), G4BohrFluctuations::Dispersion(), G4UniversalFluctuation::Dispersion(), G4mplIonisationModel::Dispersion(), G4mplIonisationWithDeltaModel::Dispersion(), G4IonFluctuations::Dispersion(), G4InitXscPAI::G4InitXscPAI(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4VEnergyLoss::GetLossWithFluct(), G4VeLowEnergyLoss::GetLossWithFluct(), G4EmCorrections::HighOrderCorrections(), G4PAIySection::Initialize(), G4EmCorrections::IonBarkasCorrection(), G4AdjointComptonModel::RapidSampleSecondaries(), G4AdjointhIonisationModel::RapidSampleSecondaries(), G4UniversalFluctuation::SampleFluctuations(), G4eBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SetupForMaterial(), and G4eBremsstrahlungRelModel::SetupForMaterial().

◆ GetElement()

◆ GetElementVector()

const G4ElementVector * G4Material::GetElementVector ( ) const
inline

Definition at line 189 of file G4Material.hh.

189{return theElementVector;}

Referenced by G4VCrossSectionHandler::ActiveElements(), AddMaterial(), G4AdjointPhotoElectricModel::AdjointCrossSection(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4RToEConvForGamma::BuildAbsorptionLengthVector(), G4AugerData::BuildAugerTransitionTable(), G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4CrossSectionHandler::BuildCrossSectionsForMaterials(), G4VRangeToEnergyConverter::BuildRangeVector(), G4Nucleus::ChooseParameters(), G4LivermoreBremsstrahlungModel::ComputeDEDXPerVolume(), G4LivermoreIonisationModel::ComputeDEDXPerVolume(), G4MuBremsstrahlungModel::ComputeDEDXPerVolume(), G4MuPairProductionModel::ComputeDEDXPerVolume(), G4eBremParametrizedModel::ComputeDEDXPerVolume(), G4eBremsstrahlungModel::ComputeDEDXPerVolume(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4eBremsstrahlungModel::CrossSectionPerVolume(), G4VEmModel::CrossSectionPerVolume(), G4AnnihiToMuPair::CrossSectionPerVolume(), G4AdjointBremsstrahlungModel::DiffCrossSectionPerVolumePrimToSecondApproximated2(), G4tgbGeometryDumper::DumpMaterial(), G4EmElementSelector::G4EmElementSelector(), G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(), G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(), G4CrossSectionDataStore::GetCrossSection(), G4HadronicProcessStore::GetElasticCrossSectionPerVolume(), G4StopElementSelector::GetElement(), G4HadronicProcessStore::GetFissionCrossSectionPerVolume(), G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(), G4QAtomicElectronScattering::GetMeanFreePath(), G4QCoherentChargeExchange::GetMeanFreePath(), G4QDiffraction::GetMeanFreePath(), G4QElastic::GetMeanFreePath(), G4QInelastic::GetMeanFreePath(), G4QIonIonElastic::GetMeanFreePath(), G4QLowEnergy::GetMeanFreePath(), G4QNGamma::GetMeanFreePath(), G4LivermoreGammaConversionModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), G4QAtomicElectronScattering::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4GoudsmitSaundersonMscModel::SampleScattering(), G4WentzelVIModel::SampleScattering(), G4WentzelVIRelModel::SampleScattering(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4CrossSectionDataStore::SampleZandA(), G4VEmModel::SelectRandomAtom(), G4PixeCrossSectionHandler::SelectRandomAtom(), G4eBremsstrahlungModel::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomElement(), G4ElementSelector::SelectZandA(), G4hQAOModel::StoppingPower(), G4VCrossSectionHandler::ValueForMaterial(), and G4PixeCrossSectionHandler::ValueForMaterial().

◆ GetFractionVector()

◆ GetIndex()

size_t G4Material::GetIndex ( ) const
inline

Definition at line 261 of file G4Material.hh.

261{return fIndexInTable;}

Referenced by G4AdjointCSManager::ComputeAdjointCS(), G4VRangeToEnergyConverter::Convert(), G4DNABornExcitationModel::CrossSectionPerVolume(), G4DNABornIonisationModel::CrossSectionPerVolume(), G4DNAChampionElasticModel::CrossSectionPerVolume(), G4DNADingfelderChargeDecreaseModel::CrossSectionPerVolume(), G4DNADingfelderChargeIncreaseModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4DNAMeltonAttachmentModel::CrossSectionPerVolume(), G4DNAMillerGreenExcitationModel::CrossSectionPerVolume(), G4DNARuddIonisationExtendedModel::CrossSectionPerVolume(), G4DNARuddIonisationModel::CrossSectionPerVolume(), G4DNASancheExcitationModel::CrossSectionPerVolume(), G4DNASancheSolvatationModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNATransformElectronModel::CrossSectionPerVolume(), G4VEmAdjointModel::DefineCurrentMaterial(), G4DNABrownianTransportation::Diffusion(), G4InitXscPAI::G4InitXscPAI(), G4PAIxSection::G4PAIxSection(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4EnergyLossTables::GetDEDX(), G4EnergyLossTables::GetDeltaLabTime(), G4EnergyLossTables::GetDeltaProperTime(), G4EnergyLossTables::GetLabTime(), G4OpRayleigh::GetMeanFreePath(), G4EnergyLossTables::GetPreciseDEDX(), G4EnergyLossTables::GetPreciseEnergyFromRange(), G4EnergyLossTables::GetPreciseRangeFromEnergy(), G4EnergyLossTables::GetProperTime(), G4EnergyLossTables::GetRange(), G4Scintillation::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4DNAMolecularMaterial::RecordMolecularMaterial(), and G4PixeCrossSectionHandler::SelectRandomAtom().

◆ GetIonisation()

◆ GetMassOfMolecule()

G4double G4Material::GetMassOfMolecule ( ) const
inline

◆ GetMatComponents()

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

Definition at line 236 of file G4Material.hh.

237 {return fMatComponents;}

Referenced by G4DNAMolecularMaterial::SearchMolecularMaterial().

◆ GetMaterial()

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

◆ GetMaterialPropertiesTable()

◆ GetMaterialTable()

◆ GetName()

const G4String & G4Material::GetName ( ) const
inline

Definition at line 177 of file G4Material.hh.

177{return fName;}

Referenced by AddMaterial(), G4GMocrenFileSceneHandler::AddSolid(), G4ErrorEnergyLoss::AlongStepDoIt(), G4IVContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4LEAlphaInelastic::ApplyYourself(), G4LEAntiKaonZeroInelastic::ApplyYourself(), G4LEAntiLambdaInelastic::ApplyYourself(), G4LEAntiNeutronInelastic::ApplyYourself(), G4LEAntiOmegaMinusInelastic::ApplyYourself(), G4LEAntiProtonInelastic::ApplyYourself(), G4LEAntiSigmaMinusInelastic::ApplyYourself(), G4LEAntiSigmaPlusInelastic::ApplyYourself(), G4LEAntiXiMinusInelastic::ApplyYourself(), G4LEAntiXiZeroInelastic::ApplyYourself(), G4LEDeuteronInelastic::ApplyYourself(), G4LEKaonMinusInelastic::ApplyYourself(), G4LEKaonPlusInelastic::ApplyYourself(), G4LEKaonZeroInelastic::ApplyYourself(), G4LELambdaInelastic::ApplyYourself(), G4LENeutronInelastic::ApplyYourself(), G4LEOmegaMinusInelastic::ApplyYourself(), G4LEPionMinusInelastic::ApplyYourself(), G4LEPionPlusInelastic::ApplyYourself(), G4LEProtonInelastic::ApplyYourself(), G4LESigmaMinusInelastic::ApplyYourself(), G4LESigmaPlusInelastic::ApplyYourself(), G4LETritonInelastic::ApplyYourself(), G4LEXiMinusInelastic::ApplyYourself(), G4LEXiZeroInelastic::ApplyYourself(), G4RPGAntiKZeroInelastic::ApplyYourself(), G4RPGAntiLambdaInelastic::ApplyYourself(), G4RPGAntiNeutronInelastic::ApplyYourself(), G4RPGAntiOmegaMinusInelastic::ApplyYourself(), G4RPGAntiProtonInelastic::ApplyYourself(), G4RPGAntiSigmaMinusInelastic::ApplyYourself(), G4RPGAntiSigmaPlusInelastic::ApplyYourself(), G4RPGAntiXiMinusInelastic::ApplyYourself(), G4RPGAntiXiZeroInelastic::ApplyYourself(), G4RPGKMinusInelastic::ApplyYourself(), G4RPGKPlusInelastic::ApplyYourself(), G4RPGKZeroInelastic::ApplyYourself(), G4RPGLambdaInelastic::ApplyYourself(), G4RPGOmegaMinusInelastic::ApplyYourself(), G4RPGSigmaMinusInelastic::ApplyYourself(), G4RPGSigmaPlusInelastic::ApplyYourself(), G4RPGXiMinusInelastic::ApplyYourself(), G4RPGXiZeroInelastic::ApplyYourself(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4AntiNeutronAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4AntiProtonAnnihilationAtRest::AtRestGetPhysicalInteractionLength(), G4KaonMinusAbsorption::AtRestGetPhysicalInteractionLength(), G4NeutronCaptureAtRest::AtRestGetPhysicalInteractionLength(), G4PionMinusAbsorptionAtRest::AtRestGetPhysicalInteractionLength(), G4IVRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4PionMinusNuclearAtRestChips::AtRestGetPhysicalInteractionLength(), G4ProtonAntiProtonAtRestChips::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4EmCorrections::BarkasCorrection(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4NeutronHPThermalScatteringData::BuildPhysicsTable(), G4ProductionCutsTable::CheckMaterialCutsCoupleInfo(), G4EmCalculator::ComputeCrossSectionPerVolume(), G4EmCalculator::ComputeDEDX(), G4EmCalculator::ComputeElectronicDEDX(), G4EmCalculator::ComputeMeanFreePath(), G4tgbVolume::ConstructG4LogVol(), G4tgbVolume::ConstructG4PhysVol(), G4VRangeToEnergyConverter::Convert(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G4PenelopeIonisationCrossSection::CrossSection(), G4PenelopeComptonModel::CrossSectionPerVolume(), G4PenelopeBremsstrahlungModel::CrossSectionPerVolume(), G4PenelopeIonisationModel::CrossSectionPerVolume(), G4PenelopeOscillatorManager::Dump(), G4EmElementSelector::Dump(), G4ProductionCutsTable::DumpCouples(), G4PenelopeRayleighModel::DumpFormFactorTable(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), G4HadronicProcess::DumpState(), G4EmModelManager::FillDEDXVector(), G4EmModelManager::FillLambdaVector(), G4EmCalculator::FindCouple(), G4EmSaturation::FindG4BirksCoefficient(), G4ElectronIonPair::FindG4MeanEnergyPerIonPair(), G4tgbMaterialMgr::FindOrBuildG4Material(), G4ForwardXrayTR::G4ForwardXrayTR(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4PenelopeOscillatorManager::GetAtomsPerMolecule(), G4CrossSectionDataStore::GetCrossSection(), G4EmCalculator::GetCrossSectionPerVolume(), G4PenelopeIonisationXSHandler::GetCrossSectionTableForCouple(), G4EmCalculator::GetCSDARange(), G4EmCalculator::GetDEDX(), G4PenelopeIonisationXSHandler::GetDensityCorrection(), G4PenelopeBremsstrahlungFS::GetEffectiveZSquared(), G4VCrossSectionDataSet::GetElementCrossSection(), G4EnergyRangeManager::GetHadronicInteraction(), G4ASTARStopping::GetIndex(), G4PSTARStopping::GetIndex(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4EmCalculator::GetKinEnergy(), GetMaterial(), G4HadronicInteraction::GetMaxEnergy(), G4PenelopeOscillatorManager::GetMeanExcitationEnergy(), G4OpRayleigh::GetMeanFreePath(), G4EmCalculator::GetMeanFreePath(), G4HadronicInteraction::GetMinEnergy(), G4PenelopeOscillatorManager::GetNumberOfZAtomsPerMolecule(), G4PenelopeOscillatorManager::GetOscillatorCompton(), G4PenelopeOscillatorManager::GetOscillatorIonisation(), G4PenelopeOscillatorManager::GetOscillatorTableCompton(), G4PenelopeOscillatorManager::GetOscillatorTableIonisation(), G4PenelopeOscillatorManager::GetPlasmaEnergySquared(), G4EmCalculator::GetRange(), G4EmCalculator::GetRangeFromRestricteDEDX(), G4PenelopeOscillatorManager::GetTotalA(), G4PenelopeOscillatorManager::GetTotalZ(), G4EmCorrections::HighOrderCorrections(), G4EmModelManager::Initialise(), G4PAIySection::Initialize(), G4GDMLWriteMaterials::MaterialWrite(), G4DNASecondOrderReaction::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4Decay::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4NistManager::PrintG4Material(), G4hImpactIonisation::PrintInfoDefinition(), GVFlashShowerParameterisation::PrintMaterial(), G4DNAMolecularMaterial::PrintNotAMolecularMaterial(), G4GDMLWriteMaterials::PropertyWrite(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4PenelopeBremsstrahlungAngular::SampleDirection(), G4PenelopeBremsstrahlungFS::SampleGammaEnergy(), G4UrbanMscModel90::SampleScattering(), G4UrbanMscModel92::SampleScattering(), G4UrbanMscModel93::SampleScattering(), G4UrbanMscModel95::SampleScattering(), G4UrbanMscModel96::SampleScattering(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4PenelopeBremsstrahlungModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4HadronicInteraction::SetMaxEnergy(), G4IonisParamMat::SetMeanExcitationEnergy(), G4HadronicInteraction::SetMinEnergy(), G4EnergySplitter::SplitEnergyInVolumes(), G4ProductionCutsTable::StoreMaterialCutsCoupleInfo(), G4GDMLRead::StripNames(), G4tgbMaterialMixtureByVolume::TransformToFractionsByWeight(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ GetNuclearInterLength()

G4double G4Material::GetNuclearInterLength ( ) const
inline

Definition at line 222 of file G4Material.hh.

222{return fNuclInterLen;}

Referenced by G4MSSteppingAction::UserSteppingAction().

◆ GetNumberOfElements()

size_t G4Material::GetNumberOfElements ( ) const
inline

Definition at line 185 of file G4Material.hh.

185{return fNumberOfElements;}

Referenced by G4VCrossSectionHandler::ActiveElements(), AddMaterial(), G4AdjointPhotoElectricModel::AdjointCrossSection(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4FissLib::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPorLEInelastic::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLFission::ApplyYourself(), G4NeutronHPThermalScattering::ApplyYourself(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4QCaptureAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4RToEConvForGamma::BuildAbsorptionLengthVector(), G4AugerData::BuildAugerTransitionTable(), G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4CrossSectionHandler::BuildCrossSectionsForMaterials(), G4NeutronHPThermalScatteringData::BuildPhysicsTable(), G4VRangeToEnergyConverter::BuildRangeVector(), G4Nucleus::ChooseParameters(), G4AdjointCSManager::ComputeAdjointCS(), G4LivermoreBremsstrahlungModel::ComputeDEDXPerVolume(), G4LivermoreIonisationModel::ComputeDEDXPerVolume(), G4MuBremsstrahlungModel::ComputeDEDXPerVolume(), G4MuPairProductionModel::ComputeDEDXPerVolume(), G4eBremParametrizedModel::ComputeDEDXPerVolume(), G4eBremsstrahlungModel::ComputeDEDXPerVolume(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4PAIySection::ComputeLowEnergyCof(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4eBremsstrahlungModel::CrossSectionPerVolume(), G4VEmModel::CrossSectionPerVolume(), G4AnnihiToMuPair::CrossSectionPerVolume(), G4AdjointBremsstrahlungModel::DiffCrossSectionPerVolumePrimToSecondApproximated2(), G4tgbGeometryDumper::DumpMaterial(), G4EmElementSelector::G4EmElementSelector(), G4Material(), G4ShellVacancy::GenerateNumberOfIonisations(), G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(), G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(), G4CrossSectionDataStore::GetCrossSection(), GVFlashShowerParameterisation::GetEffA(), GVFlashShowerParameterisation::GetEffZ(), G4HadronicProcessStore::GetElasticCrossSectionPerVolume(), G4StopElementSelector::GetElement(), G4HadronicProcessStore::GetFissionCrossSectionPerVolume(), G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(), G4NeutronIsotopeProduction::GetIsotope(), G4QAtomicElectronScattering::GetMeanFreePath(), G4QCoherentChargeExchange::GetMeanFreePath(), G4QDiffraction::GetMeanFreePath(), G4QElastic::GetMeanFreePath(), G4QInelastic::GetMeanFreePath(), G4QIonIonElastic::GetMeanFreePath(), G4QLowEnergy::GetMeanFreePath(), G4QNGamma::GetMeanFreePath(), G4KaonMinusAbsorptionAtRest::GetMeanLifeTime(), G4PiMinusAbsorptionAtRest::GetMeanLifeTime(), G4hICRU49He::HasMaterial(), G4hICRU49p::HasMaterial(), G4hSRIM2000p::HasMaterial(), G4hZiegler1977He::HasMaterial(), G4hZiegler1977p::HasMaterial(), G4hZiegler1985p::HasMaterial(), G4LivermoreGammaConversionModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4LivermoreRayleighModel::Initialise(), G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), G4QAOLowEnergyLoss::IsInCharge(), G4GDMLWriteMaterials::MaterialWrite(), G4QAtomicElectronScattering::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QNGamma::PostStepDoIt(), G4GoudsmitSaundersonMscModel::SampleScattering(), G4WentzelVIModel::SampleScattering(), G4WentzelVIRelModel::SampleScattering(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4CrossSectionDataStore::SampleZandA(), G4VEmModel::SelectRandomAtom(), G4PixeCrossSectionHandler::SelectRandomAtom(), G4eBremsstrahlungModel::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomElement(), G4ElementSelector::SelectZandA(), G4hICRU49He::StoppingPower(), G4hICRU49p::StoppingPower(), G4hQAOModel::StoppingPower(), G4hSRIM2000p::StoppingPower(), G4hZiegler1977He::StoppingPower(), G4hZiegler1977p::StoppingPower(), G4hZiegler1985p::StoppingPower(), G4VCrossSectionHandler::ValueForMaterial(), and G4PixeCrossSectionHandler::ValueForMaterial().

◆ GetNumberOfMaterials()

◆ GetPressure()

G4double G4Material::GetPressure ( ) const
inline

◆ GetRadlen()

◆ GetSandiaTable()

◆ GetState()

◆ GetTemperature()

G4double G4Material::GetTemperature ( ) const
inline

Definition at line 181 of file G4Material.hh.

181{return fTemp;}

Referenced by G4NeutronHPChannelList::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4LENDModel::ApplyYourself(), G4FissLib::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPorLEInelastic::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLFission::ApplyYourself(), G4NeutronHPThermalScattering::ApplyYourself(), G4FissionLibrary::ApplyYourself(), G4NeutronHPCaptureFS::ApplyYourself(), G4NeutronHPElasticFS::ApplyYourself(), G4NeutronHPFissionFS::ApplyYourself(), G4NeutronHPChannel::ApplyYourself(), G4NeutronHPInelasticBaseFS::BaseApply(), G4tgbMaterialMixtureByWeight::BuildG4Material(), G4NistManager::BuildMaterialWithNewDensity(), G4NeutronHPInelasticCompFS::CompositeApply(), G4tgbGeometryDumper::DumpMaterial(), G4Nucleus::G4Nucleus(), G4NeutronHPThermalScatteringData::GetCoherentCrossSection(), G4NeutronHPThermalScatteringData::GetCrossSection(), G4NeutronHPThermalScatteringData::GetIncoherentCrossSection(), G4NeutronHPThermalScatteringData::GetInelasticCrossSection(), G4LENDCrossSection::GetIsoCrossSection(), G4NeutronHPCaptureData::GetIsoCrossSection(), G4NeutronHPElasticData::GetIsoCrossSection(), G4NeutronHPFissionData::GetIsoCrossSection(), G4NeutronHPInelasticData::GetIsoCrossSection(), G4NeutronHPorLCaptureData::GetIsoCrossSection(), G4NeutronHPorLEInelasticData::GetIsoCrossSection(), G4NeutronHPorLElasticData::GetIsoCrossSection(), G4NeutronHPorLFissionData::GetIsoCrossSection(), and G4GDMLWriteMaterials::MaterialWrite().

◆ GetTotNbOfAtomsPerVolume()

◆ GetTotNbOfElectPerVolume()

G4double G4Material::GetTotNbOfElectPerVolume ( ) const
inline

Definition at line 211 of file G4Material.hh.

211{return TotNbOfElectPerVolume;}

Referenced by G4hICRU49He::StoppingPower().

◆ GetVecNbOfAtomsPerVolume()

const G4double * G4Material::GetVecNbOfAtomsPerVolume ( ) const
inline

Definition at line 205 of file G4Material.hh.

205{return VecNbOfAtomsPerVolume;}

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSection(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4FissLib::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPorLEInelastic::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLFission::ApplyYourself(), G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials(), G4Nucleus::ChooseParameters(), G4AdjointCSManager::ComputeAdjointCS(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4VEmModel::CrossSectionPerVolume(), G4AnnihiToMuPair::CrossSectionPerVolume(), G4HadronicProcessStore::GetCaptureCrossSectionPerVolume(), G4HadronicProcessStore::GetChargeExchangeCrossSectionPerVolume(), G4CrossSectionDataStore::GetCrossSection(), G4HadronicProcessStore::GetElasticCrossSectionPerVolume(), G4HadronicProcessStore::GetFissionCrossSectionPerVolume(), G4HadronicProcessStore::GetInelasticCrossSectionPerVolume(), G4QAtomicElectronScattering::GetMeanFreePath(), G4QCoherentChargeExchange::GetMeanFreePath(), G4QDiffraction::GetMeanFreePath(), G4QElastic::GetMeanFreePath(), G4QInelastic::GetMeanFreePath(), G4QIonIonElastic::GetMeanFreePath(), G4QLowEnergy::GetMeanFreePath(), G4QNGamma::GetMeanFreePath(), G4EmElementSelector::Initialise(), G4GoudsmitSaundersonMscModel::SampleScattering(), G4VCrossSectionHandler::ValueForMaterial(), and G4PixeCrossSectionHandler::ValueForMaterial().

◆ GetZ()

G4double G4Material::GetZ ( ) const

◆ operator!=()

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

Definition at line 689 of file G4Material.cc.

690{
691 return (this != (G4Material *) &right);
692}

◆ operator==()

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

Definition at line 682 of file G4Material.cc.

683{
684 return (this == (G4Material *) &right);
685}

◆ SetChemicalFormula()

void G4Material::SetChemicalFormula ( const G4String chF)
inline

Definition at line 172 of file G4Material.hh.

172{fChemicalFormula=chF;}

◆ SetMaterialPropertiesTable()

void G4Material::SetMaterialPropertiesTable ( G4MaterialPropertiesTable anMPT)
inline

Definition at line 248 of file G4Material.hh.

249 {fMaterialPropertiesTable = anMPT;}

Referenced by G4GDMLReadMaterials::PropertyRead().

◆ SetName()

void G4Material::SetName ( const G4String name)
inline

Definition at line 282 of file G4Material.hh.

282{fName=name;}

Referenced by G4GDMLRead::StripNames().

Friends And Related Function Documentation

◆ operator<< [1/3]

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

Definition at line 740 of file G4Material.cc.

741{
742 flux << &material;
743 return flux;
744}

◆ operator<< [2/3]

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

Definition at line 697 of file G4Material.cc.

698{
699 std::ios::fmtflags mode = flux.flags();
700 flux.setf(std::ios::fixed,std::ios::floatfield);
701 G4long prec = flux.precision(3);
702
703 flux
704 << " Material: " << std::setw(8) << material->fName
705 << " " << material->fChemicalFormula << " "
706 << " density: " << std::setw(6) << std::setprecision(3)
707 << G4BestUnit(material->fDensity,"Volumic Mass")
708 << " RadL: " << std::setw(7) << std::setprecision(3)
709 << G4BestUnit(material->fRadlen,"Length")
710 << " Nucl.Int.Length: " << std::setw(7) << std::setprecision(3)
711 << G4BestUnit(material->fNuclInterLen,"Length")
712 << " Imean: " << std::setw(7) << std::setprecision(3)
713 << G4BestUnit(material->GetIonisation()->GetMeanExcitationEnergy(),"Energy");
714
715 if(material->fState == kStateGas) {
716 flux
717 << " temperature: " << std::setw(6) << std::setprecision(2)
718 << (material->fTemp)/kelvin << " K"
719 << " pressure: " << std::setw(6) << std::setprecision(2)
720 << (material->fPressure)/atmosphere << " atm";
721 }
722 for (size_t i=0; i<material->fNumberOfElements; i++) {
723 flux
724 << "\n ---> " << (*(material->theElementVector))[i]
725 << "\n ElmMassFraction: "
726 << std::setw(6)<< std::setprecision(2)
727 << (material->fMassFractionVector[i])/perCent << " %"
728 << " ElmAbundance " << std::setw(6)<< std::setprecision(2)
729 << 100*(material->VecNbOfAtomsPerVolume[i])/(material->TotNbOfAtomsPerVolume)
730 << " % \n";
731 }
732 flux.precision(prec);
733 flux.setf(mode,std::ios::floatfield);
734
735 return flux;
736}
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
long G4long
Definition: G4Types.hh:68
G4double GetMeanExcitationEnergy() const
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225

◆ operator<< [3/3]

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

Definition at line 748 of file G4Material.cc.

749{
750 //Dump info for all known materials
751 flux << "\n***** Table : Nb of materials = " << MaterialTable.size()
752 << " *****\n" << G4endl;
753
754 for (size_t i=0; i<MaterialTable.size(); ++i) {
755 flux << MaterialTable[i] << G4endl << G4endl;
756 }
757
758 return flux;
759}

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