Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLWriteMaterials Class Reference

#include <G4GDMLWriteMaterials.hh>

+ Inheritance 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)
 
void SetOutputFileOverwrite (G4bool flag)
 
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 G4PhysicsFreeVector *const)
 
void PropertyConstWrite (const G4String &, const G4double, const G4MaterialPropertiesTable *)
 
- 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
 
std::vector< const G4PhysicsFreeVector * > propertyList
 
xercesc::DOMElement * materialsElement = nullptr
 
- Protected Attributes inherited from G4GDMLWriteDefine
xercesc::DOMElement * defineElement = nullptr
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc = nullptr
 
xercesc::DOMElement * extElement = nullptr
 
xercesc::DOMElement * userinfoElement = nullptr
 
G4GDMLAuxListType auxList
 
G4bool overwriteOutputFile = false
 

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 48 of file G4GDMLWriteMaterials.hh.

Constructor & Destructor Documentation

◆ G4GDMLWriteMaterials()

G4GDMLWriteMaterials::G4GDMLWriteMaterials ( )
protected

Definition at line 42 of file G4GDMLWriteMaterials.cc.

◆ ~G4GDMLWriteMaterials()

G4GDMLWriteMaterials::~G4GDMLWriteMaterials ( )
protectedvirtual

Definition at line 48 of file G4GDMLWriteMaterials.cc.

49{
50}

Member Function Documentation

◆ AddElement()

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

Definition at line 350 of file G4GDMLWriteMaterials.cc.

351{
352 for(std::size_t i = 0; i < elementList.size(); ++i) // Check if element is
353 { // already in the list!
354 if(elementList[i] == elementPtr)
355 {
356 return;
357 }
358 }
359 elementList.push_back(elementPtr);
360 ElementWrite(elementPtr);
361}
std::vector< const G4Element * > elementList
void ElementWrite(const G4Element *const)

Referenced by MaterialWrite().

◆ AddIsotope()

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

Definition at line 336 of file G4GDMLWriteMaterials.cc.

337{
338 for(std::size_t i = 0; i < isotopeList.size(); ++i) // Check if isotope is
339 { // already in the list!
340 if(isotopeList[i] == isotopePtr)
341 {
342 return;
343 }
344 }
345 isotopeList.push_back(isotopePtr);
346 IsotopeWrite(isotopePtr);
347}
void IsotopeWrite(const G4Isotope *const)
std::vector< const G4Isotope * > isotopeList

Referenced by ElementWrite().

◆ AddMaterial()

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

Definition at line 364 of file G4GDMLWriteMaterials.cc.

365{
366 for(std::size_t i = 0; i < materialList.size(); ++i) // Check if material is
367 { // already in the list!
368 if(materialList[i] == materialPtr)
369 {
370 return;
371 }
372 }
373 materialList.push_back(materialPtr);
374 MaterialWrite(materialPtr);
375}
void MaterialWrite(const G4Material *const)
std::vector< const G4Material * > materialList

Referenced by G4GDMLWriteStructure::TraverseVolumeTree().

◆ AtomWrite()

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

Definition at line 53 of file G4GDMLWriteMaterials.cc.

55{
56 xercesc::DOMElement* atomElement = NewElement("atom");
57 atomElement->setAttributeNode(NewAttribute("unit", "g/mole"));
58 atomElement->setAttributeNode(NewAttribute("value", a * mole / g));
59 element->appendChild(atomElement);
60}
xercesc::DOMElement * NewElement(const G4String &)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)

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

◆ DWrite()

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

Definition at line 63 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}

Referenced by MaterialWrite().

◆ ElementWrite()

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

Definition at line 116 of file G4GDMLWriteMaterials.cc.

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

Referenced by AddElement().

◆ IsotopeWrite()

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

Definition at line 103 of file G4GDMLWriteMaterials.cc.

104{
105 const G4String name = GenerateName(isotopePtr->GetName(), isotopePtr);
106
107 xercesc::DOMElement* isotopeElement = NewElement("isotope");
108 isotopeElement->setAttributeNode(NewAttribute("name", name));
109 isotopeElement->setAttributeNode(NewAttribute("N", isotopePtr->GetN()));
110 isotopeElement->setAttributeNode(NewAttribute("Z", isotopePtr->GetZ()));
111 materialsElement->appendChild(isotopeElement);
112 AtomWrite(isotopeElement, isotopePtr->GetA());
113}
G4int GetZ() const
Definition G4Isotope.hh:80
G4int GetN() const
Definition G4Isotope.hh:83
G4double GetA() const
Definition G4Isotope.hh:86

Referenced by AddIsotope().

◆ MaterialsWrite()

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

Implements G4GDMLWrite.

Definition at line 321 of file G4GDMLWriteMaterials.cc.

