G4GDMLReadMaterials Class Reference

#include <G4GDMLReadMaterials.hh>

Inheritance diagram for G4GDMLReadMaterials:

G4GDMLReadDefine G4GDMLRead G4GDMLReadSolids G4GDMLReadSetup G4GDMLReadParamvol G4GDMLReadStructure

Public Member Functions

G4ElementGetElement (const G4String &, G4bool verbose=true) const
G4IsotopeGetIsotope (const G4String &, G4bool verbose=true) const
G4MaterialGetMaterial (const G4String &, G4bool verbose=true) const
virtual void MaterialsRead (const xercesc::DOMElement *const)

Protected Member Functions

 G4GDMLReadMaterials ()
virtual ~G4GDMLReadMaterials ()
G4double AtomRead (const xercesc::DOMElement *const)
G4int CompositeRead (const xercesc::DOMElement *const, G4String &)
G4double DRead (const xercesc::DOMElement *const)
G4double PRead (const xercesc::DOMElement *const)
G4double TRead (const xercesc::DOMElement *const)
G4double MEERead (const xercesc::DOMElement *const)
void ElementRead (const xercesc::DOMElement *const)
G4double FractionRead (const xercesc::DOMElement *const, G4String &)
void IsotopeRead (const xercesc::DOMElement *const)
void MaterialRead (const xercesc::DOMElement *const)
void MixtureRead (const xercesc::DOMElement *const, G4Element *)
void MixtureRead (const xercesc::DOMElement *const, G4Material *)
void PropertyRead (const xercesc::DOMElement *const, G4Material *)

Detailed Description

Definition at line 51 of file G4GDMLReadMaterials.hh.


Constructor & Destructor Documentation

G4GDMLReadMaterials::G4GDMLReadMaterials (  )  [protected]

Definition at line 45 of file G4GDMLReadMaterials.cc.

00045                                          : G4GDMLReadDefine()
00046 {
00047 }

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

Definition at line 49 of file G4GDMLReadMaterials.cc.

00050 {
00051 }


Member Function Documentation

G4double G4GDMLReadMaterials::AtomRead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 54 of file G4GDMLReadMaterials.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().

Referenced by ElementRead(), IsotopeRead(), and MaterialRead().

00055 {
00056    G4double value = 0.0;
00057    G4double unit = g/mole;
00058 
00059    const xercesc::DOMNamedNodeMap* const attributes
00060          = atomElement->getAttributes();
00061    XMLSize_t attributeCount = attributes->getLength();
00062 
00063    for (XMLSize_t attribute_index=0;
00064         attribute_index<attributeCount; attribute_index++)
00065    {
00066       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00067 
00068       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00069       { continue; }
00070 
00071       const xercesc::DOMAttr* const attribute
00072             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00073       if (!attribute)
00074       {
00075         G4Exception("G4GDMLReadMaterials::AtomRead()", "InvalidRead",
00076                     FatalException, "No attribute found!");
00077         return value;
00078       }
00079       const G4String attName = Transcode(attribute->getName());
00080       const G4String attValue = Transcode(attribute->getValue());
00081 
00082       if (attName=="value") { value = eval.Evaluate(attValue); } else
00083       if (attName=="unit")  { unit = eval.Evaluate(attValue); }
00084    }
00085 
00086    return value*unit;
00087 }

G4int G4GDMLReadMaterials::CompositeRead ( const xercesc::DOMElement *  const,
G4String  
) [protected]

Definition at line 90 of file G4GDMLReadMaterials.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::EvaluateInteger(), FatalException, G4Exception(), CLHEP::detail::n, and G4GDMLRead::Transcode().

Referenced by MixtureRead().

00091 {
00092    G4int n = 0;
00093 
00094    const xercesc::DOMNamedNodeMap* const attributes
00095          = compositeElement->getAttributes();
00096    XMLSize_t attributeCount = attributes->getLength();
00097 
00098    for (XMLSize_t attribute_index=0;
00099         attribute_index<attributeCount; attribute_index++)
00100    {
00101       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00102 
00103       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00104       { continue; }
00105 
00106       const xercesc::DOMAttr* const attribute
00107             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00108       if (!attribute)
00109       {
00110         G4Exception("G4GDMLReadMaterials::CompositeRead()", "InvalidRead",
00111                     FatalException, "No attribute found!");
00112         return n;
00113       }
00114       const G4String attName = Transcode(attribute->getName());
00115       const G4String attValue = Transcode(attribute->getValue());
00116 
00117       if (attName=="n")  { n = eval.EvaluateInteger(attValue); } else
00118       if (attName=="ref") { ref = attValue; }
00119    }
00120 
00121    return n;
00122 }

