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::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.

References G4cout, G4endl, kStateGas, kStateSolid, and kStateUndefined.

00091   : fName(name)                
00092 {
00093   InitializePointers();
00094     
00095   if (density < universe_mean_density)
00096     { 
00097       G4cout << "--- Warning from G4Material::G4Material()"
00098              << " define a material with density=0 is not allowed. \n"
00099              << " The material " << name << " will be constructed with the"
00100              << " default minimal density: " << universe_mean_density/(g/cm3) 
00101              << "g/cm3" << G4endl;
00102       density = universe_mean_density;
00103     } 
00104 
00105   fDensity  = density;
00106   fState    = state;
00107   fTemp     = temp;
00108   fPressure = pressure;
00109 
00110   // Initialize theElementVector allocating one
00111   // element corresponding to this material
00112   maxNbComponents        = fNumberOfComponents = fNumberOfElements = 1;
00113   fArrayLength           = maxNbComponents;
00114   fImplicitElement       = true;
00115   theElementVector       = new G4ElementVector();
00116   theElementVector->push_back( new G4Element(name, " ", z, a));  
00117   fMassFractionVector    = new G4double[1];
00118   fMassFractionVector[0] = 1. ;
00119   fMassOfMolecule        = a/Avogadro;
00120   
00121   (*theElementVector)[0] -> increaseCountUse();
00122   
00123   if (fState == kStateUndefined)
00124     {
00125       if (fDensity > kGasThreshold) { fState = kStateSolid; }
00126       else                          { fState = kStateGas; }
00127     }
00128 
00129   ComputeDerivedQuantities();
00130 }

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.

References G4cout, G4endl, kStateGas, kStateSolid, and kStateUndefined.

00140   : fName(name)                
00141 {
00142   InitializePointers();
00143     
00144   if (density < universe_mean_density)
00145     {
00146       G4cout << "--- Warning from G4Material::G4Material()"
00147              << " define a material with density=0 is not allowed. \n"
00148              << " The material " << name << " will be constructed with the"
00149              << " default minimal density: " << universe_mean_density/(g/cm3) 
00150              << "g/cm3" << G4endl;
00151       density = universe_mean_density;
00152     }
00153         
00154   fDensity  = density;
00155   fState    = state;
00156   fTemp     = temp;
00157   fPressure = pressure;
00158     
00159   maxNbComponents     = nComponents;
00160   fArrayLength        = maxNbComponents;
00161   fNumberOfComponents = fNumberOfElements = 0;
00162   theElementVector    = new G4ElementVector();
00163   theElementVector->reserve(maxNbComponents);  
00164     
00165   if (fState == kStateUndefined) 
00166     {
00167       if (fDensity > kGasThreshold) { fState = kStateSolid; }
00168       else                          { fState = kStateGas; }
00169     }
00170 }

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.

References G4cout, G4endl, GetChemicalFormula(), GetMassOfMolecule(), GetMaterialPropertiesTable(), and GetNumberOfElements().

00179   : fName(name)                
00180 {
00181   InitializePointers();
00182     
00183   if (density < universe_mean_density)
00184     {
00185       G4cout << "--- Warning from G4Material::G4Material()"
00186              << " define a material with density=0 is not allowed. \n"
00187              << " The material " << name << " will be constructed with the"
00188              << " default minimal density: " << universe_mean_density/(g/cm3) 
00189              << "g/cm3" << G4endl;
00190       density = universe_mean_density;
00191     }
00192 
00193   fDensity  = density;
00194   fState    = state;
00195   fTemp     = temp;
00196   fPressure = pressure;
00197 
00198   fBaseMaterial = bmat;
00199   fChemicalFormula = fBaseMaterial->GetChemicalFormula();
00200   fMassOfMolecule  = fBaseMaterial->GetMassOfMolecule();
00201 
00202   fNumberOfElements = fBaseMaterial->GetNumberOfElements();     
00203   maxNbComponents = fNumberOfElements;
00204   fNumberOfComponents = fNumberOfElements;
00205 
00206   fMaterialPropertiesTable = fBaseMaterial->GetMaterialPropertiesTable();
00207 
00208   CopyPointersOfBaseMaterial();
00209 }

G4Material::~G4Material (  )  [virtual]

Definition at line 227 of file G4Material.cc.

00228 {
00229   //  G4cout << "### Destruction of material " << fName << " started" <<G4endl;
00230   if(!fBaseMaterial) {
00231     if (theElementVector)       { delete    theElementVector; }
00232     if (fMassFractionVector)    { delete [] fMassFractionVector; }
00233     if (fAtomsVector)           { delete [] fAtomsVector; }
00234     if (fSandiaTable)           { delete    fSandiaTable; }
00235   }
00236   if (fIonisation)            { delete    fIonisation; }
00237   if (VecNbOfAtomsPerVolume)  { delete [] VecNbOfAtomsPerVolume; }
00238 
00239   // Remove this material from theMaterialTable.
00240   //
00241   theMaterialTable[fIndexInTable] = 0;
00242 }

G4Material::G4Material ( __void__ &   ) 

Definition at line 216 of file G4Material.cc.

00217   : fNumberOfComponents(0), fNumberOfElements(0), theElementVector(0), 
00218     fImplicitElement(false), fMassFractionVector(0), fAtomsVector(0), 
00219     fMaterialPropertiesTable(0), fIndexInTable(0), 
00220     VecNbOfAtomsPerVolume(0), fIonisation(0), fSandiaTable(0)
00221 {
00222   InitializePointers();
00223 }


