Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GDMLWriteMaterials Class Reference

#include <G4GDMLWriteMaterials.hh>

Inheritance diagram for G4GDMLWriteMaterials:
Collaboration diagram for G4GDMLWriteMaterials:

Public Member Functions

void AddIsotope (const G4Isotope *const)
 
void AddElement (const G4Element *const)
 
void AddMaterial (const G4Material *const)
 
virtual void MaterialsWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWriteDefine
G4ThreeVector GetAngles (const G4RotationMatrix &)
 
void ScaleWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
 
void RotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void PositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void FirstrotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void FirstpositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void AddPosition (const G4String &name, const G4ThreeVector &pos)
 
virtual void DefineWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWrite
G4Transform3D Write (const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
 
void AddModule (const G4VPhysicalVolume *const topVol)
 
void AddModule (const G4int depth)
 
void AddAuxiliary (G4GDMLAuxStructType myaux)
 
virtual void SolidsWrite (xercesc::DOMElement *)=0
 
virtual void StructureWrite (xercesc::DOMElement *)=0
 
virtual G4Transform3D TraverseVolumeTree (const G4LogicalVolume *const, const G4int)=0
 
virtual void SurfacesWrite ()=0
 
virtual void SetupWrite (xercesc::DOMElement *, const G4LogicalVolume *const)=0
 
virtual void ExtensionWrite (xercesc::DOMElement *)
 
virtual void UserinfoWrite (xercesc::DOMElement *)
 
virtual void AddExtension (xercesc::DOMElement *, const G4LogicalVolume *const)
 
G4String GenerateName (const G4String &, const void *const)
 

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 Member Functions inherited from G4GDMLWriteDefine
 G4GDMLWriteDefine ()
 
virtual ~G4GDMLWriteDefine ()
 
void Scale_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Rotation_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Position_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
- Protected Member Functions inherited from G4GDMLWrite
 G4GDMLWrite ()
 
virtual ~G4GDMLWrite ()
 
VolumeMapType & VolumeMap ()
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4String &)
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4double &)
 
xercesc::DOMElement * NewElement (const G4String &)
 
G4String Modularize (const G4VPhysicalVolume *const topvol, const G4int depth)
 
void AddAuxInfo (G4GDMLAuxListType *auxInfoList, xercesc::DOMElement *element)
 
G4bool FileExists (const G4String &) const
 
PhysVolumeMapType & PvolumeMap ()
 
DepthMapType & DepthMap ()
 

Protected Attributes

std::vector< const G4Isotope * > isotopeList
 
std::vector< const G4Element * > elementList
 
std::vector< const G4Material * > materialList
 
xercesc::DOMElement * materialsElement
 
- Protected Attributes inherited from G4GDMLWriteDefine
xercesc::DOMElement * defineElement
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc
 
xercesc::DOMElement * extElement
 
xercesc::DOMElement * userinfoElement
 
XMLCh tempStr [10000]
 
G4GDMLAuxListType auxList
 

Additional Inherited Members

- Static Public Member Functions inherited from G4GDMLWrite
static void SetAddPointerToName (G4bool)
 
- Static Protected Attributes inherited from G4GDMLWriteDefine
static const G4double kRelativePrecision = DBL_EPSILON
 
static const G4double kAngularPrecision = DBL_EPSILON
 
static const G4double kLinearPrecision = DBL_EPSILON
 
- Static Protected Attributes inherited from G4GDMLWrite
static G4bool addPointerToName = true
 

Detailed Description

Definition at line 53 of file G4GDMLWriteMaterials.hh.

Constructor & Destructor Documentation

G4GDMLWriteMaterials::G4GDMLWriteMaterials ( )
protected

Definition at line 45 of file G4GDMLWriteMaterials.cc.

47 {
48 }
xercesc::DOMElement * materialsElement
G4GDMLWriteMaterials::~G4GDMLWriteMaterials ( )
protectedvirtual

Definition at line 50 of file G4GDMLWriteMaterials.cc.