G4double G4GDMLReadMaterials::DRead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 124 of file G4GDMLReadMaterials.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().

Referenced by MaterialRead().

00125 {
00126    G4double value = 0.0;
00127    G4double unit = g/cm3;
00128 
00129    const xercesc::DOMNamedNodeMap* const attributes
00130          = DElement->getAttributes();
00131    XMLSize_t attributeCount = attributes->getLength();
00132 
00133    for (XMLSize_t attribute_index=0;
00134         attribute_index<attributeCount; attribute_index++)
00135    {
00136       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00137 
00138       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00139       { continue; }
00140 
00141       const xercesc::DOMAttr* const attribute
00142             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00143       if (!attribute)
00144       {
00145         G4Exception("G4GDMLReadMaterials::DRead()", "InvalidRead",
00146                     FatalException, "No attribute found!");
00147         return value;
00148       }
00149       const G4String attName = Transcode(attribute->getName());
00150       const G4String attValue = Transcode(attribute->getValue());
00151 
00152       if (attName=="value") { value = eval.Evaluate(attValue); } else
00153       if (attName=="unit")  { unit = eval.Evaluate(attValue); }
00154    }
00155 
00156    return value*unit;
00157 }

void G4GDMLReadMaterials::ElementRead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 262 of file G4GDMLReadMaterials.cc.

References AtomRead(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), MixtureRead(), G4GDMLRead::Strip(), and G4GDMLRead::Transcode().

Referenced by MaterialsRead().

00263 {
00264    G4String name;
00265    G4String formula;
00266    G4double a = 0.0;
00267    G4double Z = 0.0;
00268 
00269    const xercesc::DOMNamedNodeMap* const attributes
00270          = elementElement->getAttributes();
00271    XMLSize_t attributeCount = attributes->getLength();
00272 
00273    for (XMLSize_t attribute_index=0;
00274         attribute_index<attributeCount; attribute_index++)
00275    {
00276       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00277 
00278       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00279       { continue; }
00280 
00281       const xercesc::DOMAttr* const attribute
00282             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00283       if (!attribute)
00284       {
00285         G4Exception("G4GDMLReadMaterials::ElementRead()", "InvalidRead",
00286                     FatalException, "No attribute found!");
00287         return;
00288       }
00289       const G4String attName = Transcode(attribute->getName());
00290       const G4String attValue = Transcode(attribute->getValue());
00291 
00292       if (attName=="name") { name = GenerateName(attValue); } else
00293       if (attName=="formula") { formula = attValue; } else
00294       if (attName=="Z") { Z = eval.Evaluate(attValue); }
00295    }
00296 
00297    G4int nComponents = 0;
00298 
00299    for (xercesc::DOMNode* iter = elementElement->getFirstChild();
00300         iter != 0; iter = iter->getNextSibling())
00301    {
00302       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00303 
00304       const xercesc::DOMElement* const child
00305             = dynamic_cast<xercesc::DOMElement*>(iter);
00306       if (!child)
00307       {
00308         G4Exception("G4GDMLReadMaterials::ElementRead()", "InvalidRead",
00309                     FatalException, "No child found!");
00310         return;
00311       }
00312       const G4String tag = Transcode(child->getTagName());
00313 
00314       if (tag=="atom") { a = AtomRead(child); }  else
00315       if (tag=="fraction") { nComponents++; }
00316    }
00317 
00318    if (nComponents>0)
00319    {
00320      MixtureRead(elementElement,
00321                  new G4Element(Strip(name),formula,nComponents));
00322    }
00323    else
00324    {
00325      new G4Element(Strip(name),formula,Z,a);
00326    }
00327 }

G4double G4GDMLReadMaterials::FractionRead ( const xercesc::DOMElement *  const,
G4String  
) [protected]

Definition at line 330 of file G4GDMLReadMaterials.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), CLHEP::detail::n, and G4GDMLRead::Transcode().

Referenced by MixtureRead().