Member Function Documentation

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

Definition at line 384 of file G4Material.cc.

References FatalException, G4cerr, G4cout, G4endl, G4Exception(), GetA(), G4Element::GetName(), and G4Element::increaseCountUse().

00385 {
00386   if(fraction < 0.0 || fraction > 1.0) {
00387     G4cout << "G4Material::AddElement ERROR for " << fName << " and " 
00388            << element->GetName() << "  mass fraction= " << fraction 
00389            << " is wrong " << G4endl;
00390     G4Exception ("G4Material::AddElement()", "mat032", FatalException,     
00391                  "Attempt to add element with wrong mass fraction");
00392   }
00393   // initialization
00394   if (fNumberOfComponents == 0) {
00395     fMassFractionVector = new G4double[fArrayLength];
00396     fAtomsVector        = new G4int   [fArrayLength];
00397   }
00398   // filling ...
00399   if (G4int(fNumberOfComponents) < maxNbComponents) {
00400     size_t el = 0;
00401     while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) { ++el; }
00402     if (el<fNumberOfElements) fMassFractionVector[el] += fraction;
00403     else {
00404       theElementVector->push_back(element); 
00405       fMassFractionVector[el] = fraction;
00406       ++fNumberOfElements;
00407       element->increaseCountUse();
00408     }
00409     ++fNumberOfComponents;  
00410   } else {
00411     G4cout << "G4Material::AddElement ERROR for " << fName << " nElement= " 
00412            <<  fNumberOfElements << G4endl;
00413     G4Exception ("G4Material::AddElement()", "mat033", FatalException, 
00414            "Attempt to add more than the declared number of elements.");
00415   }    
00416 
00417   // filled.
00418   if (G4int(fNumberOfComponents) == maxNbComponents) {
00419 
00420     size_t i=0;
00421     G4double Zmol(0.), Amol(0.);
00422     // check sum of weights -- OK?
00423     G4double wtSum(0.0);
00424     for (i=0; i<fNumberOfElements; ++i) {
00425       wtSum += fMassFractionVector[i];
00426       Zmol +=  fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
00427       Amol +=  fMassFractionVector[i]*(*theElementVector)[i]->GetA();
00428     }
00429     if (std::fabs(1.-wtSum) > perThousand) {
00430       G4cerr << "WARNING !! for " << fName << " sum of fractional masses "
00431              <<  wtSum << " is not 1 - results may be wrong" 
00432              << G4endl;
00433     }
00434     for (i=0; i<fNumberOfElements; ++i) {
00435       fAtomsVector[i] = 
00436         G4int(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA()+0.5);
00437     }
00438      
00439     ComputeDerivedQuantities();
00440   }
00441 }

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

Definition at line 341 of file G4Material.cc.

References FatalException, G4cout, G4endl, G4Exception(), and GetA().

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

00342 {   
00343   // initialization
00344   if ( fNumberOfElements == 0 ) {
00345      fAtomsVector        = new G4int   [fArrayLength];
00346      fMassFractionVector = new G4double[fArrayLength];
00347   }
00348 
00349   // filling ...
00350   if ( G4int(fNumberOfElements) < maxNbComponents ) {
00351      theElementVector->push_back(element);     
00352      fAtomsVector[fNumberOfElements] = nAtoms;
00353      fNumberOfComponents = ++fNumberOfElements;
00354      element->increaseCountUse();
00355   } else {
00356     G4cout << "G4Material::AddElement ERROR for " << fName << " nElement= " 
00357            <<  fNumberOfElements << G4endl;
00358     G4Exception ("G4Material::AddElement()", "mat031", FatalException, 
00359            "Attempt to add more than the declared number of elements.");
00360   } 
00361   // filled.
00362   if ( G4int(fNumberOfElements) == maxNbComponents ) {     
00363     // compute proportion by mass
00364     size_t i=0;
00365     G4double Amol = 0.;
00366     for (i=0; i<fNumberOfElements; ++i) {
00367       G4double w = fAtomsVector[i]*(*theElementVector)[i]->GetA(); 
00368       Amol += w;
00369       fMassFractionVector[i] = w;
00370     }
00371     for (i=0; i<fNumberOfElements; ++i) {
00372       fMassFractionVector[i] /= Amol;
00373     }
00374 
00375     fMassOfMolecule = Amol/Avogadro;
00376     ComputeDerivedQuantities();
00377   }
00378 }

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

Definition at line 447 of file G4Material.cc.

References FatalException, G4cout, G4endl, G4Exception(), GetA(), GetElementVector(), GetFractionVector(), GetName(), and GetNumberOfElements().

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