51 {
52 }

Member Function Documentation

void G4GDMLWriteMaterials::AddElement ( const G4Element * const  elementPtr)

Definition at line 296 of file G4GDMLWriteMaterials.cc.

297 {
298  for (size_t i=0;i<elementList.size();i++) // Check if element is
299  { // already in the list!
300  if (elementList[i] == elementPtr) { return; }
301  }
302  elementList.push_back(elementPtr);
303  ElementWrite(elementPtr);
304 }
std::vector< const G4Element * > elementList
void ElementWrite(const G4Element *const)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::AddIsotope ( const G4Isotope * const  isotopePtr)

Definition at line 286 of file G4GDMLWriteMaterials.cc.

287 {
288  for (size_t i=0; i<isotopeList.size(); i++) // Check if isotope is
289  { // already in the list!
290  if (isotopeList[i] == isotopePtr) { return; }
291  }
292  isotopeList.push_back(isotopePtr);
293  IsotopeWrite(isotopePtr);
294 }
void IsotopeWrite(const G4Isotope *const)
std::vector< const G4Isotope * > isotopeList

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::AddMaterial ( const G4Material * const  materialPtr)

Definition at line 306 of file G4GDMLWriteMaterials.cc.

307 {
308  for (size_t i=0;i<materialList.size();i++) // Check if material is
309  { // already in the list!
310  if (materialList[i] == materialPtr) { return; }
311  }
312  materialList.push_back(materialPtr);
313  MaterialWrite(materialPtr);
314 }
std::vector< const G4Material * > materialList
void MaterialWrite(const G4Material *const)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::AtomWrite ( xercesc::DOMElement *  element,
const G4double a 
)
protected

Definition at line 55 of file G4GDMLWriteMaterials.cc.

56 {
57  xercesc::DOMElement* atomElement = NewElement("atom");
58  atomElement->setAttributeNode(NewAttribute("unit","g/mole"));
59  atomElement->setAttributeNode(NewAttribute("value",a*mole/g));
60  element->appendChild(atomElement);
61 }
static constexpr double g
Definition: G4SIunits.hh:183
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static constexpr double mole
Definition: G4SIunits.hh:286

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::DWrite ( xercesc::DOMElement *  element,
const G4double d 
)
protected

Definition at line 64 of file G4GDMLWriteMaterials.cc.

65 {
66  xercesc::DOMElement* DElement = NewElement("D");
67  DElement->setAttributeNode(NewAttribute("unit","g/cm3"));
68  DElement->setAttributeNode(NewAttribute("value",d*cm3/g));
69  element->appendChild(DElement);
70 }
static constexpr double g
Definition: G4SIunits.hh:183
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double cm3
Definition: G4SIunits.hh:121
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::ElementWrite ( const G4Element * const  elementPtr)
protected

Definition at line 112 of file G4GDMLWriteMaterials.cc.