00331 {
00332    G4double n = 0.0;
00333 
00334    const xercesc::DOMNamedNodeMap* const attributes
00335          = fractionElement->getAttributes();
00336    XMLSize_t attributeCount = attributes->getLength();
00337 
00338    for (XMLSize_t attribute_index=0;
00339         attribute_index<attributeCount; attribute_index++)
00340    {
00341       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00342 
00343       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00344       { continue; }
00345 
00346       const xercesc::DOMAttr* const attribute
00347             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00348       if (!attribute)
00349       {
00350         G4Exception("G4GDMLReadMaterials::FractionRead()", "InvalidRead",
00351                     FatalException, "No attribute found!");
00352         return n;
00353       }
00354       const G4String attName = Transcode(attribute->getName());
00355       const G4String attValue = Transcode(attribute->getValue());
00356 
00357       if (attName=="n")   { n = eval.Evaluate(attValue); } else
00358       if (attName=="ref") { ref = attValue; }
00359    }
00360 
00361    return n;
00362 }

G4Element * G4GDMLReadMaterials::GetElement ( const G4String ,
G4bool  verbose = true 
) const

Definition at line 709 of file G4GDMLReadMaterials.cc.

References FatalException, G4NistManager::FindOrBuildElement(), G4Exception(), G4Element::GetElement(), and G4NistManager::Instance().

Referenced by MixtureRead().

00710 {
00711    G4Element* elementPtr = G4Element::GetElement(ref,false);
00712 
00713    if (!elementPtr)
00714    {
00715      elementPtr = G4NistManager::Instance()->FindOrBuildElement(ref);
00716    }
00717 
00718    if (verbose && !elementPtr)
00719    {
00720      G4String error_msg = "Referenced element '" + ref + "' was not found!";
00721      G4Exception("G4GDMLReadMaterials::GetElement()", "InvalidRead",
00722                  FatalException, error_msg);
00723    }
00724 
00725    return elementPtr;
00726 }

G4Isotope * G4GDMLReadMaterials::GetIsotope ( const G4String ,
G4bool  verbose = true 
) const

Definition at line 728 of file G4GDMLReadMaterials.cc.

References FatalException, G4Exception(), and G4Isotope::GetIsotope().

Referenced by MixtureRead().

00730 {
00731    G4Isotope* isotopePtr = G4Isotope::GetIsotope(ref,false);
00732 
00733    if (verbose && !isotopePtr)
00734    {
00735      G4String error_msg = "Referenced isotope '" + ref + "' was not found!";
00736      G4Exception("G4GDMLReadMaterials::GetIsotope()", "InvalidRead",
00737                  FatalException, error_msg);
00738    }
00739 
00740    return isotopePtr;
00741 }

G4Material * G4GDMLReadMaterials::GetMaterial ( const G4String ,
G4bool  verbose = true 
) const

Definition at line 743 of file G4GDMLReadMaterials.cc.

References FatalException, G4NistManager::FindOrBuildMaterial(), G4Exception(), G4Material::GetMaterial(), and G4NistManager::Instance().

Referenced by MixtureRead(), and G4GDMLReadStructure::VolumeRead().

00745 {
00746    G4Material *materialPtr = G4Material::GetMaterial(ref,false);
00747 
00748    if (!materialPtr)
00749    {
00750      materialPtr = G4NistManager::Instance()->FindOrBuildMaterial(ref);
00751    }
00752 
00753    if (verbose && !materialPtr)
00754    {
00755      G4String error_msg = "Referenced material '" + ref + "' was not found!";
00756      G4Exception("G4GDMLReadMaterials::GetMaterial()", "InvalidRead",
00757                  FatalException, error_msg);
00758    }
00759 
00760    return materialPtr;
00761 }

void G4GDMLReadMaterials::IsotopeRead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 365 of file G4GDMLReadMaterials.cc.