00448 {
00449   if(fraction < 0.0 || fraction > 1.0) {
00450     G4cout << "G4Material::AddMaterial ERROR for " << fName << " and " 
00451            << material->GetName() << "  mass fraction= " << fraction 
00452            << " is wrong " << G4endl;
00453     G4Exception ("G4Material::AddMaterial()", "mat034", FatalException,            
00454                  "Attempt to add material with wrong mass fraction");      
00455   }
00456   // initialization
00457   if (fNumberOfComponents == 0) {
00458     fMassFractionVector = new G4double[fArrayLength];
00459     fAtomsVector        = new G4int   [fArrayLength];
00460   }
00461 
00462   size_t nelm = material->GetNumberOfElements();
00463 
00464   // arrays should be extended
00465   if(nelm > 1) {
00466     G4int nold    = fArrayLength;
00467     fArrayLength += nelm - 1;
00468     G4double* v1 = new G4double[fArrayLength];
00469     G4int* i1    = new G4int[fArrayLength];
00470     for(G4int i=0; i<nold; ++i) {
00471       v1[i] = fMassFractionVector[i];
00472       i1[i] = fAtomsVector[i];
00473     }
00474     delete [] fAtomsVector;
00475     delete [] fMassFractionVector;
00476     fMassFractionVector = v1;
00477     fAtomsVector = i1;
00478   }
00479 
00480   // filling ...
00481   if (G4int(fNumberOfComponents) < maxNbComponents) {
00482     for (size_t elm=0; elm<nelm; ++elm)
00483       {
00484         G4Element* element = (*(material->GetElementVector()))[elm];
00485         size_t el = 0;
00486         while ((el<fNumberOfElements)&&(element!=(*theElementVector)[el])) el++;
00487         if (el < fNumberOfElements) fMassFractionVector[el] += fraction
00488                                           *(material->GetFractionVector())[elm];
00489         else {
00490           theElementVector->push_back(element); 
00491           fMassFractionVector[el] = fraction
00492                                           *(material->GetFractionVector())[elm];
00493           ++fNumberOfElements;
00494           element->increaseCountUse();
00495         }
00496       } 
00497     ++fNumberOfComponents;
00499     fMatComponents[material] = fraction;
00500       
00501   } else {
00502     G4cout << "G4Material::AddMaterial ERROR for " << fName << " nElement= " 
00503            <<  fNumberOfElements << G4endl;
00504     G4Exception ("G4Material::AddMaterial()", "mat035", FatalException, 
00505            "Attempt to add more than the declared number of components.");
00506   }    
00507 
00508   // filled.
00509   if (G4int(fNumberOfComponents) == maxNbComponents) {
00510     size_t i=0;
00511     G4double Zmol(0.), Amol(0.);
00512     // check sum of weights -- OK?
00513     G4double wtSum(0.0);
00514     for (i=0; i<fNumberOfElements; ++i) {
00515       wtSum += fMassFractionVector[i];
00516       Zmol +=  fMassFractionVector[i]*(*theElementVector)[i]->GetZ();
00517       Amol +=  fMassFractionVector[i]*(*theElementVector)[i]->GetA();
00518     }
00519     if (std::fabs(1.-wtSum) > perThousand) {
00520       G4cout << "G4Material::AddMaterial WARNING !! for " << fName 
00521              << " sum of fractional masses "
00522              <<  wtSum << " is not 1 - results may be wrong" 
00523              << G4endl;
00524     }
00525     for (i=0;i<fNumberOfElements;i++) {
00526       fAtomsVector[i] = 
00527         G4int(fMassFractionVector[i]*Amol/(*theElementVector)[i]->GetA()+0.5);
00528     }
00529      
00530     ComputeDerivedQuantities();
00531   }
00532 }

G4double G4Material::GetA (  )  const

Definition at line 617 of file G4Material.cc.

References FatalException, G4cout, G4endl, and G4Exception().

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

00618 { 
00619   if (fNumberOfElements > 1) { 
00620      G4cout << "G4Material ERROR in GetA. The material: " << fName << " is a mixture."
00621             << G4endl;
00622      G4Exception ("G4Material::GetA()", "mat037", FatalException,  
00623                   "the Atomic mass is not well defined." );
00624   } 
00625   return  (*theElementVector)[0]->GetA();      
00626 }

const G4double* G4Material::GetAtomicNumDensityVector (  )  const [inline]

Definition at line 215 of file G4Material.hh.

