G4GDMLWriteMaterials Class Reference

#include <G4GDMLWriteMaterials.hh>

Inheritance diagram for G4GDMLWriteMaterials:

G4GDMLWriteDefine G4GDMLWrite G4GDMLWriteSolids G4GDMLWriteSetup G4GDMLWriteParamvol G4GDMLWriteStructure

Public Member Functions

void AddIsotope (const G4Isotope *const)
void AddElement (const G4Element *const)
void AddMaterial (const G4Material *const)
virtual void MaterialsWrite (xercesc::DOMElement *)

Protected Member Functions

 G4GDMLWriteMaterials ()
virtual ~G4GDMLWriteMaterials ()
void AtomWrite (xercesc::DOMElement *, const G4double &)
void DWrite (xercesc::DOMElement *, const G4double &)
void PWrite (xercesc::DOMElement *, const G4double &)
void TWrite (xercesc::DOMElement *, const G4double &)
void MEEWrite (xercesc::DOMElement *, const G4double &)
void IsotopeWrite (const G4Isotope *const)
void ElementWrite (const G4Element *const)
void MaterialWrite (const G4Material *const)
void PropertyWrite (xercesc::DOMElement *, const G4Material *const)
void PropertyVectorWrite (const G4String &, const G4PhysicsOrderedFreeVector *const)

Protected Attributes

std::vector< const G4Isotope * > isotopeList
std::vector< const G4Element * > elementList
std::vector< const G4Material * > materialList
xercesc::DOMElement * materialsElement

Detailed Description

Definition at line 53 of file G4GDMLWriteMaterials.hh.


Constructor & Destructor Documentation

G4GDMLWriteMaterials::G4GDMLWriteMaterials (  )  [protected]

Definition at line 45 of file G4GDMLWriteMaterials.cc.

00046   : G4GDMLWriteDefine(), materialsElement(0)
00047 {
00048 }

G4GDMLWriteMaterials::~G4GDMLWriteMaterials (  )  [protected, virtual]

Definition at line 50 of file G4GDMLWriteMaterials.cc.

00051 {
00052 }


Member Function Documentation

void G4GDMLWriteMaterials::AddElement ( const G4Element const  ) 

Definition at line 296 of file G4GDMLWriteMaterials.cc.

References elementList, and ElementWrite().

Referenced by MaterialWrite().

00297 {
00298    for (size_t i=0;i<elementList.size();i++)     // Check if element is
00299    {                                             // already in the list!
00300       if (elementList[i] == elementPtr) { return; }
00301    }
00302    elementList.push_back(elementPtr);
00303    ElementWrite(elementPtr);
00304 }

void G4GDMLWriteMaterials::AddIsotope ( const G4Isotope const  ) 

Definition at line 286 of file G4GDMLWriteMaterials.cc.

References isotopeList, and IsotopeWrite().

Referenced by ElementWrite().

00287 {
00288    for (size_t i=0; i<isotopeList.size(); i++)   // Check if isotope is
00289    {                                             // already in the list!
00290      if (isotopeList[i] == isotopePtr)  { return; }
00291    }
00292    isotopeList.push_back(isotopePtr);
00293    IsotopeWrite(isotopePtr);
00294 }

void G4GDMLWriteMaterials::AddMaterial ( const G4Material const  ) 

Definition at line 306 of file G4GDMLWriteMaterials.cc.

References materialList, and MaterialWrite().

Referenced by G4GDMLWriteStructure::TraverseVolumeTree().

00307 {
00308    for (size_t i=0;i<materialList.size();i++)    // Check if material is
00309    {                                             // already in the list!
00310       if (materialList[i] == materialPtr)  { return; }
00311    }
00312    materialList.push_back(materialPtr);
00313    MaterialWrite(materialPtr);
00314 }

void G4GDMLWriteMaterials::AtomWrite ( xercesc::DOMElement *  ,
const G4double  
) [protected]

Definition at line 55 of file G4GDMLWriteMaterials.cc.

References G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by ElementWrite(), IsotopeWrite(), and MaterialWrite().

00056 {
00057    xercesc::DOMElement* atomElement = NewElement("atom");
00058    atomElement->setAttributeNode(NewAttribute("unit","g/mole"));
00059    atomElement->setAttributeNode(NewAttribute("value",a*mole/g));
00060    element->appendChild(atomElement);
00061 }

void G4GDMLWriteMaterials::DWrite ( xercesc::DOMElement *  ,
const G4double  
) [protected]

Definition at line 64 of file G4GDMLWriteMaterials.cc.

References G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by MaterialWrite().