References AtomRead(), G4GDMLRead::eval, G4GDMLEvaluator::EvaluateInteger(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4GDMLRead::Strip(), and G4GDMLRead::Transcode().

Referenced by MaterialsRead().

00366 {
00367    G4String name;
00368    G4int Z = 0;
00369    G4int N = 0;
00370    G4double a = 0.0;
00371 
00372    const xercesc::DOMNamedNodeMap* const attributes
00373          = isotopeElement->getAttributes();
00374    XMLSize_t attributeCount = attributes->getLength();
00375 
00376    for (XMLSize_t attribute_index=0;
00377         attribute_index<attributeCount;attribute_index++)
00378    {
00379       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00380 
00381       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00382       { continue; }
00383 
00384       const xercesc::DOMAttr* const attribute
00385             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00386       if (!attribute)
00387       {
00388         G4Exception("G4GDMLReadMaterials::IsotopeRead()", "InvalidRead",
00389                     FatalException, "No attribute found!");
00390         return;
00391       }
00392       const G4String attName = Transcode(attribute->getName());
00393       const G4String attValue = Transcode(attribute->getValue());
00394 
00395       if (attName=="name") { name = GenerateName(attValue); } else
00396       if (attName=="Z") { Z = eval.EvaluateInteger(attValue); } else
00397       if (attName=="N") { N = eval.EvaluateInteger(attValue); }
00398    }
00399 
00400    for (xercesc::DOMNode* iter = isotopeElement->getFirstChild();
00401         iter != 0; iter = iter->getNextSibling())
00402    {
00403       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00404 
00405       const xercesc::DOMElement* const child
00406             = dynamic_cast<xercesc::DOMElement*>(iter);
00407       if (!child)
00408       {
00409         G4Exception("G4GDMLReadMaterials::IsotopeRead()", "InvalidRead",
00410                     FatalException, "No child found!");
00411         return;
00412       }
00413       const G4String tag = Transcode(child->getTagName());
00414 
00415       if (tag=="atom")  { a = AtomRead(child); }
00416    }
00417 
00418    new G4Isotope(Strip(name),Z,N,a);
00419 }

void G4GDMLReadMaterials::MaterialRead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 422 of file G4GDMLReadMaterials.cc.

References AtomRead(), DRead(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4Material::GetIonisation(), G4GDMLReadDefine::GetQuantity(), kStateGas, kStateLiquid, kStateSolid, kStateUndefined, MEERead(), MixtureRead(), PRead(), PropertyRead(), G4GDMLReadDefine::RefRead(), G4IonisParamMat::SetMeanExcitationEnergy(), G4GDMLRead::Strip(), G4GDMLRead::Transcode(), and TRead().

Referenced by MaterialsRead().

00423 {
00424    G4String name;
00425    G4double Z = 0.0;
00426    G4double a = 0.0;
00427    G4double D = 0.0;
00428    G4State state = kStateUndefined;
00429    G4double T = STP_Temperature;
00430    G4double P = STP_Pressure;
00431    G4double MEE = -1.0;
00432 
00433    const xercesc::DOMNamedNodeMap* const attributes
00434          = materialElement->getAttributes();
00435    XMLSize_t attributeCount = attributes->getLength();
00436 
00437    for (XMLSize_t attribute_index=0;
00438         attribute_index<attributeCount; attribute_index++)
00439    {
00440       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00441 
00442       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00443       { continue; }
00444 
00445       const xercesc::DOMAttr* const attribute
00446             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00447       if (!attribute)
00448       {
00449         G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
00450                     FatalException, "No attribute found!");
00451         return;
00452       }
00453       const G4String attName = Transcode(attribute->getName());
00454       const G4String attValue = Transcode(attribute->getValue());
00455 
00456       if (attName=="name") { name = GenerateName(attValue); } else
00457       if (attName=="Z") { Z = eval.Evaluate(attValue); } else
00458       if (attName=="state")
00459       {
00460          if (attValue=="solid")  { state = kStateSolid;  } else
00461          if (attValue=="liquid") { state = kStateLiquid; } else
00462          if (attValue=="gas")    { state = kStateGas; }
00463       }
00464    }
00465 
00466    size_t nComponents = 0;
00467 
00468    for (xercesc::DOMNode* iter = materialElement->getFirstChild();
00469         iter != 0; iter = iter->getNextSibling())
00470    {
00471       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00472 
00473       const xercesc::DOMElement* const child
00474             = dynamic_cast<xercesc::DOMElement*>(iter);
00475       if (!child)
00476       {
00477         G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
00478                     FatalException, "No child found!");
00479         return;
00480       }
00481       const G4String tag = Transcode(child->getTagName());
00482 
00483       if (tag=="atom") { a = AtomRead(child); } else
00484       if (tag=="Dref") { D = GetQuantity(GenerateName(RefRead(child))); } else
00485       if (tag=="Pref") { P = GetQuantity(GenerateName(RefRead(child))); } else
00486       if (tag=="Tref") { T = GetQuantity(GenerateName(RefRead(child))); } else
00487       if (tag=="MEEref") { MEE = GetQuantity(GenerateName(RefRead(child))); } else
00488       if (tag=="D") { D = DRead(child); } else
00489       if (tag=="P") { P = PRead(child); } else
00490       if (tag=="T") { T = TRead(child); } else
00491       if (tag=="MEE") { MEE = MEERead(child); } else
00492       if (tag=="fraction" || tag=="composite")  { nComponents++; }
00493    }
00494 
00495    G4Material* material =  0;
00496 
00497    if (nComponents==0)
00498    {
00499      material = new G4Material(Strip(name),Z,a,D,state,T,P);
00500    }
00501    else
00502    {
00503      material = new G4Material(Strip(name),D,nComponents,state,T,P);
00504      MixtureRead(materialElement, material);
00505    }
00506    if (MEE != -1)  // ionisation potential (mean excitation energy)
00507    {
00508      material->GetIonisation()->SetMeanExcitationEnergy(MEE);
00509    }
00510 
00511    for (xercesc::DOMNode* iter = materialElement->getFirstChild();
00512         iter != 0; iter = iter->getNextSibling())
00513    {
00514       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00515 
00516       const xercesc::DOMElement* const child
00517             = dynamic_cast<xercesc::DOMElement*>(iter);
00518       if (!child)
00519       {
00520         G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
00521                     FatalException, "No child found!");
00522         return;
00523       }
00524       const G4String tag = Transcode(child->getTagName());
00525 
00526       if (tag=="property") { PropertyRead(child,material); }
00527    }
00528 }