322{
323#ifdef G4VERBOSE
324 G4cout << "G4GDML: Writing materials..." << G4endl;
325#endif
326 materialsElement = NewElement("materials");
327 element->appendChild(materialsElement);
328
329 isotopeList.clear();
330 elementList.clear();
331 materialList.clear();
332 propertyList.clear();
333}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::vector< const G4PhysicsFreeVector * > propertyList

◆ MaterialWrite()

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

Definition at line 152 of file G4GDMLWriteMaterials.cc.

153{
154 G4String state_str("undefined");
155 const G4State state = materialPtr->GetState();
156 if(state == kStateSolid)
157 {
158 state_str = "solid";
159 }
160 else if(state == kStateLiquid)
161 {
162 state_str = "liquid";
163 }
164 else if(state == kStateGas)
165 {
166 state_str = "gas";
167 }
168
169 const G4String name = GenerateName(materialPtr->GetName(), materialPtr);
170
171 xercesc::DOMElement* materialElement = NewElement("material");
172 materialElement->setAttributeNode(NewAttribute("name", name));
173 materialElement->setAttributeNode(NewAttribute("state", state_str));
174
175 // Write any property attached to the material...
176 //
177 if(materialPtr->GetMaterialPropertiesTable())
178 {
179 PropertyWrite(materialElement, materialPtr);
180 }
181
182 if(materialPtr->GetTemperature() != STP_Temperature)
183 {
184 TWrite(materialElement, materialPtr->GetTemperature());
185 }
186
187 if(materialPtr->GetPressure() != STP_Pressure)
188 {
189 PWrite(materialElement, materialPtr->GetPressure());
190 }
191
192 // Write Ionisation potential (mean excitation energy)
193 MEEWrite(materialElement,
194 materialPtr->GetIonisation()->GetMeanExcitationEnergy());
195
196 DWrite(materialElement, materialPtr->GetDensity());
197
198 const G4int NumberOfElements = (G4int)materialPtr->GetNumberOfElements();
199
200 if((NumberOfElements > 1) ||
201 (materialPtr->GetElement(0) != nullptr &&
202 materialPtr->GetElement(0)->GetNumberOfIsotopes() > 1))
203 {
204 const G4double* MassFractionVector = materialPtr->GetFractionVector();
205
206 for(G4int i = 0; i < NumberOfElements; ++i)
207 {
208 const G4String fractionref = GenerateName(
209 materialPtr->GetElement(i)->GetName(), materialPtr->GetElement(i));
210 xercesc::DOMElement* fractionElement = NewElement("fraction");
211 fractionElement->setAttributeNode(
212 NewAttribute("n", MassFractionVector[i]));
213 fractionElement->setAttributeNode(NewAttribute("ref", fractionref));
214 materialElement->appendChild(fractionElement);
215 AddElement(materialPtr->GetElement(i));
216 }
217 }
218 else
219 {
220 materialElement->setAttributeNode(NewAttribute("Z", materialPtr->GetZ()));
221 AtomWrite(materialElement, materialPtr->GetA());
222 }
223
224 // Append the material AFTER all the possible components are appended!
225 //
226 materialsElement->appendChild(materialElement);
227}
G4State
@ kStateSolid
@ kStateLiquid
@ kStateGas
void TWrite(xercesc::DOMElement *, const G4double &)
void MEEWrite(xercesc::DOMElement *, const G4double &)
void DWrite(xercesc::DOMElement *, const G4double &)
void PWrite(xercesc::DOMElement *, const G4double &)
void AddElement(const G4Element *const)
void PropertyWrite(xercesc::DOMElement *, const G4Material *const)
G4double GetMeanExcitationEnergy() const
G4double GetPressure() const
G4double GetDensity() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
G4State GetState() const
G4double GetTemperature() const
const G4Element * GetElement(G4int iel) const
G4double GetZ() const
const G4double * GetFractionVector() const
G4IonisParamMat * GetIonisation() const
G4double GetA() const
std::size_t GetNumberOfElements() const
const G4String & GetName() const

Referenced by AddMaterial().

◆ MEEWrite()

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

Definition at line 93 of file G4GDMLWriteMaterials.cc.

95{
96 xercesc::DOMElement* PElement = NewElement("MEE");
97 PElement->setAttributeNode(NewAttribute("unit", "eV"));
98 PElement->setAttributeNode(NewAttribute("value", MEE / electronvolt));
99 element->appendChild(PElement);
100}

Referenced by MaterialWrite().

◆ PropertyConstWrite()

void G4GDMLWriteMaterials::PropertyConstWrite ( const G4String & key,
const G4double pval,
const G4MaterialPropertiesTable * ptable )
protected

Definition at line 261 of file G4GDMLWriteMaterials.cc.

264{
265 const G4String matrixref = GenerateName(key, ptable);
266 xercesc::DOMElement* matrixElement = NewElement("matrix");
267 matrixElement->setAttributeNode(NewAttribute("name", matrixref));
268 matrixElement->setAttributeNode(NewAttribute("coldim", "1"));
269 std::ostringstream pvalues;
270
271 pvalues << pval;
272 matrixElement->setAttributeNode(NewAttribute("values", pvalues.str()));
273
274 defineElement->appendChild(matrixElement);
275}
xercesc::DOMElement * defineElement

