Geant4-11
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
G4ExtendedMaterial Class Reference

#include <G4ExtendedMaterial.hh>

Inheritance diagram for G4ExtendedMaterial:
G4Material

Public Member Functions

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

Static Public Member Functions

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

Private Member Functions

void ComputeDerivedQuantities ()
 
void ComputeNuclearInterLength ()
 
void ComputeRadiationLength ()
 
void CopyPointersOfBaseMaterial ()
 
void FillVectors ()
 
void InitializePointers ()
 

Private Attributes

std::vector< G4int > * fAtoms
 
G4intfAtomsVector
 
const G4MaterialfBaseMaterial
 
G4String fChemicalFormula
 
G4double fDensity
 
std::vector< const G4Element * > * fElm
 
std::vector< G4double > * fElmFrac
 
G4MaterialExtensionMap fExtensionMap
 
G4double fFreeElecDensity
 
G4int fIdxComponent
 
size_t fIndexInTable
 
G4IonisParamMatfIonisation
 
G4bool fMassFraction
 
G4doublefMassFractionVector
 
G4double fMassOfMolecule
 
std::map< G4Material *, G4doublefMatComponents
 
G4MaterialPropertiesTablefMaterialPropertiesTable
 
G4String fName
 
G4int fNbComponents
 
G4double fNuclInterLen
 
G4int fNumberOfElements
 
G4double fPressure
 
G4double fRadlen
 
G4SandiaTablefSandiaTable
 
G4State fState
 
G4double fTemp
 
G4double fTotNbOfAtomsPerVolume
 
G4double fTotNbOfElectPerVolume
 
G4doublefVecNbOfAtomsPerVolume
 
G4ElementVectortheElementVector
 

Static Private Attributes

static G4MaterialTable theMaterialTable
 

Detailed Description

Definition at line 61 of file G4ExtendedMaterial.hh.

Constructor & Destructor Documentation

◆ G4ExtendedMaterial() [1/4]

G4ExtendedMaterial::G4ExtendedMaterial ( const G4String name,
const G4Material baseMaterial 
)

Definition at line 43 of file G4ExtendedMaterial.cc.

45 : G4Material(name,baseMaterial->GetDensity(),baseMaterial,
46 baseMaterial->GetState(),baseMaterial->GetTemperature(),
47 baseMaterial->GetPressure())
48{;}
G4double GetPressure() const
Definition: G4Material.hh:179
G4double GetDensity() const
Definition: G4Material.hh:176
G4State GetState() const
Definition: G4Material.hh:177
G4double GetTemperature() const
Definition: G4Material.hh:178
G4Material(const G4String &name, G4double z, G4double a, G4double density, G4State state=kStateUndefined, G4double temp=NTP_Temperature, G4double pressure=CLHEP::STP_Pressure)
Definition: G4Material.cc:95
const char * name(G4int ptype)

◆ G4ExtendedMaterial() [2/4]

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

Definition at line 52 of file G4ExtendedMaterial.cc.

55 : G4Material(name,z,a,density,state,temp,pressure)
56{;}

◆ G4ExtendedMaterial() [3/4]

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

Definition at line 63 of file G4ExtendedMaterial.cc.

66 : G4Material(name,density,nComponents,state,temp,pressure)
67{;}

◆ G4ExtendedMaterial() [4/4]

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

Definition at line 73 of file G4ExtendedMaterial.cc.

76 : G4Material(name,density,bmat,state,temp,pressure)
77{;}

◆ ~G4ExtendedMaterial()

virtual G4ExtendedMaterial::~G4ExtendedMaterial ( )
inlinevirtual

Definition at line 103 of file G4ExtendedMaterial.hh.

103{};

Member Function Documentation

◆ AddElement() [1/2]

void G4Material::AddElement ( G4Element elm,
G4double  frac 
)
inlineinherited

Definition at line 164 of file G4Material.hh.

165 { AddElementByMassFraction(elm, frac); }
void AddElementByMassFraction(const G4Element *elm, G4double fraction)
Definition: G4Material.cc:434

References G4Material::AddElementByMassFraction().

◆ AddElement() [2/2]

void G4Material::AddElement ( G4Element elm,
G4int  nAtoms 
)
inlineinherited

◆ AddElementByMassFraction()

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

Definition at line 434 of file G4Material.cc.

435{
436 // perform checks consistency
437 if(fraction < 0.0 || fraction > 1.0) {
439 ed << "For material " << fName << " and added element "
440 << elm->GetName() << " massFraction= " << fraction
441 << " is wrong ";
442 G4Exception ("G4Material::AddElementByMassFraction()", "mat031",
443 FatalException, ed, "");
444 }
445 if(!fMassFraction) {
447 ed << "For material " << fName << " and added element "
448 << elm->GetName() << ", massFraction= " << fraction
449 << ", fIdxComponent=" << fIdxComponent
450 << " problem: cannot add by mass fraction after "
451 << "addition of elements by number of atoms";
452 G4Exception ("G4Material::AddElementByMassFraction()", "mat031",
453 FatalException, ed, "");
454 }
457 ed << "For material " << fName << " and added element "
458 << elm->GetName() << ", massFraction= " << fraction
459 << ", fIdxComponent=" << fIdxComponent
460 << "; attempt to add more than the declared number of components "
461 << fIdxComponent << " >= " << fNbComponents;
462 G4Exception ("G4Material::AddElementByMassFraction()", "mat031",
463 FatalException, ed, "");
464 }
465 if(0 == fIdxComponent) {
466 fElmFrac = new std::vector<G4double>;
467 fElm = new std::vector<const G4Element*>;
468 }
469
470 // filling
471 G4bool isAdded = false;
472 if(!fElm->empty()) {
473 for (G4int i=0; i<fNumberOfElements; ++i) {
474 if ( elm == (*fElm)[i] ) {
475 (*fElmFrac)[i] += fraction;
476 isAdded = true;
477 break;
478 }
479 }
480 }
481 if(!isAdded) {
482 fElm->push_back(elm);
483 fElmFrac->push_back(fraction);
485 }
487
488 // is filled
490}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4String & GetName() const
Definition: G4Element.hh:127
G4int fIdxComponent
Definition: G4Material.hh:351
G4int fNbComponents
Definition: G4Material.hh:350
G4bool fMassFraction
Definition: G4Material.hh:352
G4int fNumberOfElements
Definition: G4Material.hh:347
void FillVectors()
Definition: G4Material.cc:562
std::vector< G4double > * fElmFrac
Definition: G4Material.hh:356
G4String fName
Definition: G4Material.hh:362
std::vector< const G4Element * > * fElm
Definition: G4Material.hh:357

References FatalException, G4Material::fElm, G4Material::fElmFrac, G4Material::fIdxComponent, G4Material::FillVectors(), G4Material::fMassFraction, G4Material::fName, G4Material::fNbComponents, G4Material::fNumberOfElements, G4Exception(), and G4Element::GetName().

Referenced by G4Material::AddElement().

◆ AddElementByNumberOfAtoms()

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

Definition at line 360 of file G4Material.cc.