void G4GDMLReadMaterials::MaterialsRead ( const xercesc::DOMElement *  const  )  [virtual]

Implements G4GDMLRead.

Definition at line 676 of file G4GDMLReadMaterials.cc.

References G4GDMLReadDefine::DefineRead(), ElementRead(), FatalException, G4cout, G4endl, G4Exception(), IsotopeRead(), MaterialRead(), and G4GDMLRead::Transcode().

00677 {
00678    G4cout << "G4GDML: Reading materials..." << G4endl;
00679 
00680    for (xercesc::DOMNode* iter = materialsElement->getFirstChild();
00681         iter != 0; iter = iter->getNextSibling())
00682    {
00683       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00684 
00685       const xercesc::DOMElement* const child
00686             = dynamic_cast<xercesc::DOMElement*>(iter);
00687       if (!child)
00688       {
00689         G4Exception("G4GDMLReadMaterials::MaterialsRead()", "InvalidRead",
00690                     FatalException, "No child found!");
00691         return;
00692       }
00693       const G4String tag = Transcode(child->getTagName());
00694       
00695       if (tag=="define")   { DefineRead(child);  }  else 
00696       if (tag=="element")  { ElementRead(child); }  else 
00697       if (tag=="isotope")  { IsotopeRead(child); }  else 
00698       if (tag=="material") { MaterialRead(child); }
00699       else
00700       {
00701         G4String error_msg = "Unknown tag in materials: " + tag;
00702         G4Exception("G4GDMLReadMaterials::MaterialsRead()", "InvalidSetup",
00703                     FatalException, error_msg);
00704       }
00705    }
00706 }

G4double G4GDMLReadMaterials::MEERead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 227 of file G4GDMLReadMaterials.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().

Referenced by MaterialRead().

00228 {
00229    G4double value = -1;
00230    G4double unit = eV;
00231 
00232    const xercesc::DOMNamedNodeMap* const attributes = PElement->getAttributes();
00233    XMLSize_t attributeCount = attributes->getLength();
00234 
00235    for (XMLSize_t attribute_index=0;
00236         attribute_index<attributeCount; attribute_index++)
00237    {
00238       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00239 
00240       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00241       { continue; }
00242 
00243       const xercesc::DOMAttr* const attribute
00244             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
00245       if (!attribute)
00246       {
00247         G4Exception("G4GDMLReadMaterials::MEERead()", "InvalidRead",
00248                     FatalException, "No attribute found!");
00249         return value;
00250       }
00251       const G4String attName = Transcode(attribute->getName());
00252       const G4String attValue = Transcode(attribute->getValue());
00253 
00254       if (attName=="value") { value = eval.Evaluate(attValue); } else
00255       if (attName=="unit")  { unit = eval.Evaluate(attValue); }
00256    }
00257 
00258    return value*unit;
00259 }

void G4GDMLReadMaterials::MixtureRead ( const xercesc::DOMElement *  const,
G4Material  
) [protected]

Definition at line 558 of file G4GDMLReadMaterials.cc.