Referenced by G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4RToEConvForGamma::BuildAbsorptionLengthVector(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4CrossSectionHandler::BuildCrossSectionsForMaterials(), G4VRangeToEnergyConverter::BuildRangeVector(), G4MuPairProductionModel::ComputeDEDXPerVolume(), G4MuBremsstrahlungModel::ComputeDEDXPerVolume(), G4LivermoreIonisationModel::ComputeDEDXPerVolume(), G4LivermoreBremsstrahlungModel::ComputeDEDXPerVolume(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4eBremsstrahlungModel::ComputeDEDXPerVolume(), G4eBremParametrizedModel::ComputeDEDXPerVolume(), G4eBremsstrahlungModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4AdjointBremsstrahlungModel::DiffCrossSectionPerVolumePrimToSecondApproximated2(), G4StopElementSelector::GetElement(), G4ElementSelector::SelectZandA(), and G4hQAOModel::StoppingPower().

00215 {return VecNbOfAtomsPerVolume;}

const G4int* G4Material::GetAtomsVector (  )  const [inline]

Definition at line 197 of file G4Material.hh.

00197 {return fAtomsVector;}

const G4Material* G4Material::GetBaseMaterial (  )  const [inline]

Definition at line 232 of file G4Material.hh.

Referenced by G4MuElecInelasticModel::CrossSectionPerVolume(), G4MuElecElasticModel::CrossSectionPerVolume(), G4LossTableBuilder::InitialiseBaseMaterials(), and CompareMaterial::operator()().

00232 {return fBaseMaterial;}

const G4String& G4Material::GetChemicalFormula (  )  const [inline]

Definition at line 178 of file G4Material.hh.

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

00178 {return fChemicalFormula;}

G4double G4Material::GetDensity (  )  const [inline]

Definition at line 179 of file G4Material.hh.

Referenced by G4GMocrenFileSceneHandler::AddSolid(), G4NistManager::BuildMaterialWithNewDensity(), G4ProductionCutsTable::CheckMaterialInfo(), G4EmCalculator::ComputeDEDX(), G4mplIonisationWithDeltaModel::ComputeDEDXPerVolume(), G4mplIonisationModel::ComputeDEDXPerVolume(), G4EmCalculator::ComputeNuclearDEDX(), G4NistMaterialBuilder::ConstructNewGasMaterial(), G4VRangeToEnergyConverter::Convert(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G4PhysicalVolumeModel::DescribeAndDescend(), G4tgbGeometryDumper::DumpMaterial(), G4InitXscPAI::G4InitXscPAI(), G4PAIxSection::G4PAIxSection(), G4EmCalculator::GetDEDX(), G4VMscModel::GetEnergy(), G4LogicalVolume::GetMass(), G4VMscModel::GetRange(), G4LossTableBuilder::InitialiseBaseMaterials(), G4PAIySection::Initialize(), G4GDMLWriteMaterials::MaterialWrite(), G4PSDoseDeposit::ProcessHits(), G4ASCIITreeSceneHandler::RequestPrimitives(), GFlashSamplingShowerParameterisation::SetMaterial(), GFlashHomoShowerParameterisation::SetMaterial(), and G4tgbMaterialMixtureByVolume::TransformToFractionsByWeight().

00179 {return fDensity;}

G4double G4Material::GetElectronDensity (  )  const [inline]

Definition at line 216 of file G4Material.hh.

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

00216 {return TotNbOfElectPerVolume;}

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

Definition at line 201 of file G4Material.hh.

Referenced by G4NeutronHPThermalScattering::ApplyYourself(), G4NeutronHPorLFission::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLEInelastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4FissLib::ApplyYourself(), G4AdjointCSManager::ComputeAdjointCS(), G4PAIySection::ComputeLowEnergyCof(), GVFlashShowerParameterisation::GetEffA(), GVFlashShowerParameterisation::GetEffZ(), G4NeutronIsotopeProduction::GetIsotope(), G4GDMLWriteMaterials::MaterialWrite(), G4AdjointCSManager::SampleElementFromCSMatrices(), and G4VEmAdjointModel::SelectCSMatrix().

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

const G4ElementVector* G4Material::GetElementVector (  )  const [inline]

Definition at line 189 of file G4Material.hh.

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

00189 {return theElementVector;}

const G4double* G4Material::GetFractionVector (  )  const [inline]

Definition at line 193 of file G4Material.hh.

Referenced by AddMaterial(), G4QCaptureAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4tgbGeometryDumper::DumpMaterial(), GVFlashShowerParameterisation::GetEffA(), GVFlashShowerParameterisation::GetEffZ(), G4GDMLWriteMaterials::MaterialWrite(), and G4QAtomicElectronScattering::PostStepDoIt().

00193 {return fMassFractionVector;}

size_t G4Material::GetIndex (  )  const [inline]

Definition at line 261 of file G4Material.hh.

Referenced by G4AdjointCSManager::ComputeAdjointCS(), G4VRangeToEnergyConverter::Convert(), G4DNATransformElectronModel::CrossSectionPerVolume(), G4DNAScreenedRutherfordElasticModel::CrossSectionPerVolume(), G4DNASancheSolvatationModel::CrossSectionPerVolume(), G4DNASancheExcitationModel::CrossSectionPerVolume(), G4DNARuddIonisationModel::CrossSectionPerVolume(), G4DNARuddIonisationExtendedModel::CrossSectionPerVolume(), G4DNAMillerGreenExcitationModel::CrossSectionPerVolume(), G4DNAMeltonAttachmentModel::CrossSectionPerVolume(), G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume(), G4DNADingfelderChargeIncreaseModel::CrossSectionPerVolume(), G4DNADingfelderChargeDecreaseModel::CrossSectionPerVolume(), G4DNAChampionElasticModel::CrossSectionPerVolume(), G4DNABornIonisationModel::CrossSectionPerVolume(), G4DNABornExcitationModel::CrossSectionPerVolume(), 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().

00261 {return fIndexInTable;}

G4IonisParamMat* G4Material::GetIonisation (  )  const [inline]

Definition at line 225 of file G4Material.hh.

Referenced by G4EmCorrections::Bethe(), G4hImpactIonisation::BuildPhysicsTable(), G4MuBetheBlochModel::ComputeDEDXPerVolume(), G4MollerBhabhaModel::ComputeDEDXPerVolume(), G4BetheBlochModel::ComputeDEDXPerVolume(), G4EmCorrections::DensityCorrection(), G4tgbGeometryDumper::DumpMaterial(), G4ionEffectiveCharge::EffectiveCharge(), G4ElectronIonPair::FindG4MeanEnergyPerIonPair(), G4tgbMaterialMgr::FindOrBuildG4Material(), G4VEnergyLoss::GetLossWithFluct(), G4VeLowEnergyLoss::GetLossWithFluct(), G4hBetheBlochModel::LowEnergyLimit(), G4GDMLReadMaterials::MaterialRead(), G4GDMLWriteMaterials::MaterialWrite(), G4ElectronIonPair::MeanNumberOfIonsAlongStep(), G4MuBetheBlochModel::MinEnergyCut(), G4CoulombScattering::MinPrimaryEnergy(), operator<<(), G4hImpactIonisation::PrintInfoDefinition(), G4UniversalFluctuation::SampleFluctuations(), G4WentzelVIRelXSection::SetupKinematic(), G4WentzelOKandVIxSection::SetupKinematic(), and G4EmCorrections::ShellCorrectionSTD().

00225 {return fIonisation;}

G4double G4Material::GetMassOfMolecule (  )  const [inline]

Definition at line 241 of file G4Material.hh.

Referenced by G4DNASecondOrderReaction::BuildPhysicsTable(), G4Material(), and G4DNAMolecularMaterial::SearchMolecularMaterial().

00241 {return fMassOfMolecule;}

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

Definition at line 236 of file G4Material.hh.

Referenced by G4DNAMolecularMaterial::SearchMolecularMaterial().

00237                                                {return fMatComponents;}

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

Definition at line 576 of file G4Material.cc.

References G4cout, G4endl, and GetName().

Referenced by G4DNABrownianTransportation::BuildPhysicsTable(), G4ProductionCutsTable::CheckMaterialInfo(), G4EmCalculator::FindMaterial(), G3MatTable::get(), G4GDMLReadMaterials::GetMaterial(), G4DNATransformElectronModel::Initialise(), G4DNAScreenedRutherfordElasticModel::Initialise(), G4DNASancheSolvatationModel::Initialise(), G4DNASancheExcitationModel::Initialise(), G4DNARuddIonisationModel::Initialise(), G4DNARuddIonisationExtendedModel::Initialise(), G4DNAMillerGreenExcitationModel::Initialise(), G4DNAMeltonAttachmentModel::Initialise(), G4DNAEmfietzoglouExcitationModel::Initialise(), G4DNADingfelderChargeIncreaseModel::Initialise(), G4DNADingfelderChargeDecreaseModel::Initialise(), G4DNAChampionElasticModel::Initialise(), G4DNABornIonisationModel::Initialise(), G4DNABornExcitationModel::Initialise(), and G4Region::ScanVolumeTree().

00577 {  
00578   // search the material by its name 
00579   for (size_t J=0 ; J<theMaterialTable.size() ; ++J)
00580     {
00581       if (theMaterialTable[J]->GetName() == materialName)
00582         { return theMaterialTable[J]; }
00583     }
00584    
00585   // the material does not exist in the table
00586   if (warning) {
00587     G4cout << "G4Material::GetMaterial() WARNING: The material: "
00588            << materialName << " does not exist in the table. Return NULL pointer."
00589            << G4endl;
00590   }      
00591   return 0;          
00592 }

G4MaterialPropertiesTable* G4Material::GetMaterialPropertiesTable (  )  const [inline]

Definition at line 251 of file G4Material.hh.

Referenced by G4Scintillation::BuildThePhysicsTable(), G4Track::CalculateVelocityForOpticalPhoton(), G4Material(), G4OpWLS::GetMeanFreePath(), G4OpRayleigh::GetMeanFreePath(), G4OpMieHG::GetMeanFreePath(), G4OpAbsorption::GetMeanFreePath(), G4GDMLWriteMaterials::MaterialWrite(), G4Scintillation::PostStepDoIt(), G4OpWLS::PostStepDoIt(), G4OpMieHG::PostStepDoIt(), G4OpBoundaryProcess::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), G4Cerenkov::PostStepGetPhysicalInteractionLength(), G4GDMLReadMaterials::PropertyRead(), and G4GDMLWriteMaterials::PropertyWrite().

00252   {return fMaterialPropertiesTable;}

const G4MaterialTable * G4Material::GetMaterialTable (  )  [static]

Definition at line 562 of file G4Material.cc.

Referenced by G4VCrossSectionHandler::ActiveElements(), G4AugerData::BuildAugerTransitionTable(), G4AdjointCSManager::BuildCrossSectionMatrices(), G4NeutronHPThermalScatteringData::BuildPhysicsTable(), G4Scintillation::BuildThePhysicsTable(), G4VXTRenergyLoss::ComputeGasPhotoAbsCof(), G4StrawTubeXTRadiator::ComputeMediumPhotoAbsCof(), G4VXTRenergyLoss::ComputePlatePhotoAbsCof(), G4PAIPhotonModel::ComputeSandiaPhotoAbsCof(), G4PAIModel::ComputeSandiaPhotoAbsCof(), G4VRangeToEnergyConverter::Convert(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), G4NistMaterialBuilder::FindOrBuildMaterial(), G4PAIxSection::G4PAIxSection(), G4SandiaTable::G4SandiaTable(), G4DNAMolecularMaterial::GetDensityTableFor(), G4VXTRenergyLoss::GetGasCompton(), G4NistManager::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VXTRenergyLoss::GetPlateCompton(), G4PAIPhotonModel::Initialise(), G4PAIModel::Initialise(), G4ICRU73QOModel::Initialise(), G4DNAMolecularMaterial::Initialize(), G4DNAMolecularMaterial::InitializeDensity(), G4DNAMolecularMaterial::InitializeNumMolPerVol(), G4NistManager::PrintG4Material(), G4ProductionCutsTable::StoreMaterialInfo(), G4GDMLRead::StripNames(), and G4NistManager::~G4NistManager().

00563 {
00564   return &theMaterialTable;
00565 }

const G4String& G4Material::GetName (  )  const [inline]

Definition at line 177 of file G4Material.hh.

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

00177 {return fName;}

G4double G4Material::GetNuclearInterLength (  )  const [inline]

Definition at line 222 of file G4Material.hh.

Referenced by G4MSSteppingAction::UserSteppingAction().

00222 {return fNuclInterLen;}

size_t G4Material::GetNumberOfElements (  )  const [inline]

Definition at line 185 of file G4Material.hh.

Referenced by AddMaterial(), G4AdjointPhotoElectricModel::AdjointCrossSection(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4NeutronHPThermalScattering::ApplyYourself(), G4NeutronHPorLFission::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLEInelastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4FissLib::ApplyYourself(), G4QCaptureAtRest::AtRestDoIt(), G4PionMinusAbsorptionAtRest::AtRestDoIt(), G4PiMinusAbsorptionAtRest::AtRestDoIt(), G4NeutronCaptureAtRest::AtRestDoIt(), G4KaonMinusAbsorption::AtRestDoIt(), G4AntiProtonAnnihilationAtRest::AtRestDoIt(), G4AntiNeutronAnnihilationAtRest::AtRestDoIt(), G4RToEConvForGamma::BuildAbsorptionLengthVector(), G4eIonisationCrossSectionHandler::BuildCrossSectionsForMaterials(), G4CrossSectionHandler::BuildCrossSectionsForMaterials(), G4BremsstrahlungCrossSectionHandler::BuildCrossSectionsForMaterials(), G4VRangeToEnergyConverter::BuildRangeVector(), G4Nucleus::ChooseParameters(), G4AdjointCSManager::ComputeAdjointCS(), G4MuPairProductionModel::ComputeDEDXPerVolume(), G4MuBremsstrahlungModel::ComputeDEDXPerVolume(), G4LivermoreIonisationModel::ComputeDEDXPerVolume(), G4LivermoreBremsstrahlungModel::ComputeDEDXPerVolume(), G4ICRU49NuclearStoppingModel::ComputeDEDXPerVolume(), G4eBremsstrahlungRelModel::ComputeDEDXPerVolume(), G4eBremsstrahlungModel::ComputeDEDXPerVolume(), G4eBremParametrizedModel::ComputeDEDXPerVolume(), G4PAIySection::ComputeLowEnergyCof(), G4GammaConversionToMuons::ComputeMeanFreePath(), G4VEmModel::CrossSectionPerVolume(), G4eBremsstrahlungModel::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(), G4QNGamma::GetMeanFreePath(), G4QLowEnergy::GetMeanFreePath(), G4QIonIonElastic::GetMeanFreePath(), G4QInelastic::GetMeanFreePath(), G4QElastic::GetMeanFreePath(), G4QDiffraction::GetMeanFreePath(), G4QCoherentChargeExchange::GetMeanFreePath(), G4QAtomicElectronScattering::GetMeanFreePath(), G4PiMinusAbsorptionAtRest::GetMeanLifeTime(), G4KaonMinusAbsorptionAtRest::GetMeanLifeTime(), G4hZiegler1985p::HasMaterial(), G4hZiegler1977p::HasMaterial(), G4hZiegler1977He::HasMaterial(), G4hSRIM2000p::HasMaterial(), G4hICRU49p::HasMaterial(), G4hICRU49He::HasMaterial(), G4LivermoreRayleighModel::Initialise(), G4LivermorePhotoElectricModel::Initialise(), G4VAtomDeexcitation::InitialiseAtomicDeexcitation(), G4QAOLowEnergyLoss::IsInCharge(), G4GDMLWriteMaterials::MaterialWrite(), G4QNGamma::PostStepDoIt(), G4QLowEnergy::PostStepDoIt(), G4QIonIonElastic::PostStepDoIt(), G4QInelastic::PostStepDoIt(), G4QElastic::PostStepDoIt(), G4QDiffraction::PostStepDoIt(), G4QCoherentChargeExchange::PostStepDoIt(), G4QAtomicElectronScattering::PostStepDoIt(), G4WentzelVIRelModel::SampleScattering(), G4WentzelVIModel::SampleScattering(), G4GoudsmitSaundersonMscModel::SampleScattering(), G4AdjointPhotoElectricModel::SampleSecondaries(), G4CrossSectionDataStore::SampleZandA(), G4VEmModel::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomAtom(), G4PixeCrossSectionHandler::SelectRandomAtom(), G4eBremsstrahlungModel::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomElement(), G4ElementSelector::SelectZandA(), G4hZiegler1985p::StoppingPower(), G4hZiegler1977p::StoppingPower(), G4hZiegler1977He::StoppingPower(), G4hSRIM2000p::StoppingPower(), G4hQAOModel::StoppingPower(), G4hICRU49p::StoppingPower(), G4hICRU49He::StoppingPower(), G4VCrossSectionHandler::ValueForMaterial(), and G4PixeCrossSectionHandler::ValueForMaterial().

00185 {return fNumberOfElements;}

size_t G4Material::GetNumberOfMaterials (  )  [static]

Definition at line 569 of file G4Material.cc.

Referenced by G4VCrossSectionHandler::ActiveElements(), G4AugerData::BuildAugerTransitionTable(), G4NeutronHPThermalScatteringData::BuildPhysicsTable(), G4Scintillation::BuildThePhysicsTable(), G4VEnergyLoss::CopyCutVectors(), G4ElectronIonPair::DumpMeanEnergyPerIonPair(), G4VEnergyLoss::EqualCutVectors(), G4SandiaTable::G4SandiaTable(), G4PAIPhotonModel::Initialise(), and G4PAIModel::Initialise().

00570 {
00571   return theMaterialTable.size();
00572 }

G4double G4Material::GetPressure (  )  const [inline]

Definition at line 182 of file G4Material.hh.

Referenced by G4NistManager::BuildMaterialWithNewDensity(), G4tgbGeometryDumper::DumpMaterial(), and G4GDMLWriteMaterials::MaterialWrite().

00182 {return fPressure;}

G4double G4Material::GetRadlen (  )  const [inline]

Definition at line 219 of file G4Material.hh.

Referenced by G4WentzelVIRelModel::ComputeTruePathLengthLimit(), G4WentzelVIModel::ComputeTruePathLengthLimit(), G4PhysicalVolumeModel::CreateCurrentAttValues(), GFlashSamplingShowerParameterisation::SetMaterial(), GFlashHomoShowerParameterisation::SetMaterial(), G4PairProductionRelModel::SetupForMaterial(), G4eBremsstrahlungRelModel::SetupForMaterial(), and G4MSSteppingAction::UserSteppingAction().

00219 {return fRadlen;}

G4SandiaTable* G4Material::GetSandiaTable (  )  const [inline]

Definition at line 228 of file G4Material.hh.

Referenced by G4PAIModel::BuildPAIonisationTable(), G4VXTRenergyLoss::ComputeGasPhotoAbsCof(), G4StrawTubeXTRadiator::ComputeMediumPhotoAbsCof(), G4VXTRenergyLoss::ComputePlatePhotoAbsCof(), G4PEEffectModel::CrossSectionPerVolume(), G4PEEffectFluoModel::CrossSectionPerVolume(), and G4PAIySection::Initialize().

00228 {return fSandiaTable;}

G4State G4Material::GetState (  )  const [inline]

Definition at line 180 of file G4Material.hh.

Referenced by G4NistManager::BuildMaterialWithNewDensity(), G4NistMaterialBuilder::ConstructNewGasMaterial(), G4PhysicalVolumeModel::CreateCurrentAttValues(), G4tgbGeometryDumper::DumpMaterial(), G4ForwardXrayTR::GetEnergyTR(), G4hICRU49p::HasMaterial(), G4hICRU49He::HasMaterial(), G4GDMLWriteMaterials::MaterialWrite(), and G4ForwardXrayTR::PostStepDoIt().

00180 {return fState;}

G4double G4Material::GetTemperature (  )  const [inline]

Definition at line 181 of file G4Material.hh.

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

00181 {return fTemp;}

G4double G4Material::GetTotNbOfAtomsPerVolume (  )  const [inline]

Definition at line 208 of file G4Material.hh.

Referenced by G4EmCorrections::BarkasCorrection(), G4Nucleus::ChooseParameters(), G4PenelopeIonisationModel::ComputeDEDXPerVolume(), G4PenelopeBremsstrahlungModel::ComputeDEDXPerVolume(), G4MollerBhabhaModel::ComputeDEDXPerVolume(), G4PenelopeIonisationModel::CrossSectionPerVolume(), G4PenelopeComptonModel::CrossSectionPerVolume(), G4PenelopeBremsstrahlungModel::CrossSectionPerVolume(), G4MuElecInelasticModel::CrossSectionPerVolume(), G4MuElecElasticModel::CrossSectionPerVolume(), G4IonFluctuations::Dispersion(), G4EmCorrections::KShellCorrection(), G4EmCorrections::LShellCorrection(), G4EmCorrections::ShellCorrection(), and G4hICRU49He::StoppingPower().

00208 {return TotNbOfAtomsPerVolume;}

G4double G4Material::GetTotNbOfElectPerVolume (  )  const [inline]

Definition at line 211 of file G4Material.hh.

Referenced by G4IonParametrisedLossModel::DeltaRayMeanEnergyTransferRate(), and G4hICRU49He::StoppingPower().

00211 {return TotNbOfElectPerVolume;}

const G4double* G4Material::GetVecNbOfAtomsPerVolume (  )  const [inline]

Definition at line 205 of file G4Material.hh.

Referenced by G4AdjointPhotoElectricModel::AdjointCrossSection(), G4VAtomDeexcitation::AlongStepDeexcitation(), G4NeutronHPorLFission::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), G4NeutronHPorLEInelastic::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPCapture::ApplyYourself(), G4FissLib::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(), G4QNGamma::GetMeanFreePath(), G4QLowEnergy::GetMeanFreePath(), G4QIonIonElastic::GetMeanFreePath(), G4QInelastic::GetMeanFreePath(), G4QElastic::GetMeanFreePath(), G4QDiffraction::GetMeanFreePath(), G4QCoherentChargeExchange::GetMeanFreePath(), G4QAtomicElectronScattering::GetMeanFreePath(), G4EmElementSelector::Initialise(), G4GoudsmitSaundersonMscModel::SampleScattering(), G4VCrossSectionHandler::ValueForMaterial(), and G4PixeCrossSectionHandler::ValueForMaterial().

00205 {return VecNbOfAtomsPerVolume;}

G4double G4Material::GetZ (  )  const

Definition at line 604 of file G4Material.cc.

References FatalException, G4cout, G4endl, and G4Exception().

Referenced by G4tgbGeometryDumper::DumpMaterial(), GVFlashShowerParameterisation::GetEffZ(), G4PiMinusAbsorptionAtRest::GetMeanLifeTime(), G4KaonMinusAbsorptionAtRest::GetMeanLifeTime(), G4GDMLWriteMaterials::MaterialWrite(), G4VCrossSectionHandler::SelectRandomAtom(), G4PixeCrossSectionHandler::SelectRandomAtom(), G4hZiegler1985p::StoppingPower(), G4hZiegler1977p::StoppingPower(), G4hZiegler1977He::StoppingPower(), G4hSRIM2000p::StoppingPower(), G4hICRU49p::StoppingPower(), and G4hICRU49He::StoppingPower().

00605 { 
00606   if (fNumberOfElements > 1) {
00607      G4cout << "G4Material ERROR in GetZ. The material: " << fName << " is a mixture."
00608             << G4endl;
00609      G4Exception ("G4Material::GetZ()", "mat036", FatalException, 
00610                   "the Atomic number is not well defined." );
00611   } 
00612   return (*theElementVector)[0]->GetZ();      
00613 }

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

Definition at line 689 of file G4Material.cc.

00690 {
00691   return (this != (G4Material *) &right);
00692 }

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

Definition at line 682 of file G4Material.cc.

00683 {
00684   return (this == (G4Material *) &right);
00685 }

void G4Material::SetChemicalFormula ( const G4String chF  )  [inline]

Definition at line 172 of file G4Material.hh.

00172 {fChemicalFormula=chF;}

void G4Material::SetMaterialPropertiesTable ( G4MaterialPropertiesTable anMPT  )  [inline]

Definition at line 248 of file G4Material.hh.

Referenced by G4GDMLReadMaterials::PropertyRead().

00249   {fMaterialPropertiesTable = anMPT;}

void G4Material::SetName ( const G4String name  )  [inline]

Definition at line 282 of file G4Material.hh.

00282 {fName=name;}


Friends And Related Function Documentation

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

Definition at line 748 of file G4Material.cc.

00749 {
00750   //Dump info for all known materials
00751   flux << "\n***** Table : Nb of materials = " << MaterialTable.size() 
00752        << " *****\n" << G4endl;
00753         
00754   for (size_t i=0; i<MaterialTable.size(); ++i) { 
00755     flux << MaterialTable[i] << G4endl << G4endl; 
00756   }
00757 
00758   return flux;
00759 }      

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

Definition at line 740 of file G4Material.cc.

00741 {
00742   flux << &material;        
00743   return flux;
00744 }

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

Definition at line 697 of file G4Material.cc.

00698 {
00699   std::ios::fmtflags mode = flux.flags();
00700   flux.setf(std::ios::fixed,std::ios::floatfield);
00701   G4long prec = flux.precision(3);
00702   
00703   flux
00704     << " Material: "         << std::setw(8) <<  material->fName
00705     << " " << material->fChemicalFormula << " "
00706     << "  density: "         << std::setw(6) << std::setprecision(3)  
00707     << G4BestUnit(material->fDensity,"Volumic Mass") 
00708     << "  RadL: "            << std::setw(7)  << std::setprecision(3)  
00709     << G4BestUnit(material->fRadlen,"Length")
00710     << "  Nucl.Int.Length: " << std::setw(7)  << std::setprecision(3)  
00711     << G4BestUnit(material->fNuclInterLen,"Length")    
00712     << "  Imean: "           << std::setw(7)  << std::setprecision(3)  
00713     << G4BestUnit(material->GetIonisation()->GetMeanExcitationEnergy(),"Energy");
00714     
00715   if(material->fState == kStateGas) {
00716     flux
00717       << "  temperature: " << std::setw(6) << std::setprecision(2)  
00718       << (material->fTemp)/kelvin << " K"
00719       << "  pressure: "    << std::setw(6) << std::setprecision(2)   
00720       << (material->fPressure)/atmosphere << " atm";
00721   }
00722   for (size_t i=0; i<material->fNumberOfElements; i++) {
00723     flux 
00724       << "\n   ---> " << (*(material->theElementVector))[i] 
00725       << "\n          ElmMassFraction: " 
00726       << std::setw(6)<< std::setprecision(2) 
00727       << (material->fMassFractionVector[i])/perCent << " %" 
00728       << "  ElmAbundance "     << std::setw(6)<< std::setprecision(2) 
00729       << 100*(material->VecNbOfAtomsPerVolume[i])/(material->TotNbOfAtomsPerVolume)
00730       << " % \n";
00731   }
00732   flux.precision(prec);    
00733   flux.setf(mode,std::ios::floatfield);
00734             
00735   return flux;
00736 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:28 2013 for Geant4 by  doxygen 1.4.7