361{
362 // perform checks consistency
363 if(0 == fIdxComponent) {
364 fMassFraction = false;
365 fAtoms = new std::vector<G4int>;
366 fElm = new std::vector<const G4Element*>;
367 }
370 ed << "For material " << fName << " and added element "
371 << elm->GetName() << " with Natoms=" << nAtoms
372 << " wrong attempt to add more than the declared number of elements "
373 << fIdxComponent << " >= " << fNbComponents;
374 G4Exception ("G4Material::AddElementByNumberOfAtoms()", "mat031",
375 FatalException, ed, "");
376 }
377 if(fMassFraction) {
379 ed << "For material " << fName << " and added element "
380 << elm->GetName() << " with Natoms=" << nAtoms
381 << " problem: cannot add by number of atoms after "
382 << "addition of elements by mass fraction";
383 G4Exception ("G4Material::AddElementByNumberOfAtoms()", "mat031",
384 FatalException, ed, "");
385 }
386 // filling
387 G4bool isAdded = false;
388 if(!fElm->empty()) {
389 for (G4int i=0; i<fNumberOfElements; ++i) {
390 if ( elm == (*fElm)[i] ) {
391 (*fAtoms)[i] += nAtoms;
392 isAdded = true;
393 break;
394 }
395 }
396 }
397 if(!isAdded) {
398 fElm->push_back(elm);
399 fAtoms->push_back(nAtoms);
401 }
403
404 // is filled - complete composition of atoms
410
411 G4double Amol = 0.;
412 for (G4int i=0; i<fNumberOfElements; ++i) {
413 theElementVector->push_back((*fElm)[i]);
414 fAtomsVector[i] = (*fAtoms)[i];
415 G4double w = fAtomsVector[i]*(*fElm)[i]->GetA();
416 Amol += w;
417 fMassFractionVector[i] = w;
418 }
419 for (G4int i=0; i<fNumberOfElements; ++i) {
420 fMassFractionVector[i] /= Amol;
421 }
422 delete fAtoms;
423 delete fElm;
426 }
427}
std::vector< const G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:83
void ComputeDerivedQuantities()
Definition: G4Material.cc:288
G4double * fMassFractionVector
Definition: G4Material.hh:328
std::vector< G4int > * fAtoms
Definition: G4Material.hh:355
G4ElementVector * theElementVector
Definition: G4Material.hh:326
G4double fMassOfMolecule
Definition: G4Material.hh:343
G4int * fAtomsVector
Definition: G4Material.hh:327
float Avogadro
Definition: hepunit.py:252

References source.hepunit::Avogadro, G4Material::ComputeDerivedQuantities(), FatalException, G4Material::fAtoms, G4Material::fAtomsVector, G4Material::fElm, G4Material::fIdxComponent, G4Material::fMassFraction, G4Material::fMassFractionVector, G4Material::fMassOfMolecule, G4Material::fName, G4Material::fNbComponents, G4Material::fNumberOfElements, G4Exception(), G4Element::GetName(), and G4Material::theElementVector.

Referenced by G4Material::AddElement().

◆ AddMaterial()

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

Definition at line 496 of file G4Material.cc.

497{
498 if(fraction < 0.0 || fraction > 1.0) {
500 ed << "For material " << fName << " and added material "
501 << material->GetName() << ", massFraction= " << fraction
502 << " is wrong ";
503 G4Exception ("G4Material::AddMaterial()", "mat031", FatalException,
504 ed, "");
505 }
506 if(!fMassFraction) {
508 ed << "For material " << fName << " and added material "
509 << material->GetName() << ", massFraction= " << fraction
510 << ", fIdxComponent=" << fIdxComponent
511 << " problem: cannot add by mass fraction after "
512 << "addition of elements by number of atoms";
513 G4Exception ("G4Material::AddMaterial()", "mat031", FatalException,
514 ed, "");
515 }
518 ed << "For material " << fName << " and added material "
519 << material->GetName() << ", massFraction= " << fraction
520 << "; attempt to add more than the declared number of components "
521 << fIdxComponent << " >= " << fNbComponents;
522 G4Exception ("G4Material::AddMaterial()", "mat031", FatalException,
523 ed, "");
524 }
525
526 if(0 == fIdxComponent) {
527 fElmFrac = new std::vector<G4double>;
528 fElm = new std::vector<const G4Element*>;
529 }
530
531 // filling
532 G4int nelm = material->GetNumberOfElements();
533 for(G4int j=0; j<nelm; ++j) {
534 auto elm = material->GetElement(j);
535 auto frac = material->GetFractionVector();
536 G4bool isAdded = false;
537 if(!fElm->empty()) {
538 for (G4int i=0; i<fNumberOfElements; ++i) {
539 if ( elm == (*fElm)[i] ) {
540 (*fElmFrac)[i] += fraction*frac[j];
541 isAdded = true;
542 break;
543 }
544 }
545 }
546 if(!isAdded) {
547 fElm->push_back(elm);
548 fElmFrac->push_back(fraction*frac[j]);
550 }
551 }
552
553 fMatComponents[material] = fraction;
555
556 // is filled
558}
std::map< G4Material *, G4double > fMatComponents
Definition: G4Material.hh:360
string material
Definition: eplot.py:19

References FatalException, G4Material::fElm, G4Material::fElmFrac, G4Material::fIdxComponent, G4Material::FillVectors(), G4Material::fMassFraction, G4Material::fMatComponents, G4Material::fName, G4Material::fNbComponents, G4Material::fNumberOfElements, G4Exception(), and eplot::material.

Referenced by G4tgbMaterialMixtureByVolume::BuildG4Material(), G4tgbMaterialMixtureByWeight::BuildG4Material(), DetectorConstruction::DefineMaterials(), ExN03DetectorConstruction::DefineMaterials(), and export_G4Material().

◆ begin()

G4MaterialExtensionMap::const_iterator G4ExtendedMaterial::begin ( ) const
inline

Definition at line 124 of file G4ExtendedMaterial.hh.

124{ return fExtensionMap.begin(); }
G4MaterialExtensionMap fExtensionMap

References fExtensionMap.

◆ cbegin()

G4MaterialExtensionMap::const_iterator G4ExtendedMaterial::cbegin ( ) const
inline

Definition at line 125 of file G4ExtendedMaterial.hh.

125{ return fExtensionMap.cbegin(); }

References fExtensionMap.

◆ cend()

G4MaterialExtensionMap::const_iterator G4ExtendedMaterial::cend ( ) const
inline

Definition at line 127 of file G4ExtendedMaterial.hh.

127{ return fExtensionMap.cend(); }

References fExtensionMap.

◆ ComputeDensityEffectOnFly()

void G4Material::ComputeDensityEffectOnFly ( G4bool  val)
inherited

Definition at line 658 of file G4Material.cc.

659{
660#ifdef G4MULTITHREADED
661 G4MUTEXLOCK(&materialMutex);
662#endif
663 if (nullptr == fIonisation) { fIonisation = new G4IonisParamMat(this); }
665#ifdef G4MULTITHREADED
666 G4MUTEXUNLOCK(&materialMutex);
667#endif
668}
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
void ComputeDensityEffectOnFly(G4bool)
G4IonisParamMat * fIonisation
Definition: G4Material.hh:331

References G4IonisParamMat::ComputeDensityEffectOnFly(), G4Material::fIonisation, G4MUTEXLOCK, and G4MUTEXUNLOCK.

Referenced by G4Material::CopyPointersOfBaseMaterial(), and G4NistManager::SetDensityEffectCalculatorFlag().

◆ ComputeDerivedQuantities()

void G4Material::ComputeDerivedQuantities ( )
privateinherited

Definition at line 288 of file G4Material.cc.

289{
290 // Header routine to compute various properties of material.
291 //
292 // Number of atoms per volume (per element), total nb of electrons per volume
293 G4double Zi, Ai;
295 delete [] fVecNbOfAtomsPerVolume;
298 fFreeElecDensity = 0.;
299 const G4double elecTh = 15.*CLHEP::eV; // threshold for conductivity e-
300 for (G4int i=0; i<fNumberOfElements; ++i) {
301 Zi = (*theElementVector)[i]->GetZ();
302 Ai = (*theElementVector)[i]->GetA();
306 if(fState != kStateGas) {
309 }
310 }
311
314
315 if (!fIonisation) { fIonisation = new G4IonisParamMat(this); }
316 if (!fSandiaTable){ fSandiaTable = new G4SandiaTable(this); }
317}
@ kStateGas
Definition: G4Material.hh:111
static G4int GetNumberOfFreeElectrons(G4int Z, G4double th)
G4double fTotNbOfElectPerVolume
Definition: G4Material.hh:340
G4double fFreeElecDensity
Definition: G4Material.hh:335
G4State fState
Definition: G4Material.hh:345
G4double fDensity
Definition: G4Material.hh:334
G4SandiaTable * fSandiaTable
Definition: G4Material.hh:332
G4double fTotNbOfAtomsPerVolume
Definition: G4Material.hh:339
G4double * fVecNbOfAtomsPerVolume
Definition: G4Material.hh:329
void ComputeNuclearInterLength()
Definition: G4Material.cc:612
void ComputeRadiationLength()
Definition: G4Material.cc:601
static constexpr double eV