References G4Material::AddElement(), G4Material::AddMaterial(), CompositeRead(), FatalException, FractionRead(), G4Exception(), G4GDMLRead::GenerateName(), GetElement(), GetMaterial(), CLHEP::detail::n, and G4GDMLRead::Transcode().

00560 {
00561    for (xercesc::DOMNode* iter = mixtureElement->getFirstChild();
00562         iter != 0; iter = iter->getNextSibling())
00563    {
00564       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00565 
00566       const xercesc::DOMElement* const child
00567             = dynamic_cast<xercesc::DOMElement*>(iter);
00568       if (!child)
00569       {
00570         G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidRead",
00571                     FatalException, "No child found!");
00572         return;
00573       }
00574       const G4String tag = Transcode(child->getTagName());
00575 
00576       if (tag=="fraction")
00577       {
00578          G4String ref;
00579          G4double n = FractionRead(child,ref);
00580          
00581          G4Material *materialPtr = GetMaterial(GenerateName(ref,true), false);
00582          G4Element *elementPtr = GetElement(GenerateName(ref,true), false);
00583 
00584          if (materialPtr != 0) { material->AddMaterial(materialPtr,n); } else
00585          if (elementPtr != 0)  { material->AddElement(elementPtr,n); }
00586 
00587          if ((materialPtr == 0) && (elementPtr == 0))
00588          {
00589             G4String error_msg = "Referenced material/element '"
00590                                + GenerateName(ref,true) + "' was not found!";
00591             G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidSetup",
00592                         FatalException, error_msg);   
00593          }
00594       } 
00595       else if (tag=="composite")
00596       {
00597          G4String ref;
00598          G4int n = CompositeRead(child,ref);
00599 
00600          G4Element *elementPtr = GetElement(GenerateName(ref,true));
00601          material->AddElement(elementPtr,n);
00602       }
00603    }
00604 }

void G4GDMLReadMaterials::MixtureRead ( const xercesc::DOMElement *  const,
G4Element  
) [protected]

Definition at line 531 of file G4GDMLReadMaterials.cc.

References G4Element::AddIsotope(), FatalException, FractionRead(), G4Exception(), G4GDMLRead::GenerateName(), GetIsotope(), CLHEP::detail::n, and G4GDMLRead::Transcode().

Referenced by ElementRead(), and MaterialRead().

00532 {
00533    for (xercesc::DOMNode* iter = mixtureElement->getFirstChild();
00534         iter != 0; iter = iter->getNextSibling())
00535    {
00536       if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
00537 
00538       const xercesc::DOMElement* const child
00539             = dynamic_cast<xercesc::DOMElement*>(iter);
00540       if (!child)
00541       {
00542         G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidRead",
00543                     FatalException, "No child found!");
00544         return;
00545       }
00546       const G4String tag = Transcode(child->getTagName());
00547 
00548       if (tag=="fraction")
00549       {
00550          G4String ref;
00551          G4double n = FractionRead(child,ref);
00552          element->AddIsotope(GetIsotope(GenerateName(ref,true)),n);
00553       }
00554    }
00555 }

G4double G4GDMLReadMaterials::PRead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 159 of file G4GDMLReadMaterials.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), pascal, and G4GDMLRead::Transcode().

Referenced by MaterialRead().

00160 {
00161    G4double value = STP_Pressure;
00162    G4double unit = pascal;
00163 
00164    const xercesc::DOMNamedNodeMap* const attributes = PElement->getAttributes();
00165    XMLSize_t attributeCount = attributes->getLength();
00166 
00167    for (XMLSize_t attribute_index=0;
00168         attribute_index<attributeCount; attribute_index++)
00169    {
00170       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00171 
00172       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00173       { continue; }
00174 
00175       const xercesc::DOMAttr* const attribute
00176             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00177       if (!attribute)
00178       {
00179         G4Exception("G4GDMLReadMaterials::PRead()", "InvalidRead",
00180                     FatalException, "No attribute found!");
00181         return value;
00182       }
00183       const G4String attName = Transcode(attribute->getName());
00184       const G4String attValue = Transcode(attribute->getValue());
00185 
00186       if (attName=="value") { value = eval.Evaluate(attValue); } else
00187       if (attName=="unit")  { unit = eval.Evaluate(attValue); }
00188    }
00189 
00190    return value*unit;
00191 }

void G4GDMLReadMaterials::PropertyRead ( const xercesc::DOMElement *  const,
G4Material  
) [protected]