Referenced by PropertyWrite().

◆ PropertyVectorWrite()

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

Definition at line 230 of file G4GDMLWriteMaterials.cc.

232{
233 for(std::size_t i = 0; i < propertyList.size(); ++i) // Check if property is
234 { // already in the list!
235 if(propertyList[i] == pvec)
236 {
237 return;
238 }
239 }
240 propertyList.push_back(pvec);
241
242 const G4String matrixref = GenerateName(key, pvec);
243 xercesc::DOMElement* matrixElement = NewElement("matrix");
244 matrixElement->setAttributeNode(NewAttribute("name", matrixref));
245 matrixElement->setAttributeNode(NewAttribute("coldim", "2"));
246 std::ostringstream pvalues;
247 for(std::size_t i = 0; i < pvec->GetVectorLength(); ++i)
248 {
249 if(i != 0)
250 {
251 pvalues << " ";
252 }
253 pvalues << pvec->Energy(i) << " " << (*pvec)[i];
254 }
255 matrixElement->setAttributeNode(NewAttribute("values", pvalues.str()));
256
257 defineElement->appendChild(matrixElement);
258}
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const

Referenced by PropertyWrite(), and G4GDMLWriteSolids::PropertyWrite().

◆ PropertyWrite()

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

Definition at line 278 of file G4GDMLWriteMaterials.cc.

280{
281 xercesc::DOMElement* propElement;
283
284 auto pvec = ptable->GetProperties();
285 auto cvec = ptable->GetConstProperties();
286
287 for(size_t i = 0; i < pvec.size(); ++i)
288 {
289 if (pvec[i] != nullptr)
290 {
291 propElement = NewElement("property");
292 propElement->setAttributeNode(
293 NewAttribute("name", ptable->GetMaterialPropertyNames()[i]));
294 propElement->setAttributeNode(NewAttribute(
295 "ref", GenerateName(ptable->GetMaterialPropertyNames()[i],
296 pvec[i])));
298 pvec[i]);
299 matElement->appendChild(propElement);
300 }
301 }
302
303 for(size_t i = 0; i < cvec.size(); ++i)
304 {
305 if (cvec[i].second == true)
306 {
307 propElement = NewElement("property");
308 propElement->setAttributeNode(NewAttribute(
309 "name", ptable->GetMaterialConstPropertyNames()[i]));
310 propElement->setAttributeNode(NewAttribute(
312 ptable)));
314 cvec[i].first, ptable);
315 matElement->appendChild(propElement);
316 }
317 }
318}
void PropertyVectorWrite(const G4String &, const G4PhysicsFreeVector *const)
void PropertyConstWrite(const G4String &, const G4double, const G4MaterialPropertiesTable *)
const std::vector< G4String > & GetMaterialConstPropertyNames() const
const std::vector< std::pair< G4double, G4bool > > & GetConstProperties() const
const std::vector< G4MaterialPropertyVector * > & GetProperties() const
const std::vector< G4String > & GetMaterialPropertyNames() const

Referenced by MaterialWrite().

◆ PWrite()

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

Definition at line 73 of file G4GDMLWriteMaterials.cc.

75{
76 xercesc::DOMElement* PElement = NewElement("P");
77 PElement->setAttributeNode(NewAttribute("unit", "pascal"));
78 PElement->setAttributeNode(NewAttribute("value", P / hep_pascal));
79 element->appendChild(PElement);
80}

Referenced by MaterialWrite().

◆ TWrite()

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

Definition at line 83 of file G4GDMLWriteMaterials.cc.

85{
86 xercesc::DOMElement* TElement = NewElement("T");
87 TElement->setAttributeNode(NewAttribute("unit", "K"));
88 TElement->setAttributeNode(NewAttribute("value", T / kelvin));
89 element->appendChild(TElement);
90}

Referenced by MaterialWrite().

Member Data Documentation

◆ elementList

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

Definition at line 80 of file G4GDMLWriteMaterials.hh.

Referenced by AddElement(), and MaterialsWrite().

◆ isotopeList

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

Definition at line 79 of file G4GDMLWriteMaterials.hh.

Referenced by AddIsotope(), and MaterialsWrite().

◆ materialList

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

Definition at line 81 of file G4GDMLWriteMaterials.hh.

Referenced by AddMaterial(), and MaterialsWrite().

◆ materialsElement

xercesc::DOMElement* G4GDMLWriteMaterials::materialsElement = nullptr
protected

Definition at line 83 of file G4GDMLWriteMaterials.hh.

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

◆ propertyList

std::vector<const G4PhysicsFreeVector*> G4GDMLWriteMaterials::propertyList
protected

Definition at line 82 of file G4GDMLWriteMaterials.hh.

Referenced by MaterialsWrite(), and PropertyVectorWrite().


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