00065 {
00066    xercesc::DOMElement* DElement = NewElement("D");
00067    DElement->setAttributeNode(NewAttribute("unit","g/cm3"));
00068    DElement->setAttributeNode(NewAttribute("value",d*cm3/g));
00069    element->appendChild(DElement);
00070 }

void G4GDMLWriteMaterials::ElementWrite ( const G4Element const  )  [protected]

Definition at line 112 of file G4GDMLWriteMaterials.cc.

References AddIsotope(), AtomWrite(), G4GDMLWrite::GenerateName(), G4Element::GetA(), G4Element::GetIsotope(), G4Isotope::GetName(), G4Element::GetName(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZ(), materialsElement, G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by AddElement().

00113 {
00114    const G4String name = GenerateName(elementPtr->GetName(),elementPtr);
00115 
00116    xercesc::DOMElement* elementElement = NewElement("element");
00117    elementElement->setAttributeNode(NewAttribute("name",name));
00118 
00119    const size_t NumberOfIsotopes = elementPtr->GetNumberOfIsotopes();
00120 
00121    if (NumberOfIsotopes>0)
00122    {
00123       const G4double* RelativeAbundanceVector =
00124             elementPtr->GetRelativeAbundanceVector();             
00125       for (size_t i=0;i<NumberOfIsotopes;i++)
00126       {
00127          G4String fractionref = GenerateName(elementPtr->GetIsotope(i)->GetName(),
00128                                              elementPtr->GetIsotope(i));
00129          xercesc::DOMElement* fractionElement = NewElement("fraction");
00130          fractionElement->setAttributeNode(NewAttribute("n",
00131                                            RelativeAbundanceVector[i]));
00132          fractionElement->setAttributeNode(NewAttribute("ref",fractionref));
00133          elementElement->appendChild(fractionElement);
00134          AddIsotope(elementPtr->GetIsotope(i));
00135       }
00136    }
00137    else
00138    {
00139       elementElement->setAttributeNode(NewAttribute("Z",elementPtr->GetZ()));
00140       AtomWrite(elementElement,elementPtr->GetA());
00141    }
00142 
00143    materialsElement->appendChild(elementElement);
00144      // Append the element AFTER all the possible components are appended!
00145 }

void G4GDMLWriteMaterials::IsotopeWrite ( const G4Isotope const  )  [protected]

Definition at line 100 of file G4GDMLWriteMaterials.cc.

References AtomWrite(), G4GDMLWrite::GenerateName(), G4Isotope::GetA(), G4Isotope::GetN(), G4Isotope::GetName(), G4Isotope::GetZ(), materialsElement, G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by AddIsotope().

00101 {
00102    const G4String name = GenerateName(isotopePtr->GetName(),isotopePtr);
00103 
00104    xercesc::DOMElement* isotopeElement = NewElement("isotope");
00105    isotopeElement->setAttributeNode(NewAttribute("name",name));
00106    isotopeElement->setAttributeNode(NewAttribute("N",isotopePtr->GetN()));
00107    isotopeElement->setAttributeNode(NewAttribute("Z",isotopePtr->GetZ()));
00108    materialsElement->appendChild(isotopeElement);
00109    AtomWrite(isotopeElement,isotopePtr->GetA());
00110 }

void G4GDMLWriteMaterials::MaterialsWrite ( xercesc::DOMElement *   )  [virtual]

Implements G4GDMLWrite.

Definition at line 274 of file G4GDMLWriteMaterials.cc.

References elementList, G4cout, G4endl, isotopeList, materialList, materialsElement, and G4GDMLWrite::NewElement().

00275 {
00276    G4cout << "G4GDML: Writing materials..." << G4endl;
00277 
00278    materialsElement = NewElement("materials");
00279    element->appendChild(materialsElement);
00280 
00281    isotopeList.clear();
00282    elementList.clear();
00283    materialList.clear();
00284 }

void G4GDMLWriteMaterials::MaterialWrite ( const G4Material const  )  [protected]

Definition at line 147 of file G4GDMLWriteMaterials.cc.

References AddElement(), AtomWrite(), DWrite(), G4GDMLWrite::GenerateName(), G4Material::GetA(), G4Material::GetDensity(), G4Material::GetElement(), G4Material::GetFractionVector(), G4Material::GetIonisation(), G4Material::GetMaterialPropertiesTable(), G4IonisParamMat::GetMeanExcitationEnergy(), G4Element::GetName(), G4Material::GetName(), G4Material::GetNumberOfElements(), G4Element::GetNumberOfIsotopes(), G4Material::GetPressure(), G4Material::GetState(), G4Material::GetTemperature(), G4Material::GetZ(), kStateGas, kStateLiquid, kStateSolid, materialsElement, MEEWrite(), G4GDMLWrite::NewAttribute(), G4GDMLWrite::NewElement(), PropertyWrite(), PWrite(), and TWrite().

Referenced by AddMaterial().

00148 {
00149    G4String state_str("undefined");
00150    const G4State state = materialPtr->GetState();
00151    if (state==kStateSolid) { state_str = "solid"; } else
00152    if (state==kStateLiquid) { state_str = "liquid"; } else
00153    if (state==kStateGas) { state_str = "gas"; }
00154 
00155    const G4String name = GenerateName(materialPtr->GetName(), materialPtr);
00156 
00157    xercesc::DOMElement* materialElement = NewElement("material");
00158    materialElement->setAttributeNode(NewAttribute("name",name));
00159    materialElement->setAttributeNode(NewAttribute("state",state_str));
00160 
00161    // Write any property attached to the material...
00162    //
00163    if (materialPtr->GetMaterialPropertiesTable())
00164    {
00165      PropertyWrite(materialElement, materialPtr);
00166    }
00167 
00168    if (materialPtr->GetTemperature() != STP_Temperature)
00169      { TWrite(materialElement,materialPtr->GetTemperature()); }
00170    if (materialPtr->GetPressure() != STP_Pressure)
00171      { PWrite(materialElement,materialPtr->GetPressure()); }
00172 
00173    // Write Ionisation potential (mean excitation energy)
00174    MEEWrite(materialElement,materialPtr->GetIonisation()->GetMeanExcitationEnergy());
00175    
00176    DWrite(materialElement,materialPtr->GetDensity());
00177   
00178    const size_t NumberOfElements = materialPtr->GetNumberOfElements();
00179 
00180    if ( (NumberOfElements>1)
00181       || ( materialPtr->GetElement(0)
00182         && materialPtr->GetElement(0)->GetNumberOfIsotopes()>1 ) )
00183    {
00184       const G4double* MassFractionVector = materialPtr->GetFractionVector();
00185 
00186       for (size_t i=0;i<NumberOfElements;i++)
00187       {
00188          const G4String fractionref =
00189                         GenerateName(materialPtr->GetElement(i)->GetName(),
00190                                      materialPtr->GetElement(i));
00191          xercesc::DOMElement* fractionElement = NewElement("fraction");
00192          fractionElement->setAttributeNode(NewAttribute("n",
00193                                            MassFractionVector[i]));
00194          fractionElement->setAttributeNode(NewAttribute("ref",fractionref));
00195          materialElement->appendChild(fractionElement);
00196          AddElement(materialPtr->GetElement(i));
00197       }
00198    }
00199    else
00200    {
00201       materialElement->setAttributeNode(NewAttribute("Z",materialPtr->GetZ()));
00202       AtomWrite(materialElement,materialPtr->GetA());
00203    }
00204 
00205    // Append the material AFTER all the possible components are appended!
00206    //
00207    materialsElement->appendChild(materialElement);
00208 }

void G4GDMLWriteMaterials::MEEWrite ( xercesc::DOMElement *  ,
const G4double  
) [protected]

Definition at line 91 of file G4GDMLWriteMaterials.cc.

References G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by MaterialWrite().

00092 {
00093    xercesc::DOMElement* PElement = NewElement("MEE");
00094    PElement->setAttributeNode(NewAttribute("unit","eV"));
00095    PElement->setAttributeNode(NewAttribute("value",MEE/electronvolt));
00096    element->appendChild(PElement);
00097 }

void G4GDMLWriteMaterials::PropertyVectorWrite ( const G4String ,
const G4PhysicsOrderedFreeVector const 
) [protected]

Definition at line 210 of file G4GDMLWriteMaterials.cc.

References G4GDMLWriteDefine::defineElement, G4PhysicsVector::Energy(), G4GDMLWrite::GenerateName(), G4PhysicsVector::GetVectorLength(), G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by PropertyWrite().

00212 {
00213    const G4String matrixref = GenerateName(key, pvec);
00214    xercesc::DOMElement* matrixElement = NewElement("matrix");
00215    matrixElement->setAttributeNode(NewAttribute("name", matrixref));
00216    matrixElement->setAttributeNode(NewAttribute("coldim", "2"));
00217    std::ostringstream pvalues;
00218    for (size_t i=0; i<pvec->GetVectorLength(); i++)
00219    {
00220        if (i!=0)  { pvalues << " "; }
00221        pvalues << pvec->Energy(i) << " " << (*pvec)[i];
00222    }
00223    matrixElement->setAttributeNode(NewAttribute("values", pvalues.str()));
00224 
00225    defineElement->appendChild(matrixElement);
00226 }

void G4GDMLWriteMaterials::PropertyWrite ( xercesc::DOMElement *  ,
const G4Material const 
) [protected]

Definition at line 228 of file G4GDMLWriteMaterials.cc.

References G4GDMLWriteDefine::defineElement, G4Exception(), G4GDMLWrite::GenerateName(), G4Material::GetMaterialPropertiesTable(), G4Material::GetName(), G4MaterialPropertiesTable::GetPropertiesCMap(), G4MaterialPropertiesTable::GetPropertiesMap(), JustWarning, G4GDMLWrite::NewAttribute(), G4GDMLWrite::NewElement(), and PropertyVectorWrite().

Referenced by MaterialWrite().

00230 {
00231    xercesc::DOMElement* propElement;
00232    G4MaterialPropertiesTable* ptable = mat->GetMaterialPropertiesTable();
00233    const std::map< G4String, G4PhysicsOrderedFreeVector*,
00234                  std::less<G4String> >* pmap = ptable->GetPropertiesMap();
00235    const std::map< G4String, G4double,
00236                  std::less<G4String> >* cmap = ptable->GetPropertiesCMap();
00237    std::map< G4String, G4PhysicsOrderedFreeVector*,
00238                  std::less<G4String> >::const_iterator mpos;
00239    std::map< G4String, G4double,
00240                  std::less<G4String> >::const_iterator cpos;
00241    for (mpos=pmap->begin(); mpos!=pmap->end(); mpos++)
00242    {
00243       propElement = NewElement("property");
00244       propElement->setAttributeNode(NewAttribute("name", mpos->first));
00245       propElement->setAttributeNode(NewAttribute("ref",
00246                                     GenerateName(mpos->first, mpos->second)));
00247       if (mpos->second)
00248       {
00249          PropertyVectorWrite(mpos->first, mpos->second);
00250          matElement->appendChild(propElement);
00251       }
00252       else
00253       {
00254          G4String warn_message = "Null pointer for material property -"
00255                   + mpos->first + "- of material -" + mat->GetName() + "- !";
00256          G4Exception("G4GDMLWriteMaterials::PropertyWrite()", "NullPointer",
00257                      JustWarning, warn_message);
00258          continue;
00259       }
00260    }
00261    for (cpos=cmap->begin(); cpos!=cmap->end(); cpos++)
00262    {
00263       propElement = NewElement("property");
00264       propElement->setAttributeNode(NewAttribute("name", cpos->first));
00265       propElement->setAttributeNode(NewAttribute("ref", cpos->first));
00266       xercesc::DOMElement* constElement = NewElement("constant");
00267       constElement->setAttributeNode(NewAttribute("name", cpos->first));
00268       constElement->setAttributeNode(NewAttribute("value", cpos->second));
00269       defineElement->appendChild(constElement);
00270       matElement->appendChild(propElement);
00271    }
00272 }

void G4GDMLWriteMaterials::PWrite ( xercesc::DOMElement *  ,
const G4double  
) [protected]

Definition at line 73 of file G4GDMLWriteMaterials.cc.

References G4GDMLWrite::NewAttribute(), G4GDMLWrite::NewElement(), and pascal.

Referenced by MaterialWrite().

00074 {
00075    xercesc::DOMElement* PElement = NewElement("P");
00076    PElement->setAttributeNode(NewAttribute("unit","pascal"));
00077    PElement->setAttributeNode(NewAttribute("value",P/pascal));
00078    element->appendChild(PElement);
00079 }

void G4GDMLWriteMaterials::TWrite ( xercesc::DOMElement *  ,
const G4double  
) [protected]

Definition at line 82 of file G4GDMLWriteMaterials.cc.

References G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by MaterialWrite().

00083 {
00084    xercesc::DOMElement* TElement = NewElement("T");
00085    TElement->setAttributeNode(NewAttribute("unit","K"));
00086    TElement->setAttributeNode(NewAttribute("value",T/kelvin));
00087    element->appendChild(TElement);
00088 }


Field Documentation

std::vector<const G4Element*> G4GDMLWriteMaterials::elementList [protected]

Definition at line 83 of file G4GDMLWriteMaterials.hh.

Referenced by AddElement(), and MaterialsWrite().

std::vector<const G4Isotope*> G4GDMLWriteMaterials::isotopeList [protected]

Definition at line 82 of file G4GDMLWriteMaterials.hh.

Referenced by AddIsotope(), and MaterialsWrite().

std::vector<const G4Material*> G4GDMLWriteMaterials::materialList [protected]

Definition at line 84 of file G4GDMLWriteMaterials.hh.

Referenced by AddMaterial(), and MaterialsWrite().

xercesc::DOMElement* G4GDMLWriteMaterials::materialsElement [protected]

Definition at line 85 of file G4GDMLWriteMaterials.hh.

Referenced by ElementWrite(), IsotopeWrite(), MaterialsWrite(), and MaterialWrite().


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