Definition at line 607 of file G4GDMLReadMaterials.cc.

References G4MaterialPropertiesTable::AddConstProperty(), G4MaterialPropertiesTable::AddProperty(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4GDMLMatrix::Get(), G4GDMLMatrix::GetCols(), G4Material::GetMaterialPropertiesTable(), G4GDMLReadDefine::GetMatrix(), G4GDMLMatrix::GetRows(), G4PhysicsOrderedFreeVector::InsertValues(), G4Material::SetMaterialPropertiesTable(), G4GDMLRead::Strip(), and G4GDMLRead::Transcode().

Referenced by MaterialRead().

00609 {
00610    G4String name;
00611    G4String ref;
00612    G4GDMLMatrix matrix;
00613 
00614    const xercesc::DOMNamedNodeMap* const attributes
00615          = propertyElement->getAttributes();
00616    XMLSize_t attributeCount = attributes->getLength();
00617 
00618    for (XMLSize_t attribute_index=0;
00619         attribute_index<attributeCount; attribute_index++)
00620    {
00621       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00622 
00623       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00624       { continue; }
00625 
00626       const xercesc::DOMAttr* const attribute
00627             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00628       if (!attribute)
00629       {
00630         G4Exception("G4GDMLReadMaterials::PropertyRead()", "InvalidRead",
00631                     FatalException, "No attribute found!");
00632         return;
00633       }
00634       const G4String attName = Transcode(attribute->getName());
00635       const G4String attValue = Transcode(attribute->getValue());
00636 
00637       if (attName=="name") { name = GenerateName(attValue); } else
00638       if (attName=="ref")  { matrix = GetMatrix(ref=attValue); }
00639    }
00640 
00641    /*
00642    if (matrix.GetCols() != 2)
00643    {
00644      G4String error_msg = "Referenced matrix '" + ref
00645             + "' should have \n two columns as a property table for material: "
00646             + material->GetName();
00647      G4Exception("G4GDMLReadMaterials::PropertyRead()", "InvalidRead",
00648                  FatalException, error_msg);
00649    }
00650    */
00651 
00652    if (matrix.GetRows() == 0) { return; }
00653 
00654    G4MaterialPropertiesTable* matprop=material->GetMaterialPropertiesTable();
00655    if (!matprop)
00656    {
00657      matprop = new G4MaterialPropertiesTable();
00658      material->SetMaterialPropertiesTable(matprop);
00659    }
00660    if (matrix.GetCols() == 1)  // constant property assumed
00661    {
00662      matprop->AddConstProperty(Strip(name), matrix.Get(0,0));
00663    }
00664    else  // build the material properties vector
00665    {
00666      G4MaterialPropertyVector* propvect = new G4MaterialPropertyVector();
00667      for (size_t i=0; i<matrix.GetRows(); i++)
00668      {
00669        propvect->InsertValues(matrix.Get(i,0),matrix.Get(i,1));
00670      }
00671      matprop->AddProperty(Strip(name),propvect);
00672    }
00673 }

G4double G4GDMLReadMaterials::TRead ( const xercesc::DOMElement *  const  )  [protected]

Definition at line 193 of file G4GDMLReadMaterials.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), and G4GDMLRead::Transcode().

Referenced by MaterialRead().

00194 {
00195    G4double value = STP_Temperature;
00196    G4double unit = kelvin;
00197 
00198    const xercesc::DOMNamedNodeMap* const attributes = TElement->getAttributes();
00199    XMLSize_t attributeCount = attributes->getLength();
00200 
00201    for (XMLSize_t attribute_index=0;
00202         attribute_index<attributeCount; attribute_index++)
00203    {
00204       xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
00205 
00206       if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
00207       { continue; }
00208 
00209       const xercesc::DOMAttr* const attribute
00210             = dynamic_cast<xercesc::DOMAttr*>(attribute_node);   
00211       if (!attribute)
00212       {
00213         G4Exception("G4GDMLReadMaterials::TRead()", "InvalidRead",
00214                     FatalException, "No attribute found!");
00215         return value;
00216       }
00217       const G4String attName = Transcode(attribute->getName());
00218       const G4String attValue = Transcode(attribute->getValue());
00219 
00220       if (attName=="value") { value = eval.Evaluate(attValue); } else
00221       if (attName=="unit")  { unit = eval.Evaluate(attValue); }
00222    }
00223 
00224    return value*unit;
00225 }


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