113 {
114  const G4String name = GenerateName(elementPtr->GetName(),elementPtr);
115 
116  xercesc::DOMElement* elementElement = NewElement("element");
117  elementElement->setAttributeNode(NewAttribute("name",name));
118 
119  const size_t NumberOfIsotopes = elementPtr->GetNumberOfIsotopes();
120 
121  if (NumberOfIsotopes>0)
122  {
123  const G4double* RelativeAbundanceVector =
124  elementPtr->GetRelativeAbundanceVector();
125  for (size_t i=0;i<NumberOfIsotopes;i++)
126  {
127  G4String fractionref = GenerateName(elementPtr->GetIsotope(i)->GetName(),
128  elementPtr->GetIsotope(i));
129  xercesc::DOMElement* fractionElement = NewElement("fraction");
130  fractionElement->setAttributeNode(NewAttribute("n",
131  RelativeAbundanceVector[i]));
132  fractionElement->setAttributeNode(NewAttribute("ref",fractionref));
133  elementElement->appendChild(fractionElement);
134  AddIsotope(elementPtr->GetIsotope(i));
135  }
136  }
137  else
138  {
139  elementElement->setAttributeNode(NewAttribute("Z",elementPtr->GetZ()));
140  AtomWrite(elementElement,elementPtr->GetA());
141  }
142 
143  materialsElement->appendChild(elementElement);
144  // Append the element AFTER all the possible components are appended!
145 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const XML_Char * name
Definition: expat.h:151
const G4String & GetName() const
Definition: G4Isotope.hh:88
xercesc::DOMElement * materialsElement
G4double GetZ() const
Definition: G4Element.hh:131
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetA() const
Definition: G4Element.hh:139
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
void AtomWrite(xercesc::DOMElement *, const G4double &)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
void AddIsotope(const G4Isotope *const)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::IsotopeWrite ( const G4Isotope * const  isotopePtr)
protected

Definition at line 100 of file G4GDMLWriteMaterials.cc.

101 {
102  const G4String name = GenerateName(isotopePtr->GetName(),isotopePtr);
103 
104  xercesc::DOMElement* isotopeElement = NewElement("isotope");
105  isotopeElement->setAttributeNode(NewAttribute("name",name));
106  isotopeElement->setAttributeNode(NewAttribute("N",isotopePtr->GetN()));
107  isotopeElement->setAttributeNode(NewAttribute("Z",isotopePtr->GetZ()));
108  materialsElement->appendChild(isotopeElement);
109  AtomWrite(isotopeElement,isotopePtr->GetA());
110 }
const XML_Char * name
Definition: expat.h:151
G4double GetA() const
Definition: G4Isotope.hh:97
const G4String & GetName() const
Definition: G4Isotope.hh:88
xercesc::DOMElement * materialsElement
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4int GetN() const
Definition: G4Isotope.hh:94
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4int GetZ() const
Definition: G4Isotope.hh:91
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
void AtomWrite(xercesc::DOMElement *, const G4double &)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::MaterialsWrite ( xercesc::DOMElement *  element)
virtual

Implements G4GDMLWrite.

Definition at line 274 of file G4GDMLWriteMaterials.cc.

275 {
276  G4cout << "G4GDML: Writing materials..." << G4endl;
277 
278  materialsElement = NewElement("materials");
279  element->appendChild(materialsElement);
280 
281  isotopeList.clear();
282  elementList.clear();
283  materialList.clear();
284 }
xercesc::DOMElement * materialsElement
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
std::vector< const G4Element * > elementList
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
std::vector< const G4Material * > materialList
std::vector< const G4Isotope * > isotopeList

Here is the call graph for this function:

void G4GDMLWriteMaterials::MaterialWrite ( const G4Material * const  materialPtr)
protected

Definition at line 147 of file G4GDMLWriteMaterials.cc.

148 {
149  G4String state_str("undefined");
150  const G4State state = materialPtr->GetState();
151  if (state==kStateSolid) { state_str = "solid"; } else
152  if (state==kStateLiquid) { state_str = "liquid"; } else
153  if (state==kStateGas) { state_str = "gas"; }
154 
155  const G4String name = GenerateName(materialPtr->GetName(), materialPtr);
156 
157  xercesc::DOMElement* materialElement = NewElement("material");
158  materialElement->setAttributeNode(NewAttribute("name",name));
159  materialElement->setAttributeNode(NewAttribute("state",state_str));
160 
161  // Write any property attached to the material...
162  //
163  if (materialPtr->GetMaterialPropertiesTable())
164  {
165  PropertyWrite(materialElement, materialPtr);
166  }
167 
168  if (materialPtr->GetTemperature() != STP_Temperature)
169  { TWrite(materialElement,materialPtr->GetTemperature()); }
170  if (materialPtr->GetPressure() != STP_Pressure)
171  { PWrite(materialElement,materialPtr->GetPressure()); }
172 
173  // Write Ionisation potential (mean excitation energy)
174  MEEWrite(materialElement,materialPtr->GetIonisation()->GetMeanExcitationEnergy());
175 
176  DWrite(materialElement,materialPtr->GetDensity());
177 
178  const size_t NumberOfElements = materialPtr->GetNumberOfElements();
179 
180  if ( (NumberOfElements>1)
181  || ( materialPtr->GetElement(0)
182  && materialPtr->GetElement(0)->GetNumberOfIsotopes()>1 ) )
183  {
184  const G4double* MassFractionVector = materialPtr->GetFractionVector();
185 
186  for (size_t i=0;i<NumberOfElements;i++)
187  {
188  const G4String fractionref =
189  GenerateName(materialPtr->GetElement(i)->GetName(),
190  materialPtr->GetElement(i));
191  xercesc::DOMElement* fractionElement = NewElement("fraction");
192  fractionElement->setAttributeNode(NewAttribute("n",
193  MassFractionVector[i]));
194  fractionElement->setAttributeNode(NewAttribute("ref",fractionref));
195  materialElement->appendChild(fractionElement);
196  AddElement(materialPtr->GetElement(i));
197  }
198  }
199  else
200  {
201  materialElement->setAttributeNode(NewAttribute("Z",materialPtr->GetZ()));
202  AtomWrite(materialElement,materialPtr->GetA());
203  }
204 
205  // Append the material AFTER all the possible components are appended!
206  //
207  materialsElement->appendChild(materialElement);
208 }
G4double GetPressure() const
Definition: G4Material.hh:183
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const XML_Char * name
Definition: expat.h:151
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double GetZ() const
Definition: G4Material.cc:623
void MEEWrite(xercesc::DOMElement *, const G4double &)
void DWrite(xercesc::DOMElement *, const G4double &)
G4State
Definition: G4Material.hh:114
const G4String & GetName() const
Definition: G4Material.hh:178
xercesc::DOMElement * materialsElement
G4double GetDensity() const
Definition: G4Material.hh:180
static constexpr double STP_Pressure
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
void TWrite(xercesc::DOMElement *, const G4double &)
static constexpr double STP_Temperature
void PWrite(xercesc::DOMElement *, const G4double &)
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetA() const
Definition: G4Material.cc:636
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
void AddElement(const G4Element *const)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
void AtomWrite(xercesc::DOMElement *, const G4double &)
G4double GetMeanExcitationEnergy() const
G4double GetTemperature() const
Definition: G4Material.hh:182
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4State GetState() const
Definition: G4Material.hh:181
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
const G4double * GetFractionVector() const
Definition: G4Material.hh:194
void PropertyWrite(xercesc::DOMElement *, const G4Material *const)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::MEEWrite ( xercesc::DOMElement *  element,
const G4double MEE 
)
protected

Definition at line 91 of file G4GDMLWriteMaterials.cc.

92 {
93  xercesc::DOMElement* PElement = NewElement("MEE");
94  PElement->setAttributeNode(NewAttribute("unit","eV"));
95  PElement->setAttributeNode(NewAttribute("value",MEE/electronvolt));
96  element->appendChild(PElement);
97 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double electronvolt
Definition: G4SIunits.hh:206
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::PropertyVectorWrite ( const G4String key,
const G4PhysicsOrderedFreeVector * const  pvec 
)
protected

Definition at line 210 of file G4GDMLWriteMaterials.cc.

212 {
213  const G4String matrixref = GenerateName(key, pvec);
214  xercesc::DOMElement* matrixElement = NewElement("matrix");
215  matrixElement->setAttributeNode(NewAttribute("name", matrixref));
216  matrixElement->setAttributeNode(NewAttribute("coldim", "2"));
217  std::ostringstream pvalues;
218  for (size_t i=0; i<pvec->GetVectorLength(); i++)
219  {
220  if (i!=0) { pvalues << " "; }
221  pvalues << pvec->Energy(i) << " " << (*pvec)[i];
222  }
223  matrixElement->setAttributeNode(NewAttribute("values", pvalues.str()));
224 
225  defineElement->appendChild(matrixElement);
226 }
xercesc::DOMElement * defineElement
size_t GetVectorLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double Energy(size_t index) const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::PropertyWrite ( xercesc::DOMElement *  matElement,
const G4Material * const  mat 
)
protected

Definition at line 228 of file G4GDMLWriteMaterials.cc.

230 {
231  xercesc::DOMElement* propElement;
233  const std::map< G4String, G4PhysicsOrderedFreeVector*,
234  std::less<G4String> >* pmap = ptable->GetPropertiesMap();
235  const std::map< G4String, G4double,
236  std::less<G4String> >* cmap = ptable->GetPropertiesCMap();
237  std::map< G4String, G4PhysicsOrderedFreeVector*,
238  std::less<G4String> >::const_iterator mpos;
239  std::map< G4String, G4double,
240  std::less<G4String> >::const_iterator cpos;
241  for (mpos=pmap->begin(); mpos!=pmap->end(); mpos++)
242  {
243  propElement = NewElement("property");
244  propElement->setAttributeNode(NewAttribute("name", mpos->first));
245  propElement->setAttributeNode(NewAttribute("ref",
246  GenerateName(mpos->first, mpos->second)));
247  if (mpos->second)
248  {
249  PropertyVectorWrite(mpos->first, mpos->second);
250  matElement->appendChild(propElement);
251  }
252  else
253  {
254  G4String warn_message = "Null pointer for material property -"
255  + mpos->first + "- of material -" + mat->GetName() + "- !";
256  G4Exception("G4GDMLWriteMaterials::PropertyWrite()", "NullPointer",
257  JustWarning, warn_message);
258  continue;
259  }
260  }
261  for (cpos=cmap->begin(); cpos!=cmap->end(); cpos++)
262  {
263  propElement = NewElement("property");
264  propElement->setAttributeNode(NewAttribute("name", cpos->first));
265  propElement->setAttributeNode(NewAttribute("ref", cpos->first));
266  xercesc::DOMElement* constElement = NewElement("constant");
267  constElement->setAttributeNode(NewAttribute("name", cpos->first));
268  constElement->setAttributeNode(NewAttribute("value", cpos->second));
269  defineElement->appendChild(constElement);
270  matElement->appendChild(propElement);
271  }
272 }
const G4String & GetName() const
Definition: G4Material.hh:178
xercesc::DOMElement * defineElement
const std::map< G4String, G4MaterialPropertyVector *, std::less< G4String > > * GetPropertiesMap() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
const std::map< G4String, G4double, std::less< G4String > > * GetPropertiesCMap() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void PropertyVectorWrite(const G4String &, const G4PhysicsOrderedFreeVector *const)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::PWrite ( xercesc::DOMElement *  element,
const G4double P 
)
protected

Definition at line 73 of file G4GDMLWriteMaterials.cc.

74 {
75  xercesc::DOMElement* PElement = NewElement("P");
76  PElement->setAttributeNode(NewAttribute("unit","pascal"));
77  PElement->setAttributeNode(NewAttribute("value",P/hep_pascal));
78  element->appendChild(PElement);
79 }
static constexpr double hep_pascal
Definition: G4SIunits.hh:235
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static double P[]
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteMaterials::TWrite ( xercesc::DOMElement *  element,
const G4double T 
)
protected

Definition at line 82 of file G4GDMLWriteMaterials.cc.

83 {
84  xercesc::DOMElement* TElement = NewElement("T");
85  TElement->setAttributeNode(NewAttribute("unit","K"));
86  TElement->setAttributeNode(NewAttribute("value",T/kelvin));
87  element->appendChild(TElement);
88 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double kelvin
Definition: G4SIunits.hh:281
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

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

Definition at line 83 of file G4GDMLWriteMaterials.hh.

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

Definition at line 82 of file G4GDMLWriteMaterials.hh.

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

Definition at line 84 of file G4GDMLWriteMaterials.hh.

xercesc::DOMElement* G4GDMLWriteMaterials::materialsElement
protected

Definition at line 85 of file G4GDMLWriteMaterials.hh.


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