Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4Material Class Reference

#include <G4Material.hh>

Inheritance diagram for G4Material:

Public Member Functions

 G4Material (const G4String &name, G4double z, G4double a, G4double density, G4State state=kStateUndefined, G4double temp=NTP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
 G4Material (const G4String &name, G4double density, G4int nComponents, G4State state=kStateUndefined, G4double temp=NTP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
 G4Material (const G4String &name, G4double density, const G4Material *baseMaterial, G4State state=kStateUndefined, G4double temp=NTP_Temperature, G4double pressure=CLHEP::STP_Pressure)
 
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
 
const std::map< G4Material
*, G4double > & 
GetMatComponents () const
 
G4double GetMassOfMolecule () const
 
G4double GetZ () const
 
G4double GetA () const
 
void SetMaterialPropertiesTable (G4MaterialPropertiesTable *anMPT)
 
G4MaterialPropertiesTableGetMaterialPropertiesTable () const
 
size_t GetIndex () const
 
 G4Material (__void__ &)
 
void SetName (const G4String &name)
 
virtual G4bool IsExtended () const
 

Static Public Member Functions

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

Friends

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

Detailed Description

Definition at line 120 of file G4Material.hh.

Constructor & Destructor Documentation

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

Definition at line 90 of file G4Material.cc.

93  : fName(name)
94 {
95  InitializePointers();
96 
97  if (density < universe_mean_density)
98  {
99  G4cout << " G4Material WARNING:"
100  << " define a material with density=0 is not allowed. \n"
101  << " The material " << name << " will be constructed with the"
102  << " default minimal density: " << universe_mean_density/(g/cm3)
103  << "g/cm3" << G4endl;
104  density = universe_mean_density;
105  }
106 
107  fDensity = density;
108  fState = state;
109  fTemp = temp;
110  fPressure = pressure;
111 
112  // Initialize theElementVector allocating one
113  // element corresponding to this material
114  maxNbComponents = fNumberOfComponents = fNumberOfElements = 1;
115  fArrayLength = maxNbComponents;
116  theElementVector = new G4ElementVector();
117 
118  const std::vector<G4String> elmnames =
120  G4String enam, snam;
121  G4int iz = G4lrint(z);
122  if(iz < (G4int)elmnames.size()) {
123  snam = elmnames[iz];
124  enam = snam;
125  } else {
126  enam = "ELM_" + name;
127  snam = name;
128  }
129  theElementVector->push_back(new G4Element(enam, snam, z, a));
130 
131  fMassFractionVector = new G4double[1];
132  fMassFractionVector[0] = 1. ;
133  fMassOfMolecule = a/CLHEP::Avogadro;
134 
135  if (fState == kStateUndefined)
136  {
137  if (fDensity > kGasThreshold) { fState = kStateSolid; }
138  else { fState = kStateGas; }
139  }
140 
141  ComputeDerivedQuantities();
142 }
const XML_Char * name
Definition: expat.h:151
std::vector< G4Element * > G4ElementVector
int kGasThreshold
Definition: hepunit.py:304
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
int universe_mean_density
Definition: hepunit.py:307
const std::vector< G4String > & GetNistElementNames() const
int G4int
Definition: G4Types.hh:78
static constexpr double Avogadro
static G4NistManager * Instance()
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
static constexpr double cm3
Definition: G4SIunits.hh:121
int G4lrint(double ad)
Definition: templates.hh:163
tuple z
Definition: test.py:28
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Definition at line 149 of file G4Material.cc.

152  : fName(name)
153 {
154  InitializePointers();
155 
156  if (density < universe_mean_density)
157  {
158  G4cout << "--- Warning from G4Material::G4Material()"
159  << " define a material with density=0 is not allowed. \n"
160  << " The material " << name << " will be constructed with the"
161  << " default minimal density: " << universe_mean_density/(g/cm3)
162  << "g/cm3" << G4endl;
163  density = universe_mean_density;
164  }
165 
166  fDensity = density;
167  fState = state;
168  fTemp = temp;
169  fPressure = pressure;
170 
171  maxNbComponents = nComponents;
172  fArrayLength = maxNbComponents;
173  fNumberOfComponents = fNumberOfElements = 0;
174  theElementVector = new G4ElementVector();
175  theElementVector->reserve(maxNbComponents);
176 
177  if (fState == kStateUndefined)
178  {
179  if (fDensity > kGasThreshold) { fState = kStateSolid; }
180  else { fState = kStateGas; }
181  }
182 }
std::vector< G4Element * > G4ElementVector
int kGasThreshold
Definition: hepunit.py:304
int universe_mean_density
Definition: hepunit.py:307
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
static constexpr double cm3
Definition: G4SIunits.hh:121
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

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

Definition at line 188 of file G4Material.cc.

191  : fName(name)
192 {
193  InitializePointers();
194 
195  if (density < universe_mean_density)
196  {
197  G4cout << "--- Warning from G4Material::G4Material()"
198  << " define a material with density=0 is not allowed. \n"
199  << " The material " << name << " will be constructed with the"
200  << " default minimal density: " << universe_mean_density/(g/cm3)
201  << "g/cm3" << G4endl;
202  density = universe_mean_density;
203  }
204 
205  fDensity = density;
206  fState = state;
207  fTemp = temp;
208  fPressure = pressure;
209 
210  fBaseMaterial = bmat;
211  fChemicalFormula = fBaseMaterial->GetChemicalFormula();
212  fMassOfMolecule = fBaseMaterial->GetMassOfMolecule();
213 
214  fNumberOfElements = fBaseMaterial->GetNumberOfElements();
215  maxNbComponents = fNumberOfElements;
216  fNumberOfComponents = fNumberOfElements;
217 
218  CopyPointersOfBaseMaterial();
219 }
int universe_mean_density
Definition: hepunit.py:307
const G4String & GetChemicalFormula() const
Definition: G4Material.hh:179
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
static constexpr double cm3
Definition: G4SIunits.hh:121
G4double GetMassOfMolecule() const
Definition: G4Material.hh:242
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:186

Here is the call graph for this function:

G4Material::~G4Material ( )
virtual

Definition at line 234 of file G4Material.cc.

235 {
236  // G4cout << "### Destruction of material " << fName << " started" <<G4endl;
237  if(fBaseMaterial == nullptr) {
238  delete theElementVector;
239  delete fSandiaTable;
240  //delete fMaterialPropertiesTable;
241  if(fMassFractionVector) { delete [] fMassFractionVector; }
242  if(fAtomsVector) { delete [] fAtomsVector; }
243  }
244  delete fIonisation;
245  if(VecNbOfAtomsPerVolume) { delete [] VecNbOfAtomsPerVolume; }
246 
247  // Remove this material from theMaterialTable.
248  //
249  theMaterialTable[fIndexInTable] = nullptr;
250 }
G4Material::G4Material ( __void__ &  )

Definition at line 226 of file G4Material.cc.

227  : fName("")
228 {
229  InitializePointers();
230 }

Member Function Documentation

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

Definition at line 362 of file G4Material.cc.

363 {
364  // initialization
365  if ( fNumberOfElements == 0 ) {
366  fAtomsVector = new G4int [fArrayLength];
367  fMassFractionVector = new G4double[fArrayLength];
368  }
369 
370  // filling ...
371  if ( fNumberOfElements < maxNbComponents ) {
372  theElementVector->push_back(element);
373  fAtomsVector[fNumberOfElements] = nAtoms;
374  fNumberOfComponents = ++fNumberOfElements;
375  } else {
376  G4cout << "G4Material::AddElement ERROR for " << fName << " nElement= "
377  << fNumberOfElements << G4endl;
378  G4Exception ("G4Material::AddElement()", "mat031", FatalException,
379  "Attempt to add more than the declared number of elements.");
380  }
381  // filled.
382  if ( fNumberOfElements == maxNbComponents ) {
383  // compute proportion by mass
384  G4int i=0;
385  G4double Amol = 0.;
386  for (i=0; i<fNumberOfElements; ++i) {
387  G4double w = fAtomsVector[i]*(*theElementVector)[i]->GetA();
388  Amol += w;
389  fMassFractionVector[i] = w;
390  }
391  for (i=0; i<fNumberOfElements; ++i) {
392  fMassFractionVector[i] /= Amol;
393  }
394 
395  fMassOfMolecule = Amol/Avogadro;
396  ComputeDerivedQuantities();
397  }
398 }
int G4int
Definition: G4Types.hh:78
float Avogadro
Definition: hepunit.py:253
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 404 of file G4Material.cc.

405 {
406  if(fraction < 0.0 || fraction > 1.0) {
407  G4cout << "G4Material::AddElement ERROR for " << fName << " and "
408  << element->GetName() << " mass fraction= " << fraction
409  << " is wrong " << G4endl;
410  G4Exception ("G4Material::AddElement()", "mat032", FatalException,
411  "Attempt to add element with wrong mass fraction");
412  }
413  // initialization
414  if (fNumberOfComponents == 0) {
415  fMassFractionVector = new G4double[fArrayLength];
416  fAtomsVector = new G4int [fArrayLength];
417  }
418  // filling ...
419  if (fNumberOfComponents < maxNbComponents) {
420  G4int el = 0;
421  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
422  while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) { ++el; }
423  if (el<fNumberOfElements) fMassFractionVector[el] += fraction;
424  else {
425  theElementVector->push_back(element);
426  fMassFractionVector[el] = fraction;
427  ++fNumberOfElements;
428  }
429  ++fNumberOfComponents;
430  } else {
431  G4cout << "G4Material::AddElement ERROR for " << fName << " nElement= "
432  << fNumberOfElements << G4endl;
433  G4Exception ("G4Material::AddElement()", "mat033", FatalException,
434  "Attempt to add more than the declared number of elements.");
435  }
436 
437  // filled.
438  if (fNumberOfComponents == maxNbComponents) {
439 
440  G4int i=0;
441  G4double Zmol(0.), Amol(0.);
442  // check sum of weights -- OK?
443  G4double wtSum(0.0);
444  for (i=0; i<fNumberOfElements; ++i) {
445  wtSum += fMassFractionVector[i];
446  Zmol += fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
447  Amol += fMassFractionVector[i]*(*theElementVector)[i]->GetA();
448  }
449  if (std::fabs(1.-wtSum) > perThousand) {
450  G4cerr << "WARNING !! for " << fName << " sum of fractional masses "
451  << wtSum << " is not 1 - results may be wrong"
452  << G4endl;
453  }
454  for (i=0; i<fNumberOfElements; ++i) {
455  fAtomsVector[i] =
456  G4lrint(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA());
457  }
458 
459  ComputeDerivedQuantities();
460  }
461 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetA() const
Definition: G4Material.cc:636
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
static constexpr double perThousand
Definition: G4SIunits.hh:333
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
G4GLOB_DLL std::ostream G4cerr

Here is the call graph for this function:

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

store massFraction of material component

Definition at line 467 of file G4Material.cc.

468 {
469  if(fraction < 0.0 || fraction > 1.0) {
470  G4cout << "G4Material::AddMaterial ERROR for " << fName << " and "
471  << material->GetName() << " mass fraction= " << fraction
472  << " is wrong ";
473  G4Exception ("G4Material::AddMaterial()", "mat034", FatalException,
474  "Attempt to add material with wrong mass fraction");
475  }
476  // initialization
477  if (fNumberOfComponents == 0) {
478  fMassFractionVector = new G4double[fArrayLength];
479  fAtomsVector = new G4int [fArrayLength];
480  }
481 
482  G4int nelm = material->GetNumberOfElements();
483 
484  // arrays should be extended
485  if(nelm > 1) {
486  G4int nold = fArrayLength;
487  fArrayLength += nelm - 1;
488  G4double* v1 = new G4double[fArrayLength];
489  G4int* i1 = new G4int[fArrayLength];
490  for(G4int i=0; i<nold; ++i) {
491  v1[i] = fMassFractionVector[i];
492  i1[i] = fAtomsVector[i];
493  }
494  delete [] fAtomsVector;
495  delete [] fMassFractionVector;
496  fMassFractionVector = v1;
497  fAtomsVector = i1;
498  }
499 
500  // filling ...
501  if (fNumberOfComponents < maxNbComponents) {
502  for (G4int elm=0; elm<nelm; ++elm)
503  {
504  G4Element* element = (*(material->GetElementVector()))[elm];
505  G4int el = 0;
506  // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
507  while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) el++;
508  if (el < fNumberOfElements) fMassFractionVector[el] += fraction
509  *(material->GetFractionVector())[elm];
510  else {
511  theElementVector->push_back(element);
512  fMassFractionVector[el] = fraction
513  *(material->GetFractionVector())[elm];
514  ++fNumberOfElements;
515  }
516  }
517  ++fNumberOfComponents;
519  fMatComponents[material] = fraction;
520 
521  } else {
522  G4cout << "G4Material::AddMaterial ERROR for " << fName << " nElement= "
523  << fNumberOfElements << G4endl;
524  G4Exception ("G4Material::AddMaterial()", "mat035", FatalException,
525  "Attempt to add more than the declared number of components.");
526  }
527 
528  // filled.
529  if (fNumberOfComponents == maxNbComponents) {
530  G4int i=0;
531  G4double Zmol(0.), Amol(0.);
532  // check sum of weights -- OK?
533  G4double wtSum(0.0);
534  for (i=0; i<fNumberOfElements; ++i) {
535  wtSum += fMassFractionVector[i];
536  Zmol += fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
537  Amol += fMassFractionVector[i]*(*theElementVector)[i]->GetA();
538  }
539  if (std::fabs(1.-wtSum) > perThousand) {
540  G4cout << "G4Material::AddMaterial WARNING !! for " << fName
541  << " sum of fractional masses "
542  << wtSum << " is not 1 - results may be wrong"
543  << G4endl;
544  }
545  for (i=0; i<fNumberOfElements; ++i) {
546  fAtomsVector[i] =
547  G4lrint(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA());
548  }
549 
550  ComputeDerivedQuantities();
551  }
552 }
const G4String & GetName() const
Definition: G4Material.hh:178
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetA() const
Definition: G4Material.cc:636
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
static constexpr double perThousand
Definition: G4SIunits.hh:333
double G4double
Definition: G4Types.hh:76
const G4double * GetFractionVector() const
Definition: G4Material.hh:194

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4Material::GetA ( ) const

Definition at line 636 of file G4Material.cc.

637 {
638  if (fNumberOfElements > 1) {
639  G4cout << "G4Material ERROR in GetA. The material: " << fName
640  << " is a mixture.";
641  G4Exception ("G4Material::GetA()", "mat037", FatalException,
642  "the Atomic mass is not well defined." );
643  }
644  return (*theElementVector)[0]->GetA();
645 }
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

const G4double* G4Material::GetAtomicNumDensityVector ( ) const
inline

Definition at line 216 of file G4Material.hh.

216 {return VecNbOfAtomsPerVolume;}

Here is the caller graph for this function:

const G4int* G4Material::GetAtomsVector ( ) const
inline

Definition at line 198 of file G4Material.hh.

198 {return fAtomsVector;}

Here is the caller graph for this function:

const G4Material* G4Material::GetBaseMaterial ( ) const
inline

Definition at line 233 of file G4Material.hh.

233 {return fBaseMaterial;}

Here is the caller graph for this function:

const G4String& G4Material::GetChemicalFormula ( ) const
inline

Definition at line 179 of file G4Material.hh.

179 {return fChemicalFormula;}

Here is the caller graph for this function:

G4double G4Material::GetDensity ( ) const
inline

Definition at line 180 of file G4Material.hh.

180 {return fDensity;}
G4double G4Material::GetElectronDensity ( ) const
inline

Definition at line 217 of file G4Material.hh.

217 {return TotNbOfElectPerVolume;}

Here is the caller graph for this function:

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

Definition at line 202 of file G4Material.hh.

202 {return (*theElementVector)[iel];}

Here is the caller graph for this function:

const G4ElementVector* G4Material::GetElementVector ( ) const
inline

Definition at line 190 of file G4Material.hh.

190 {return theElementVector;}
const G4double* G4Material::GetFractionVector ( ) const
inline

Definition at line 194 of file G4Material.hh.

194 {return fMassFractionVector;}

Here is the caller graph for this function:

size_t G4Material::GetIndex ( ) const
inline

Definition at line 262 of file G4Material.hh.

262 {return fIndexInTable;}
G4IonisParamMat* G4Material::GetIonisation ( ) const
inline

Definition at line 226 of file G4Material.hh.

226 {return fIonisation;}

Here is the caller graph for this function:

G4double G4Material::GetMassOfMolecule ( ) const
inline

Definition at line 242 of file G4Material.hh.

242 {return fMassOfMolecule;}

Here is the caller graph for this function:

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

Definition at line 237 of file G4Material.hh.

238  {return fMatComponents;}

Here is the caller graph for this function:

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

Definition at line 602 of file G4Material.cc.

603 {
604  // search the material by its name
605  for (size_t J=0 ; J<theMaterialTable.size() ; ++J)
606  {
607  if (theMaterialTable[J]->GetName() == materialName)
608  { return theMaterialTable[J]; }
609  }
610 
611  // the material does not exist in the table
612  if (warning) {
613  G4cout << "G4Material::GetMaterial() WARNING: The material: "
614  << materialName
615  << " does not exist in the table. Return NULL pointer."
616  << G4endl;
617  }
618  return nullptr;
619 }
const G4String & GetName() const
Definition: G4Material.hh:178
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4MaterialPropertiesTable* G4Material::GetMaterialPropertiesTable ( ) const
inline

Definition at line 252 of file G4Material.hh.

253  {return fMaterialPropertiesTable;}

Here is the caller graph for this function:

G4MaterialTable * G4Material::GetMaterialTable ( )
static

Definition at line 587 of file G4Material.cc.

588 {
589  return &theMaterialTable;
590 }
const G4String& G4Material::GetName ( void  ) const
inline

Definition at line 178 of file G4Material.hh.

178 {return fName;}
G4double G4Material::GetNuclearInterLength ( ) const
inline

Definition at line 223 of file G4Material.hh.

223 {return fNuclInterLen;}

Here is the caller graph for this function:

size_t G4Material::GetNumberOfElements ( ) const
inline

Definition at line 186 of file G4Material.hh.

186 {return fNumberOfElements;}
size_t G4Material::GetNumberOfMaterials ( )
static

Definition at line 594 of file G4Material.cc.

595 {
596  return theMaterialTable.size();
597 }

Here is the caller graph for this function:

G4double G4Material::GetPressure ( ) const
inline

Definition at line 183 of file G4Material.hh.

183 {return fPressure;}

Here is the caller graph for this function:

G4double G4Material::GetRadlen ( ) const
inline

Definition at line 220 of file G4Material.hh.

220 {return fRadlen;}

Here is the caller graph for this function:

G4SandiaTable* G4Material::GetSandiaTable ( ) const
inline

Definition at line 229 of file G4Material.hh.

229 {return fSandiaTable;}

Here is the caller graph for this function:

G4State G4Material::GetState ( ) const
inline

Definition at line 181 of file G4Material.hh.

181 {return fState;}

Here is the caller graph for this function:

G4double G4Material::GetTemperature ( ) const
inline

Definition at line 182 of file G4Material.hh.

182 {return fTemp;}

Here is the caller graph for this function:

G4double G4Material::GetTotNbOfAtomsPerVolume ( ) const
inline

Definition at line 209 of file G4Material.hh.

209 {return TotNbOfAtomsPerVolume;}

Here is the caller graph for this function:

G4double G4Material::GetTotNbOfElectPerVolume ( ) const
inline

Definition at line 212 of file G4Material.hh.

212 {return TotNbOfElectPerVolume;}

Here is the caller graph for this function:

const G4double* G4Material::GetVecNbOfAtomsPerVolume ( ) const
inline

Definition at line 206 of file G4Material.hh.

206 {return VecNbOfAtomsPerVolume;}

Here is the caller graph for this function:

G4double G4Material::GetZ ( ) const

Definition at line 623 of file G4Material.cc.

624 {
625  if (fNumberOfElements > 1) {
626  G4cout << "G4Material ERROR in GetZ. The material: " << fName
627  << " is a mixture.";
628  G4Exception ("G4Material::GetZ()", "mat036", FatalException,
629  "the Atomic number is not well defined." );
630  }
631  return (*theElementVector)[0]->GetZ();
632 }
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4Material::IsExtended ( ) const
virtual

Reimplemented in G4ExtendedMaterial.

Definition at line 724 of file G4Material.cc.

725 { return false; }

Here is the caller graph for this function:

void G4Material::SetChemicalFormula ( const G4String chF)
inline

Definition at line 173 of file G4Material.hh.

173 {fChemicalFormula=chF;}

Here is the caller graph for this function:

void G4Material::SetMaterialPropertiesTable ( G4MaterialPropertiesTable anMPT)
inline

Definition at line 249 of file G4Material.hh.

250  {fMaterialPropertiesTable = anMPT;}

Here is the caller graph for this function:

void G4Material::SetName ( const G4String name)
inline

Definition at line 279 of file G4Material.hh.

279 {fName=name;}
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

Friends And Related Function Documentation

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

Definition at line 650 of file G4Material.cc.

651 {
652  std::ios::fmtflags mode = flux.flags();
653  flux.setf(std::ios::fixed,std::ios::floatfield);
654  G4long prec = flux.precision(3);
655 
656  flux
657  << " Material: " << std::setw(8) << material->fName
658  << " " << material->fChemicalFormula << " "
659  << " density: " << std::setw(6) << std::setprecision(3)
660  << G4BestUnit(material->fDensity,"Volumic Mass")
661  << " RadL: " << std::setw(7) << std::setprecision(3)
662  << G4BestUnit(material->fRadlen,"Length")
663  << " Nucl.Int.Length: " << std::setw(7) << std::setprecision(3)
664  << G4BestUnit(material->fNuclInterLen,"Length")
665  << "\n" << std::setw(30)
666  << " Imean: " << std::setw(7) << std::setprecision(3)
668  "Energy");
669 
670  if(material->fState == kStateGas) {
671  flux
672  << " temperature: " << std::setw(6) << std::setprecision(2)
673  << (material->fTemp)/kelvin << " K"
674  << " pressure: " << std::setw(6) << std::setprecision(2)
675  << (material->fPressure)/atmosphere << " atm";
676  }
677  flux << "\n";
678 
679  for (G4int i=0; i<material->fNumberOfElements; i++) {
680  flux
681  << "\n ---> " << (*(material->theElementVector))[i]
682  << "\n ElmMassFraction: "
683  << std::setw(6)<< std::setprecision(2)
684  << (material->fMassFractionVector[i])/perCent << " %"
685  << " ElmAbundance " << std::setw(6)<< std::setprecision(2)
686  << 100*(material->VecNbOfAtomsPerVolume[i])
687  /(material->TotNbOfAtomsPerVolume)
688  << " % \n";
689  }
690  flux.precision(prec);
691  flux.setf(mode,std::ios::floatfield);
692 
693  if(material->IsExtended())
694  { static_cast<const G4ExtendedMaterial*>(material)->Print(flux); }
695 
696  return flux;
697 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
static constexpr double perCent
Definition: G4SIunits.hh:332
long G4long
Definition: G4Types.hh:80
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
virtual G4bool IsExtended() const
Definition: G4Material.cc:724
static const double prec
Definition: RanecuEngine.cc:58
void Print(const std::vector< T > &data)
Definition: DicomRun.hh:109
static constexpr double kelvin
Definition: G4SIunits.hh:281
G4double GetMeanExcitationEnergy() const
static constexpr double atmosphere
Definition: G4SIunits.hh:237
std::ostream& operator<< ( std::ostream &  flux,
const G4Material material 
)
friend

Definition at line 701 of file G4Material.cc.

702 {
703  flux << &material;
704  return flux;
705 }
string material
Definition: eplot.py:19
std::ostream& operator<< ( std::ostream &  flux,
G4MaterialTable  MaterialTable 
)
friend

Definition at line 709 of file G4Material.cc.

710 {
711  //Dump info for all known materials
712  flux << "\n***** Table : Nb of materials = " << MaterialTable.size()
713  << " *****\n" << G4endl;
714 
715  for (size_t i=0; i<MaterialTable.size(); ++i) {
716  flux << MaterialTable[i] << G4endl << G4endl;
717  }
718 
719  return flux;
720 }
#define G4endl
Definition: G4ios.hh:61

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