References source.hepunit::Avogadro, G4Material::ComputeNuclearInterLength(), G4Material::ComputeRadiationLength(), CLHEP::eV, G4Material::fDensity, G4Material::fFreeElecDensity, G4Material::fIonisation, G4Material::fMassFractionVector, G4Material::fNumberOfElements, G4Material::fSandiaTable, G4Material::fState, G4Material::fTotNbOfAtomsPerVolume, G4Material::fTotNbOfElectPerVolume, G4Material::fVecNbOfAtomsPerVolume, G4AtomicShells::GetNumberOfFreeElectrons(), and kStateGas.

Referenced by G4Material::AddElementByNumberOfAtoms(), G4Material::FillVectors(), and G4Material::G4Material().

◆ ComputeNuclearInterLength()

void G4Material::ComputeNuclearInterLength ( )
privateinherited

Definition at line 612 of file G4Material.cc.

613{
614 const G4double lambda0 = 35*CLHEP::g/CLHEP::cm2;
615 const G4double twothird = 2.0/3.0;
616 G4double NILinv = 0.0;
617 for (G4int i=0; i<fNumberOfElements; ++i) {
618 G4int Z = (*theElementVector)[i]->GetZasInt();
619 G4double A = (*theElementVector)[i]->GetN();
620 if(1 == Z) {
621 NILinv += fVecNbOfAtomsPerVolume[i]*A;
622 } else {
623 NILinv += fVecNbOfAtomsPerVolume[i]*G4Exp(twothird*G4Log(A));
624 }
625 }
626 NILinv *= amu/lambda0;
627 fNuclInterLen = (NILinv <= 0.0 ? DBL_MAX : 1./NILinv);
628}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
const G4int Z[17]
const G4double A[17]
G4double fNuclInterLen
Definition: G4Material.hh:342
static constexpr double g
static constexpr double cm2
#define DBL_MAX
Definition: templates.hh:62

References A, source.hepunit::amu, CLHEP::cm2, DBL_MAX, G4Material::fNuclInterLen, G4Material::fNumberOfElements, G4Material::fVecNbOfAtomsPerVolume, CLHEP::g, G4Exp(), G4Log(), and Z.

Referenced by G4Material::ComputeDerivedQuantities().

◆ ComputeRadiationLength()

void G4Material::ComputeRadiationLength ( )
privateinherited

Definition at line 601 of file G4Material.cc.

602{
603 G4double radinv = 0.0 ;
604 for (G4int i=0; i<fNumberOfElements; ++i) {
605 radinv += fVecNbOfAtomsPerVolume[i]*((*theElementVector)[i]->GetfRadTsai());
606 }
607 fRadlen = (radinv <= 0.0 ? DBL_MAX : 1./radinv);
608}
G4double fRadlen
Definition: G4Material.hh:341

References DBL_MAX, G4Material::fNumberOfElements, G4Material::fRadlen, and G4Material::fVecNbOfAtomsPerVolume.

Referenced by G4Material::ComputeDerivedQuantities().

◆ CopyPointersOfBaseMaterial()

void G4Material::CopyPointersOfBaseMaterial ( )
privateinherited

Definition at line 321 of file G4Material.cc.

322{
327
329
335
339 for (G4int i=0; i<fNumberOfElements; ++i) {
340 fVecNbOfAtomsPerVolume[i] = factor*v[i];
341 }
342 fRadlen = fBaseMaterial->GetRadlen()/factor;
344
345 if(!fIonisation) { fIonisation = new G4IonisParamMat(this); }
349 }
350
353}
@ kStateUndefined
Definition: G4Material.hh:111
G4double GetMeanExcitationEnergy() const
void SetMeanExcitationEnergy(G4double value)
G4DensityEffectCalculator * GetDensityEffectCalculator()
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
const G4Material * fBaseMaterial
Definition: G4Material.hh:318
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
const G4double * GetFractionVector() const
Definition: G4Material.hh:190
G4double GetTotNbOfElectPerVolume() const
Definition: G4Material.hh:208
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:222
G4double GetFreeElectronDensity() const
Definition: G4Material.hh:175
const G4int * GetAtomsVector() const
Definition: G4Material.hh:194
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:225
G4double GetRadlen() const
Definition: G4Material.hh:216
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
void ComputeDensityEffectOnFly(G4bool)
Definition: G4Material.cc:658
G4double GetNuclearInterLength() const
Definition: G4Material.hh:219
G4MaterialPropertiesTable * fMaterialPropertiesTable
Definition: G4Material.hh:319

References G4Material::ComputeDensityEffectOnFly(), G4Material::fAtomsVector, G4Material::fBaseMaterial, G4Material::fDensity, G4Material::fFreeElecDensity, G4Material::fIonisation, G4Material::fMassFractionVector, G4Material::fMaterialPropertiesTable, G4Material::fNuclInterLen, G4Material::fNumberOfElements, G4Material::fRadlen, G4Material::fSandiaTable, G4Material::fState, G4Material::fTotNbOfAtomsPerVolume, G4Material::fTotNbOfElectPerVolume, G4Material::fVecNbOfAtomsPerVolume, G4Material::GetAtomsVector(), G4Material::GetDensity(), G4IonisParamMat::GetDensityEffectCalculator(), G4Material::GetElementVector(), G4Material::GetFractionVector(), G4Material::GetFreeElectronDensity(), G4Material::GetIonisation(), G4Material::GetMaterialPropertiesTable(), G4IonisParamMat::GetMeanExcitationEnergy(), G4Material::GetNuclearInterLength(), G4Material::GetRadlen(), G4Material::GetSandiaTable(), G4Material::GetState(), G4Material::GetTotNbOfAtomsPerVolume(), G4Material::GetTotNbOfElectPerVolume(), G4Material::GetVecNbOfAtomsPerVolume(), kStateUndefined, G4IonisParamMat::SetMeanExcitationEnergy(), and G4Material::theElementVector.

Referenced by G4Material::G4Material().

◆ end()

G4MaterialExtensionMap::const_iterator G4ExtendedMaterial::end ( ) const
inline

Definition at line 126 of file G4ExtendedMaterial.hh.

126{ return fExtensionMap.end(); }

References fExtensionMap.

◆ FillVectors()

void G4Material::FillVectors ( )
privateinherited

Definition at line 562 of file G4Material.cc.

563{
564 // there are material components
569 G4double wtSum(0.0);
570 for (G4int i=0; i<fNumberOfElements; ++i) {
571 theElementVector->push_back((*fElm)[i]);
572 fMassFractionVector[i] = (*fElmFrac)[i];
573 wtSum += fMassFractionVector[i];
574 }
575 delete fElmFrac;
576 delete fElm;
577
578 // check sum of weights -- OK?
579 if (std::abs(1.-wtSum) > perThousand) {
581 ed << "For material " << fName << " sum of fractional masses "
582 << wtSum << " is not 1 - results may be wrong";
583 G4Exception ("G4Material::FillVectors()", "mat031", JustWarning,
584 ed, "");
585 }
586 G4double coeff = (wtSum > 0.0) ? 1./wtSum : 1.0;
587 G4double Amol(0.);
588 for (G4int i=0; i<fNumberOfElements; ++i) {
589 fMassFractionVector[i] *= coeff;
590 Amol += fMassFractionVector[i]*(*theElementVector)[i]->GetA();
591 }
592 for (G4int i=0; i<fNumberOfElements; ++i) {
593 fAtomsVector[i] =
595 }
597}
@ JustWarning
static constexpr double perThousand
Definition: G4SIunits.hh:326
G4double GetA() const
Definition: G4Material.cc:750
int G4lrint(double ad)
Definition: templates.hh:134

References G4Material::ComputeDerivedQuantities(), G4Material::fAtomsVector, G4Material::fElm, G4Material::fElmFrac, G4Material::fMassFractionVector, G4Material::fName, G4Material::fNumberOfElements, G4Exception(), G4lrint(), G4Material::GetA(), JustWarning, perThousand, and G4Material::theElementVector.

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

◆ GetA()

G4double G4Material::GetA ( ) const
inherited

Definition at line 750 of file G4Material.cc.

751{
752 if (fNumberOfElements > 1) {
754 ed << "For material " << fName << " ERROR in GetA() - Nelm="
755 << fNumberOfElements << " > 1, which is not allowed";
756 G4Exception ("G4Material::GetA()", "mat036", FatalException,
757 ed, "");
758 }
759 return (*theElementVector)[0]->GetA();
760}

References FatalException, G4Material::fName, G4Material::fNumberOfElements, G4Exception(), and G4Material::theElementVector.

Referenced by G4tgbGeometryDumper::DumpMaterial(), export_G4Material(), G4Material::FillVectors(), GVFlashShowerParameterisation::GetEffA(), G4Material::GetMaterial(), and G4GDMLWriteMaterials::MaterialWrite().

◆ GetAtomicNumDensityVector()

const G4double * G4Material::GetAtomicNumDensityVector ( ) const
inlineinherited

◆ GetAtomsVector()

const G4int * G4Material::GetAtomsVector ( ) const
inlineinherited

Definition at line 194 of file G4Material.hh.

194{return fAtomsVector;}

References G4Material::fAtomsVector.

Referenced by G4Material::CopyPointersOfBaseMaterial(), and G4VLEPTSModel::ReadParam().

◆ GetBaseMaterial()

const G4Material * G4Material::GetBaseMaterial ( ) const
inlineinherited

◆ GetChemicalFormula()

const G4String & G4Material::GetChemicalFormula ( ) const
inlineinherited

◆ GetDensity()

G4double G4Material::GetDensity ( ) const
inlineinherited

Definition at line 176 of file G4Material.hh.

176{return fDensity;}

References G4Material::fDensity.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4NistManager::BuildMaterialWithNewDensity(), G4VLEPTSModel::BuildMeanFreePathTable(), G4ProductionCutsTable::CheckMaterialInfo(), G4AdjointPrimaryGenerator::ComputeAccumulatedDepthVectorAlongBackRay(), G4EmCalculator::ComputeDEDX(), G4IonisParamMat::ComputeDensityEffectParameters(), G4TablesForExtrapolator::ComputeElectronDEDX(), G4SandiaTable::ComputeMatSandiaMatrixPAI(), G4eDPWAElasticDCS::ComputeMParams(), G4TablesForExtrapolator::ComputeMuonDEDX(), G4EmCalculator::ComputeNuclearDEDX(), G4TablesForExtrapolator::ComputeProtonDEDX(), G4NistMaterialBuilder::ConstructNewGasMaterial(), G4Material::CopyPointersOfBaseMaterial(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G4PhysicalVolumeModel::DescribeAndDescend(), G4AtimaFluctuations::Dispersion(), G4EmSaturation::DumpBirksCoefficients(), G4tgbGeometryDumper::DumpMaterial(), export_G4Material(), G4InitXscPAI::G4InitXscPAI(), G4PAIxSection::G4PAIxSection(), G4EmCalculator::GetDEDX(), G4VMscModel::GetEnergy(), G4LogicalVolume::GetMass(), G4Material::GetMaterial(), G4VMscModel::GetRange(), G4GSPWACorrections::InitDataMaterial(), G4DNAMolecularMaterial::InitializeDensity(), G4GSMottCorrection::InitMCDataMaterial(), G4GoudsmitSaundersonTable::InitMoliereMSCParams(), G4GDMLWriteMaterials::MaterialWrite(), G4PSDoseDeposit::ProcessHits(), G4PhysicalVolumeMassScene::ProcessVolume(), G4ErrorFreeTrajState::PropagateErrorIoni(), G4ErrorFreeTrajState::PropagateErrorMSC(), G4ASCIITreeSceneHandler::RequestPrimitives(), G4IonisParamMat::SetDensityEffectParameters(), GFlashHomoShowerParameterisation::SetMaterial(), GFlashSamplingShowerParameterisation::SetMaterial(), and G4tgbMaterialMixtureByVolume::TransformToFractionsByWeight().

◆ GetElectronDensity()

G4double G4Material::GetElectronDensity ( ) const
inlineinherited

◆ GetElement()

const G4Element * G4Material::GetElement ( G4int  iel) const
inlineinherited

◆ GetElementVector()

const G4ElementVector * G4Material::GetElementVector ( ) const
inlineinherited

Definition at line 186 of file G4Material.hh.

186{return theElementVector;}

References G4Material::theElementVector.

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSection(), G4SBBremTable::BuildSamplingTables(), G4Nucleus::ChooseParameters(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4IonisParamMat::ComputeDensityEffectParameters(), G4WentzelVIRelModel::ComputeEffectiveMass(), G4IonisParamMat::ComputeFluctModel(), G4IonisParamMat::ComputeIonParameters(), G4SandiaTable::ComputeMatSandiaMatrix(), G4SandiaTable::ComputeMatSandiaMatrixPAI(), G4SandiaTable::ComputeMatTable(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4IonisParamMat::ComputeMeanParameters(), G4eDPWAElasticDCS::ComputeMParams(), G4WentzelVIModel::ComputeSecondMoment(), G4CrystalExtension::ComputeStructureFactor(), G4CrystalExtension::ComputeStructureFactorGeometrical(), G4WentzelVIModel::ComputeTransportXSectionPerVolume(), G4Material::CopyPointersOfBaseMaterial(), G4tgbGeometryDumper::DumpMaterial(), export_G4Material(), G4EmElementSelector::G4EmElementSelector(), G4HadElementSelector::G4HadElementSelector(), G4CrystalExtension::GetAtomPos(), G4CrossSectionDataStore::GetCrossSection(), G4IonICRU73Data::GetDEDX(), G4MuNeutrinoNucleusTotXsc::GetElementCrossSection(), G4GSPWACorrections::InitDataMaterial(), G4GSPWACorrections::InitDataPerElement(), G4LivermoreGammaConversion5DModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4EmSaturation::InitialiseBirksCoefficient(), G4ElasticHadrNucleusHE::InitialiseModel(), G4GSMottCorrection::InitMCDataMaterial(), G4GSMottCorrection::InitMCDataPerElement(), G4GoudsmitSaundersonTable::InitMoliereMSCParams(), G4IonICRU73Data::ReadElementData(), G4WentzelVIModel::SampleScattering(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4GammaConversionToMuons::SelectRandomAtom(), G4ElementSelector::SelectZandA(), and G4EmCorrections::SetupKinematics().

◆ GetFractionVector()

const G4double * G4Material::GetFractionVector ( ) const
inlineinherited

◆ GetFreeElectronDensity()

G4double G4Material::GetFreeElectronDensity ( ) const
inlineinherited

◆ GetIndex()

size_t G4Material::GetIndex ( ) const
inlineinherited

Definition at line 256 of file G4Material.hh.

256{return fIndexInTable;}
size_t fIndexInTable
Definition: G4Material.hh:346

References G4Material::fIndexInTable.

Referenced by G4NistMaterialBuilder::BuildMaterial(), G4VLEPTSModel::BuildMeanFreePathTable(), G4VLEPTSModel::BuildPhysicsTable(), G4AdjointCSManager::ComputeAdjointCS(), G4EnergyLossForExtrapolator::ComputeDEDX(), G4EnergyLossForExtrapolator::ComputeEnergy(), G4EnergyLossForExtrapolator::ComputeRange(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::CrossSectionPerVolume(), export_G4Material(), G4DNAElectronHoleRecombination::FindReactant(), G4InitXscPAI::G4InitXscPAI(), G4PAIxSection::G4PAIxSection(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4IonICRU73Data::GetDEDX(), G4EnergyLossTables::GetDEDX(), G4EnergyLossTables::GetDeltaLabTime(), G4EnergyLossTables::GetDeltaProperTime(), G4EnergyLossTables::GetLabTime(), G4VLEPTSModel::GetMeanFreePath(), G4OpRayleigh::GetMeanFreePath(), G4DNADummyModel::GetNumMoleculePerVolumeUnitForMaterial(), G4DNAModelInterface::GetNumMoleculePerVolumeUnitForMaterial(), G4DNAModelInterface::GetNumMolPerVolUnitForComponentInComposite(), G4EnergyLossTables::GetPreciseDEDX(), G4EnergyLossTables::GetPreciseEnergyFromRange(), G4EnergyLossTables::GetPreciseRangeFromEnergy(), G4EnergyLossTables::GetProperTime(), G4EnergyLossTables::GetRange(), G4GoudsmitSaundersonMscModel::GetTransportMeanFreePath(), G4GoudsmitSaundersonMscModel::GetTransportMeanFreePathOnly(), G4GSPWACorrections::InitDataMaterial(), G4GSPWACorrections::InitDataPerMaterials(), G4IonICRU73Data::Initialise(), G4EmSaturation::InitialiseBirksCoefficient(), G4GSMottCorrection::InitMCDataMaterial(), G4GSMottCorrection::InitMCDataPerMaterials(), G4GoudsmitSaundersonTable::InitMoliereMSCParams(), G4GoudsmitSaundersonTable::InitSCPCorrection(), G4Scintillation::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpWLS2::PostStepDoIt(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4DNAMolecularMaterial::RecordMolecularMaterial(), G4VEmAdjointModel::SelectCSMatrix(), G4EnergyLossForExtrapolator::SetupKinematics(), G4EnergyLossForExtrapolator::TrueStepLength(), and G4EmSaturation::VisibleEnergyDeposition().

◆ GetIonisation()

G4IonisParamMat * G4Material::GetIonisation ( ) const
inlineinherited

Definition at line 222 of file G4Material.hh.

222{return fIonisation;}

References G4Material::fIonisation.

Referenced by G4EmCorrections::Bethe(), G4NistMaterialBuilder::BuildMaterial(), G4DensityEffectCalculator::ComputeDensityCorrection(), G4UrbanAdjointMscModel::ComputeTruePathLengthLimit(), G4GoudsmitSaundersonMscModel::ComputeTruePathLengthLimit(), G4Material::CopyPointersOfBaseMaterial(), DetectorConstruction::DefineMaterials(), ExN03DetectorConstruction::DefineMaterials(), G4EmCorrections::DensityCorrection(), G4AtimaFluctuations::Dispersion(), G4EmSaturation::DumpBirksCoefficients(), G4tgbGeometryDumper::DumpMaterial(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), export_G4Material(), G4ElectronIonPair::FindG4MeanEnergyPerIonPair(), G4tgbMaterialMgr::FindOrBuildG4Material(), G4DensityEffectCalculator::G4DensityEffectCalculator(), G4EmSaturation::InitialiseBirksCoefficient(), G4GoudsmitSaundersonTable::InitSCPCorrection(), G4eDPWAElasticDCS::InitSCPCorrection(), G4GDMLWriteMaterials::MaterialWrite(), G4IonParametrisedLossModel::MinEnergyCut(), G4MuBetheBlochModel::MinEnergyCut(), G4mplIonisationWithDeltaModel::MinEnergyCut(), G4AtimaEnergyLossModel::MinEnergyCut(), G4BetheBlochModel::MinEnergyCut(), G4BraggIonModel::MinEnergyCut(), G4BraggModel::MinEnergyCut(), G4LindhardSorensenIonModel::MinEnergyCut(), G4CoulombScattering::MinPrimaryEnergy(), G4IonFluctuations::RelativisticFactor(), G4IonisParamMat::SetDensityEffectParameters(), G4WentzelOKandVIxSection::SetupKinematic(), G4WentzelVIRelXSection::SetupKinematic(), G4EmCorrections::ShellCorrectionSTD(), and G4EmSaturation::VisibleEnergyDeposition().

◆ GetMassOfMolecule()

G4double G4Material::GetMassOfMolecule ( ) const
inlineinherited

◆ GetMatComponents()

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

◆ GetMaterial() [1/3]

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

Definition at line 686 of file G4Material.cc.

687{
688 // search the material by its name
689 for (size_t j=0; j<theMaterialTable.size(); ++j) {
690 if (theMaterialTable[j]->GetName() == materialName) {
691 return theMaterialTable[j];
692 }
693 }
694
695 // the material does not exist in the table
696 if (warn) {
697 G4cout << "G4Material::GetMaterial() WARNING: The material: "
698 << materialName
699 << " does not exist in the table. Return NULL pointer."
700 << G4endl;
701 }
702 return nullptr;
703}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4MaterialTable theMaterialTable
Definition: G4Material.hh:316
const G4String & GetName() const
Definition: G4Material.hh:173

References G4cout, G4endl, G4Material::GetName(), and G4Material::theMaterialTable.

Referenced by G4DNABrownianTransportation::BuildPhysicsTable(), G4LossTableBuilder::BuildTableForModel(), G4ProductionCutsTable::CheckMaterialInfo(), MyDetectorConstruction::Construct(), QDetectorConstruction::Construct(), demo::ConstructGeom(), Lesson1::ConstructGeom(), run::ConstructGeom(), test::ConstructGeom(), test_voxel::ConstructGeom(), G4EzWorld::CreateWorld(), G4DNAELSEPAElasticModel::CrossSectionPerVolume(), G4DNADummyModel::CrossSectionPerVolume(), G4EmCalculator::FindMaterial(), G4GDMLReadMaterials::GetMaterial(), G4DNAChampionElasticModel::Initialise(), G4DNACPA100ElasticModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4TDNAOneStepThermalizationModel< MODEL >::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNATransformElectronModel::Initialise(), G4DNAUeharaScreenedRutherfordElasticModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4DNABornExcitationModel1::Initialise(), G4DNABornExcitationModel2::Initialise(), G4DNABornIonisationModel1::Initialise(), G4DNABornIonisationModel2::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNAEmfietzoglouIonisationModel::Initialise(), G4DNARelativisticIonisationModel::Initialise(), G4DNADummyModel::Initialise(), G4DNAIonElasticModel::Initialise(), ExN03DetectorConstruction::SetAbsorberMaterial(), ExN03DetectorConstruction::SetGapMaterial(), G4DNAMolecularMaterial::SetMolecularConfiguration(), and G4ProductionCutsTable::UpdateCoupleTable().

◆ GetMaterial() [2/3]

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

Definition at line 707 of file G4Material.cc.

708{
709 // search the material by its name
710 for (size_t j=0; j<theMaterialTable.size(); ++j) {
711 G4Material* mat = theMaterialTable[j];
712 if (1 == mat->GetNumberOfElements() &&
713 z == mat->GetZ() && a == mat->GetA() && dens == mat->GetDensity()) {
714 return mat;
715 }
716 }
717 return nullptr;
718}
G4double GetZ() const
Definition: G4Material.cc:736
size_t GetNumberOfElements() const
Definition: G4Material.hh:182

References G4Material::GetA(), G4Material::GetDensity(), G4Material::GetNumberOfElements(), G4Material::GetZ(), and G4Material::theMaterialTable.

◆ GetMaterial() [3/3]

G4Material * G4Material::GetMaterial ( size_t  nComp,
G4double  dens 
)
staticinherited

Definition at line 722 of file G4Material.cc.

723{
724 // search the material by its name
725 for (size_t j=0; j<theMaterialTable.size(); ++j) {
726 G4Material* mat = theMaterialTable[j];
727 if (nComp == mat->GetNumberOfElements() && dens == mat->GetDensity()) {
728 return mat;
729 }
730 }
731 return nullptr;
732}

References G4Material::GetDensity(), G4Material::GetNumberOfElements(), and G4Material::theMaterialTable.

◆ GetMaterialPropertiesTable()

G4MaterialPropertiesTable * G4Material::GetMaterialPropertiesTable ( ) const
inlineinherited

◆ GetMaterialTable()

G4MaterialTable * G4Material::GetMaterialTable ( )
staticinherited

Definition at line 672 of file G4Material.cc.

673{
674 return &theMaterialTable;
675}

References G4Material::theMaterialTable.

Referenced by G4VCrossSectionHandler::ActiveElements(), G4PixeCrossSectionHandler::ActiveElements(), G4AdjointCSManager::BuildCrossSectionMatrices(), G4PixeCrossSectionHandler::BuildCrossSectionsForMaterials(), G4DNAModelInterface::BuildMaterialMolPerVolTable(), G4DNAModelInterface::BuildMaterialParticleModelTable(), G4NistMaterialBuilder::BuildNistMaterial(), G4ParticleHPThermalScattering::buildPhysicsTable(), G4ParticleHPThermalScatteringData::BuildPhysicsTable(), G4VLEPTSModel::BuildPhysicsTable(), G4Cerenkov::BuildPhysicsTable(), G4Scintillation::BuildPhysicsTable(), G4OpRayleigh::BuildPhysicsTable(), G4OpWLS::BuildPhysicsTable(), G4OpWLS2::BuildPhysicsTable(), G4TablesForExtrapolator::ComputeElectronDEDX(), G4VXTRenergyLoss::ComputeGasPhotoAbsCof(), G4PAIxSection::ComputeLowEnergyCof(), G4StrawTubeXTRadiator::ComputeMediumPhotoAbsCof(), G4TablesForExtrapolator::ComputeMuonDEDX(), G4VXTRenergyLoss::ComputePlatePhotoAbsCof(), G4TablesForExtrapolator::ComputeProtonDEDX(), G4TablesForExtrapolator::ComputeTrasportXS(), DetectorConstruction::DefineMaterials(), ExN03DetectorConstruction::DefineMaterials(), G4EmSaturation::DumpBirksCoefficients(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), export_G4Material(), G4NistMaterialBuilder::FindMaterial(), G4PAIxSection::G4PAIxSection(), G4SandiaTable::G4SandiaTable(), G4DNAMolecularMaterial::GetDensityTableFor(), G4VXTRenergyLoss::GetGasCompton(), G4NistManager::GetMaterial(), getMaterialTable(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VXTRenergyLoss::GetPlateCompton(), G4TablesForExtrapolator::Initialisation(), G4ICRU90StoppingData::Initialise(), G4ASTARStopping::Initialise(), G4PSTARStopping::Initialise(), G4LEPTSElasticModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4eeToTwoGammaModel::Initialise(), G4ICRU73QOModel::Initialise(), G4PEEffectFluoModel::Initialise(), G4HadronXSDataTable::Initialise(), G4EmSaturation::InitialiseG4Saturation(), G4DNAMolecularMaterial::Initialize(), G4DNAMolecularMaterial::InitializeDensity(), G4GoudsmitSaundersonTable::InitMoliereMSCParams(), G4VDNAModel::IsMaterialDefine(), G4eIonisationParameters::LoadData(), G4NistManager::PrintG4Material(), G4ProductionCutsTable::StoreMaterialInfo(), G4GDMLRead::StripNames(), and G4NistManager::~G4NistManager().

◆ GetName()

const G4String & G4Material::GetName ( ) const
inlineinherited

Definition at line 173 of file G4Material.hh.

173{return fName;}

References G4Material::fName.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4VtkSceneHandler::AddSolid(), G4ErrorEnergyLoss::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4ParticleHPInelastic::ApplyYourself(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4EmCorrections::BarkasCorrection(), G4DNAModelInterface::BuildMaterialMolPerVolTable(), G4DNAModelInterface::BuildMaterialParticleModelTable(), G4VLEPTSModel::BuildPhysicsTable(), G4PenelopeIonisationXSHandler::BuildXSTable(), G4PenelopeBremsstrahlungModel::BuildXSTable(), G4ProductionCutsTable::CheckMaterialCutsCoupleInfo(), G4EmCalculator::ComputeCrossSectionPerVolume(), G4DensityEffectCalculator::ComputeDensityCorrection(), G4IonisParamMat::ComputeDensityEffectParameters(), G4TablesForExtrapolator::ComputeElectronDEDX(), G4SandiaTable::ComputeMatSandiaMatrix(), G4SandiaTable::ComputeMatSandiaMatrixPAI(), G4SandiaTable::ComputeMatTable(), G4EmCalculator::ComputeMeanFreePath(), G4TablesForExtrapolator::ComputeMuonDEDX(), G4TablesForExtrapolator::ComputeProtonDEDX(), G4TablesForExtrapolator::ComputeTrasportXS(), ExN03DetectorConstruction::ConstructCalorimeter(), G4tgbVolume::ConstructG4LogVol(), G4tgbVolume::ConstructG4PhysVol(), DetectorConstruction::ConstructVolumes(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G4DNAModelInterface::CrossSectionPerVolume(), G4EmElementSelector::Dump(), G4EmSaturation::DumpBirksCoefficients(), G4ProductionCutsTable::DumpCouples(), G4PenelopeRayleighModel::DumpFormFactorTable(), G4PenelopeRayleighModelMI::DumpFormFactorTable(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), G4HadronicProcess::DumpState(), G4MuonicAtomDecay::DumpState(), G4ExceptionHandler::DumpTrackInfo(), G4EmCorrections::EffectiveChargeCorrection(), export_G4Material(), G4DensityEffectCalculator::FermiDeltaCalculation(), G4EmModelManager::FillDEDXVector(), G4EmModelManager::FillLambdaVector(), G4EmSaturation::FindG4BirksCoefficient(), G4ElectronIonPair::FindG4MeanEnergyPerIonPair(), G4VEmProcess::FindLambdaMax(), G4IonisParamMat::FindMeanExcitationEnergy(), G4tgbMaterialMgr::FindOrBuildG4Material(), G4PenelopeIonisationCrossSection::FindShellIDIndex(), G4ForwardXrayTR::G4ForwardXrayTR(), G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4PenelopeOscillatorManager::GetAtomsPerMolecule(), G4CrossSectionDataStore::GetCrossSection(), G4EmCalculator::GetCrossSectionPerVolume(), G4PenelopeBremsstrahlungModel::GetCrossSectionTableForCouple(), G4EmCalculator::GetCSDARange(), G4EmCalculator::GetDEDX(), G4PenelopeIonisationXSHandler::GetDensityCorrection(), G4VCrossSectionDataSet::GetElementCrossSection(), G4PenelopeRayleighModel::GetFSquared(), G4PenelopeRayleighModelMI::GetFSquared(), G4ESTARStopping::GetIndex(), G4CrossSectionDataStore::GetIsoCrossSection(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4EmCalculator::GetKinEnergy(), G4LatticeManager::GetLattice(), G4Material::GetMaterial(), G4PenelopeOscillatorManager::GetMeanExcitationEnergy(), G4EmCalculator::GetMeanFreePath(), G4PenelopeOscillatorManager::GetNumberOfZAtomsPerMolecule(), G4DNAModelInterface::GetNumMoleculePerVolumeUnitForMaterial(), G4DNAModelInterface::GetNumMolPerVolUnitForComponentInComposite(), G4PenelopeOscillatorManager::GetOscillatorTableCompton(), G4PenelopeOscillatorManager::GetOscillatorTableIonisation(), G4PenelopeOscillatorManager::GetPlasmaEnergySquared(), G4PenelopeRayleighModel::GetPMaxTable(), G4PenelopeRayleighModelMI::GetPMaxTable(), G4EmCalculator::GetRangeFromRestricteDEDX(), G4PenelopeOscillatorManager::GetTotalA(), G4PenelopeOscillatorManager::GetTotalZ(), G4EmCorrections::HighOrderCorrections(), G4ICRU90StoppingData::Initialise(), G4ASTARStopping::Initialise(), G4IonICRU73Data::Initialise(), G4PSTARStopping::Initialise(), G4LEPTSElasticModel::Initialise(), G4PAIModel::Initialise(), G4PAIPhotModel::Initialise(), G4EmSaturation::InitialiseBirksCoefficient(), G4Material::InitializePointers(), G4PenelopeRayleighModel::InitializeSamplingAlgorithm(), G4PenelopeRayleighModelMI::InitializeSamplingAlgorithm(), G4LatticeManager::LoadLattice(), G4GDMLWriteMaterials::MaterialWrite(), G4DNASecondOrderReaction::PostStepDoIt(), G4MicroElecSurface::PostStepDoIt(), G4VXTRenergyLoss::PostStepDoIt(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4Decay::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), ExN03DetectorConstruction::PrintCalorParameters(), G4NistManager::PrintG4Material(), GVFlashShowerParameterisation::PrintMaterial(), G4DNAMolecularMaterial::PrintNotAMolecularMaterial(), DetectorConstruction::PrintParameters(), G4ErrorFreeTrajState::PropagateErrorMSC(), RegisterExtension(), G4ASCIITreeSceneHandler::RequestPrimitives(), RetrieveExtension(), G4VLEPTSModel::SampleEnergyLoss(), G4PenelopeBremsstrahlungFS::SampleGammaEnergy(), G4LivermoreComptonModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4DNAModelInterface::SampleSecondaries(), G4IonisParamMat::SetMeanExcitationEnergy(), G4EmCalculator::SetupMaterial(), G4EnergySplitter::SplitEnergyInVolumes(), G4ProductionCutsTable::StoreMaterialCutsCoupleInfo(), G4GDMLRead::StripNames(), G4ParallelWorldProcess::SwitchMaterial(), G4tgbMaterialMixtureByVolume::TransformToFractionsByWeight(), G4GDMLWriteStructure::TraverseVolumeTree(), and G4EmCalculator::UpdateParticle().

◆ GetNuclearInterLength()

G4double G4Material::GetNuclearInterLength ( ) const
inlineinherited

◆ GetNumberOfElements()

size_t G4Material::GetNumberOfElements ( ) const
inlineinherited

Definition at line 182 of file G4Material.hh.

182{return fNumberOfElements;}

References G4Material::fNumberOfElements.

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSection(), G4FissLib::ApplyYourself(), G4ParticleHPCapture::ApplyYourself(), G4ParticleHPFission::ApplyYourself(), G4ParticleHPInelastic::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4ParticleHPElastic::ApplyYourself(), G4ErrorFreeTrajState::CalculateEffectiveZandA(), G4EmCalculator::CheckMaterial(), G4Nucleus::ChooseParameters(), G4AdjointCSManager::ComputeAdjointCS(), G4CrossSectionDataStore::ComputeCrossSection(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4IonisParamMat::ComputeDensityEffectOnFly(), G4IonisParamMat::ComputeDensityEffectParameters(), G4WentzelVIRelModel::ComputeEffectiveMass(), G4IonisParamMat::ComputeFluctModel(), G4IonisParamMat::ComputeIonParameters(), G4SandiaTable::ComputeMatSandiaMatrix(), G4SandiaTable::ComputeMatSandiaMatrixPAI(), G4SandiaTable::ComputeMatTable(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4IonisParamMat::ComputeMeanParameters(), G4eDPWAElasticDCS::ComputeMParams(), G4WentzelVIModel::ComputeSecondMoment(), G4WentzelVIModel::ComputeTransportXSectionPerVolume(), G4VEmModel::CrossSectionPerVolume(), G4tgbGeometryDumper::DumpMaterial(), G4DensityEffectCalculator::G4DensityEffectCalculator(), G4EmElementSelector::G4EmElementSelector(), G4HadElementSelector::G4HadElementSelector(), G4Material::G4Material(), G4CrossSectionDataStore::GetCrossSection(), G4IonICRU73Data::GetDEDX(), GVFlashShowerParameterisation::GetEffA(), GVFlashShowerParameterisation::GetEffZ(), G4Material::GetMaterial(), G4GSPWACorrections::InitDataMaterial(), G4IonICRU73Data::Initialise(), G4LivermoreGammaConversion5DModel::Initialise(), G4eDPWACoulombScatteringModel::Initialise(), G4WentzelVIModel::Initialise(), G4HadronXSDataTable::Initialise(), G4EmSaturation::InitialiseBirksCoefficient(), G4GSMottCorrection::InitMCDataMaterial(), G4GoudsmitSaundersonTable::InitMoliereMSCParams(), G4GDMLWriteMaterials::MaterialWrite(), G4IonICRU73Data::ReadElementData(), G4VLEPTSModel::ReadParam(), G4WentzelVIModel::SampleScattering(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4CrossSectionDataStore::SampleZandA(), G4GammaConversionToMuons::SelectRandomAtom(), G4VEmModel::SelectRandomAtom(), G4VEmModel::SelectRandomAtomNumber(), G4ElementSelector::SelectZandA(), and G4EmCorrections::SetupKinematics().

◆ GetNumberOfExtensions()

G4int G4ExtendedMaterial::GetNumberOfExtensions ( ) const
inline

Definition at line 119 of file G4ExtendedMaterial.hh.

120 { return G4int(fExtensionMap.size()); }

References fExtensionMap.

◆ GetNumberOfMaterials()

size_t G4Material::GetNumberOfMaterials ( )
staticinherited

◆ GetPressure()

G4double G4Material::GetPressure ( ) const
inlineinherited

◆ GetRadlen()

G4double G4Material::GetRadlen ( ) const
inlineinherited

◆ GetSandiaTable()

G4SandiaTable * G4Material::GetSandiaTable ( ) const
inlineinherited

◆ GetState()

G4State G4Material::GetState ( ) const
inlineinherited

◆ GetTemperature()

G4double G4Material::GetTemperature ( ) const
inlineinherited

Definition at line 178 of file G4Material.hh.

178{return fTemp;}
G4double fTemp
Definition: G4Material.hh:336

References G4Material::fTemp.

Referenced by G4ParticleHPChannelList::ApplyYourself(), G4FissLib::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4LENDModel::ApplyYourself(), G4ParticleHPCapture::ApplyYourself(), G4ParticleHPFission::ApplyYourself(), G4ParticleHPInelastic::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4ParticleHPElastic::ApplyYourself(), G4FissionLibrary::ApplyYourself(), G4ParticleHPCaptureFS::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4ParticleHPFissionFS::ApplyYourself(), G4ParticleHPChannel::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4tgbMaterialMixtureByWeight::BuildG4Material(), G4NistManager::BuildMaterialWithNewDensity(), G4ParticleHPInelasticCompFS::CompositeApply(), G4IonisParamMat::ComputeDensityEffectParameters(), G4DNABrownianTransportation::ComputeStep(), G4NistMaterialBuilder::ConstructNewGasMaterial(), G4tgbGeometryDumper::DumpMaterial(), export_G4Material(), G4DNAElectronHoleRecombination::FindReactant(), G4Nucleus::G4Nucleus(), G4ParticleHPThermalScatteringData::GetCoherentCrossSection(), G4ParticleHPThermalScatteringData::GetCrossSection(), G4ParticleHPThermalScatteringData::GetIncoherentCrossSection(), G4ParticleHPThermalScatteringData::GetInelasticCrossSection(), and G4GDMLWriteMaterials::MaterialWrite().

◆ GetTotNbOfAtomsPerVolume()

G4double G4Material::GetTotNbOfAtomsPerVolume ( ) const
inlineinherited

◆ GetTotNbOfElectPerVolume()

G4double G4Material::GetTotNbOfElectPerVolume ( ) const
inlineinherited

◆ GetVecNbOfAtomsPerVolume()

const G4double * G4Material::GetVecNbOfAtomsPerVolume ( ) const
inlineinherited

◆ GetZ()

G4double G4Material::GetZ ( ) const
inherited

Definition at line 736 of file G4Material.cc.

737{
738 if (fNumberOfElements > 1) {
740 ed << "For material " << fName << " ERROR in GetZ() - Nelm="
741 << fNumberOfElements << " > 1, which is not allowed";
742 G4Exception ("G4Material::GetZ()", "mat036", FatalException,
743 ed, "");
744 }
745 return (*theElementVector)[0]->GetZ();
746}

References FatalException, G4Material::fName, G4Material::fNumberOfElements, G4Exception(), and G4Material::theElementVector.

Referenced by G4tgbGeometryDumper::DumpMaterial(), export_G4Material(), GVFlashShowerParameterisation::GetEffZ(), G4Material::GetMaterial(), and G4GDMLWriteMaterials::MaterialWrite().

◆ InitializePointers()

void G4Material::InitializePointers ( )
privateinherited

Definition at line 247 of file G4Material.cc.

248{
249 fBaseMaterial = nullptr;
250 fMaterialPropertiesTable = nullptr;
251 theElementVector = nullptr;
252 fAtomsVector = nullptr;
253 fMassFractionVector = nullptr;
254 fVecNbOfAtomsPerVolume = nullptr;
255
256 fIonisation = nullptr;
257 fSandiaTable = nullptr;
258
263
265
267
268 fMassFraction = true;
269
270 fChemicalFormula = "";
271
272 // initilized data members
273
274 // Store in the static Table of Materials
276 for(size_t i=0; i<fIndexInTable; ++i) {
277 if(theMaterialTable[i]->GetName() == fName) {
278 G4cout << "G4Material WARNING: duplicate name of material "
279 << fName << G4endl;
280 break;
281 }
282 }
283 theMaterialTable.push_back(this);
284}

References G4Material::fAtomsVector, G4Material::fBaseMaterial, G4Material::fChemicalFormula, G4Material::fDensity, G4Material::fFreeElecDensity, G4Material::fIdxComponent, G4Material::fIndexInTable, G4Material::fIonisation, G4Material::fMassFraction, G4Material::fMassFractionVector, G4Material::fMassOfMolecule, G4Material::fMaterialPropertiesTable, G4Material::fName, G4Material::fNbComponents, G4Material::fNuclInterLen, G4Material::fNumberOfElements, G4Material::fPressure, G4Material::fRadlen, G4Material::fSandiaTable, G4Material::fState, G4Material::fTemp, G4Material::fTotNbOfAtomsPerVolume, G4Material::fTotNbOfElectPerVolume, G4Material::fVecNbOfAtomsPerVolume, G4cout, G4endl, G4Material::GetName(), kStateUndefined, G4Material::theElementVector, and G4Material::theMaterialTable.

Referenced by G4Material::G4Material().

◆ IsExtended()

G4bool G4ExtendedMaterial::IsExtended ( ) const
virtual

Reimplemented from G4Material.

Definition at line 119 of file G4ExtendedMaterial.cc.

120{ return true; }

◆ operator!=()

G4bool G4Material::operator!= ( const G4Material ) const
deleteinherited

◆ operator==()

G4bool G4Material::operator== ( const G4Material ) const
deleteinherited

◆ Print()

void G4ExtendedMaterial::Print ( std::ostream &  flux) const

Definition at line 124 of file G4ExtendedMaterial.cc.

125{
126 flux << "\n Registered material extensions :\n";
127 auto iter = fExtensionMap.begin();
128 for(;iter!=fExtensionMap.end();iter++)
129 { flux << " " << iter->first << "\n"; }
130}

References fExtensionMap.

◆ RegisterExtension()

void G4ExtendedMaterial::RegisterExtension ( std::unique_ptr< G4VMaterialExtension extension)

Definition at line 84 of file G4ExtendedMaterial.cc.

85{
86 auto iter = fExtensionMap.find(extension->GetName());
87 if(iter!=fExtensionMap.end())
88 {
90 msg << "G4ExtendedMaterial <"<<GetName()<<"> already has extension for "
91 << extension->GetName()
92 << ". Extension is replaced.";
93 G4Exception("G4ExtendedMaterial::RegisterExtension(...)","MatExt001",JustWarning,msg);
94 }
95 fExtensionMap.insert(std::make_pair(extension->GetName(),std::move(extension)));
96}

References fExtensionMap, G4Exception(), G4Material::GetName(), and JustWarning.

◆ RetrieveExtension()

G4VMaterialExtension * G4ExtendedMaterial::RetrieveExtension ( const G4String name)

Definition at line 102 of file G4ExtendedMaterial.cc.

103{
104 const auto iter = fExtensionMap.find(name);
105 if(iter!=fExtensionMap.end())
106 { return iter->second.get(); }
107 else
108 {
110 msg << "G4ExtendedMAterial <"<<GetName()<<"> cannot find extension for "
111 << name;
112 G4Exception("G4ExtendedMaterial::RetreiveExtension(...)","MatExt002",JustWarning,msg);
113 return nullptr;
114 }
115}

References fExtensionMap, G4Exception(), G4Material::GetName(), JustWarning, and G4InuclParticleNames::name().

Referenced by G4PenelopeRayleighModelMI::CrossSectionPerVolume(), G4Channeling::GetMatData(), and G4PenelopeRayleighModelMI::GetMIData().

◆ SetChemicalFormula()

void G4Material::SetChemicalFormula ( const G4String chF)
inherited

Definition at line 632 of file G4Material.cc.

633{
634#ifdef G4MULTITHREADED
635 G4MUTEXLOCK(&materialMutex);
636#endif
637 fChemicalFormula = chF;
638#ifdef G4MULTITHREADED
639 G4MUTEXUNLOCK(&materialMutex);
640#endif
641}

References G4Material::fChemicalFormula, G4MUTEXLOCK, and G4MUTEXUNLOCK.

Referenced by G4NistMaterialBuilder::BuildMaterial(), and export_G4Material().

◆ SetFreeElectronDensity()

void G4Material::SetFreeElectronDensity ( G4double  val)
inherited

Definition at line 645 of file G4Material.cc.

646{
647#ifdef G4MULTITHREADED
648 G4MUTEXLOCK(&materialMutex);
649#endif
650 if(val >= 0.) { fFreeElecDensity = val; }
651#ifdef G4MULTITHREADED
652 G4MUTEXUNLOCK(&materialMutex);
653#endif
654}

References G4Material::fFreeElecDensity, G4MUTEXLOCK, and G4MUTEXUNLOCK.

◆ SetMaterialPropertiesTable()

void G4Material::SetMaterialPropertiesTable ( G4MaterialPropertiesTable anMPT)
inherited

Definition at line 839 of file G4Material.cc.

840{
841 if(nullptr != anMPT && fMaterialPropertiesTable != anMPT) {
842#ifdef G4MULTITHREADED
843 G4MUTEXLOCK(&materialMutex);
844 if(fMaterialPropertiesTable != anMPT) {
845#endif
848#ifdef G4MULTITHREADED
849 }
850 G4MUTEXUNLOCK(&materialMutex);
851#endif
852 }
853}

References G4Material::fMaterialPropertiesTable, G4MUTEXLOCK, and G4MUTEXUNLOCK.

Referenced by export_G4Material().

◆ SetName()

void G4Material::SetName ( const G4String name)
inlineinherited

Definition at line 285 of file G4Material.hh.

285{fName=name;}

References G4Material::fName, and G4InuclParticleNames::name().

Referenced by export_G4Material(), and G4GDMLRead::StripNames().

Field Documentation

◆ fAtoms

std::vector<G4int>* G4Material::fAtoms
privateinherited

Definition at line 355 of file G4Material.hh.

Referenced by G4Material::AddElementByNumberOfAtoms().

◆ fAtomsVector

G4int* G4Material::fAtomsVector
privateinherited

◆ fBaseMaterial

const G4Material* G4Material::fBaseMaterial
privateinherited

◆ fChemicalFormula

G4String G4Material::fChemicalFormula
privateinherited

◆ fDensity

G4double G4Material::fDensity
privateinherited

◆ fElm

std::vector<const G4Element*>* G4Material::fElm
privateinherited

◆ fElmFrac

std::vector<G4double>* G4Material::fElmFrac
privateinherited

◆ fExtensionMap

G4MaterialExtensionMap G4ExtendedMaterial::fExtensionMap
private

◆ fFreeElecDensity

G4double G4Material::fFreeElecDensity
privateinherited

◆ fIdxComponent

G4int G4Material::fIdxComponent
privateinherited

◆ fIndexInTable

size_t G4Material::fIndexInTable
privateinherited

◆ fIonisation

G4IonisParamMat* G4Material::fIonisation
privateinherited

◆ fMassFraction

G4bool G4Material::fMassFraction
privateinherited

◆ fMassFractionVector

G4double* G4Material::fMassFractionVector
privateinherited

◆ fMassOfMolecule

G4double G4Material::fMassOfMolecule
privateinherited

◆ fMatComponents

std::map<G4Material*, G4double> G4Material::fMatComponents
privateinherited

Definition at line 360 of file G4Material.hh.

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

◆ fMaterialPropertiesTable

G4MaterialPropertiesTable* G4Material::fMaterialPropertiesTable
privateinherited

◆ fName

G4String G4Material::fName
privateinherited

◆ fNbComponents

G4int G4Material::fNbComponents
privateinherited

◆ fNuclInterLen

G4double G4Material::fNuclInterLen
privateinherited

◆ fNumberOfElements

G4int G4Material::fNumberOfElements
privateinherited

◆ fPressure

G4double G4Material::fPressure
privateinherited

◆ fRadlen

G4double G4Material::fRadlen
privateinherited

◆ fSandiaTable

G4SandiaTable* G4Material::fSandiaTable
privateinherited

◆ fState

G4State G4Material::fState
privateinherited

◆ fTemp

G4double G4Material::fTemp
privateinherited

◆ fTotNbOfAtomsPerVolume

G4double G4Material::fTotNbOfAtomsPerVolume
privateinherited

◆ fTotNbOfElectPerVolume

G4double G4Material::fTotNbOfElectPerVolume
privateinherited

◆ fVecNbOfAtomsPerVolume

G4double* G4Material::fVecNbOfAtomsPerVolume
privateinherited

◆ theElementVector

G4ElementVector* G4Material::theElementVector
privateinherited

◆ theMaterialTable

G4MaterialTable G4Material::theMaterialTable
staticprivateinherited

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