Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLWriteMaterials.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4GDMLWriteMaterials implementation
27//
28// Author: Zoltan Torzsok, November 2007
29// --------------------------------------------------------------------
30
31#include <sstream>
32
34
36#include "G4SystemOfUnits.hh"
37#include "G4Element.hh"
38#include "G4Isotope.hh"
39#include "G4Material.hh"
40
41// --------------------------------------------------------------------
44{
45}
46
47// --------------------------------------------------------------------
49{
50}
51
52// --------------------------------------------------------------------
53void G4GDMLWriteMaterials::AtomWrite(xercesc::DOMElement* element,
54 const G4double& a)
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}
61
62// --------------------------------------------------------------------
63void G4GDMLWriteMaterials::DWrite(xercesc::DOMElement* element,
64 const G4double& d)
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}
71
72// --------------------------------------------------------------------
73void G4GDMLWriteMaterials::PWrite(xercesc::DOMElement* element,
74 const G4double& P)
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}
81
82// --------------------------------------------------------------------
83void G4GDMLWriteMaterials::TWrite(xercesc::DOMElement* element,
84 const G4double& T)
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}
91
92// --------------------------------------------------------------------
93void G4GDMLWriteMaterials::MEEWrite(xercesc::DOMElement* element,
94 const G4double& MEE)
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}
101
102// --------------------------------------------------------------------
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}
114
115// --------------------------------------------------------------------
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 std::size_t NumberOfIsotopes = elementPtr->GetNumberOfIsotopes();
124
125 if(NumberOfIsotopes > 0)
126 {
127 const G4double* RelativeAbundanceVector
128 = elementPtr->GetRelativeAbundanceVector();
129 for(std::size_t 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}
150
151// --------------------------------------------------------------------
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 std::size_t NumberOfElements = 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(std::size_t 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}
228
229// --------------------------------------------------------------------
231 const G4String& key, const G4PhysicsOrderedFreeVector* const pvec)
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}
259
260// --------------------------------------------------------------------
262 const G4String& key, const G4double pval,
263 const G4MaterialPropertiesTable* ptable)
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}
276
277// --------------------------------------------------------------------
278void G4GDMLWriteMaterials::PropertyWrite(xercesc::DOMElement* matElement,
279 const G4Material* const mat)
280{
281 xercesc::DOMElement* propElement;
283
284 auto pmap = ptable->GetPropertyMap();
285 auto cmap = ptable->GetConstPropertyMap();
286
287 for(auto mpos = pmap->cbegin(); mpos != pmap->cend(); ++mpos)
288 {
289 propElement = NewElement("property");
290 propElement->setAttributeNode(
291 NewAttribute("name", ptable->GetMaterialPropertyNames()[mpos->first]));
292 propElement->setAttributeNode(NewAttribute(
293 "ref", GenerateName(ptable->GetMaterialPropertyNames()[mpos->first],
294 mpos->second)));
295 if(mpos->second)
296 {
297 PropertyVectorWrite(ptable->GetMaterialPropertyNames()[mpos->first],
298 mpos->second);
299 matElement->appendChild(propElement);
300 }
301 else
302 {
303 G4String warn_message = "Null pointer for material property -" +
304 ptable->GetMaterialPropertyNames()[mpos->first] +
305 "- of material -" + mat->GetName() + "- !";
306 G4Exception("G4GDMLWriteMaterials::PropertyWrite()", "NullPointer",
307 JustWarning, warn_message);
308 continue;
309 }
310 }
311
312 for(auto cpos = cmap->cbegin(); cpos != cmap->cend(); ++cpos)
313 {
314 propElement = NewElement("property");
315 propElement->setAttributeNode(NewAttribute(
316 "name", ptable->GetMaterialConstPropertyNames()[cpos->first]));
317 propElement->setAttributeNode(NewAttribute(
318 "ref", GenerateName(ptable->GetMaterialConstPropertyNames()[cpos->first],
319 ptable)));
320
322 cpos->second, ptable);
323
324 matElement->appendChild(propElement);
325 }
326}
327
328// --------------------------------------------------------------------
329void G4GDMLWriteMaterials::MaterialsWrite(xercesc::DOMElement* element)
330{
331#ifdef G4VERBOSE
332 G4cout << "G4GDML: Writing materials..." << G4endl;
333#endif
334 materialsElement = NewElement("materials");
335 element->appendChild(materialsElement);
336
337 isotopeList.clear();
338 elementList.clear();
339 materialList.clear();
340 propertyList.clear();
341}
342
343// --------------------------------------------------------------------
344void G4GDMLWriteMaterials::AddIsotope(const G4Isotope* const isotopePtr)
345{
346 for(std::size_t i = 0; i < isotopeList.size(); ++i) // Check if isotope is
347 { // already in the list!
348 if(isotopeList[i] == isotopePtr)
349 {
350 return;
351 }
352 }
353 isotopeList.push_back(isotopePtr);
354 IsotopeWrite(isotopePtr);
355}
356
357// --------------------------------------------------------------------
358void G4GDMLWriteMaterials::AddElement(const G4Element* const elementPtr)
359{
360 for(std::size_t i = 0; i < elementList.size(); ++i) // Check if element is
361 { // already in the list!
362 if(elementList[i] == elementPtr)
363 {
364 return;
365 }
366 }
367 elementList.push_back(elementPtr);
368 ElementWrite(elementPtr);
369}
370
371// --------------------------------------------------------------------
372void G4GDMLWriteMaterials::AddMaterial(const G4Material* const materialPtr)
373{
374 for(std::size_t i = 0; i < materialList.size(); ++i) // Check if material is
375 { // already in the list!
376 if(materialList[i] == materialPtr)
377 {
378 return;
379 }
380 }
381 materialList.push_back(materialPtr);
382 MaterialWrite(materialPtr);
383}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4State
Definition: G4Material.hh:111
@ kStateSolid
Definition: G4Material.hh:111
@ kStateLiquid
Definition: G4Material.hh:111
@ kStateGas
Definition: G4Material.hh:111
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:130
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4double GetA() const
Definition: G4Element.hh:138
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
const G4String & GetName() const
Definition: G4Element.hh:126
xercesc::DOMElement * defineElement
std::vector< const G4PhysicsOrderedFreeVector * > propertyList
void TWrite(xercesc::DOMElement *, const G4double &)
void MEEWrite(xercesc::DOMElement *, const G4double &)
void MaterialWrite(const G4Material *const)
void DWrite(xercesc::DOMElement *, const G4double &)
void AtomWrite(xercesc::DOMElement *, const G4double &)
virtual void MaterialsWrite(xercesc::DOMElement *)
void AddIsotope(const G4Isotope *const)
void PropertyVectorWrite(const G4String &, const G4PhysicsOrderedFreeVector *const)
void AddMaterial(const G4Material *const)
std::vector< const G4Element * > elementList
xercesc::DOMElement * materialsElement
std::vector< const G4Material * > materialList
void PWrite(xercesc::DOMElement *, const G4double &)
void ElementWrite(const G4Element *const)
void AddElement(const G4Element *const)
void PropertyConstWrite(const G4String &, const G4double, const G4MaterialPropertiesTable *)
void PropertyWrite(xercesc::DOMElement *, const G4Material *const)
void IsotopeWrite(const G4Isotope *const)
std::vector< const G4Isotope * > isotopeList
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:181
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:132
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:155
G4double GetMeanExcitationEnergy() const
G4int GetZ() const
Definition: G4Isotope.hh:90
G4int GetN() const
Definition: G4Isotope.hh:93
const G4String & GetName() const
Definition: G4Isotope.hh:87
G4double GetA() const
Definition: G4Isotope.hh:96
std::vector< G4String > GetMaterialPropertyNames() const
std::vector< G4String > GetMaterialConstPropertyNames() const
const std::map< G4int, G4double, std::less< G4int > > * GetConstPropertyMap() const
const std::map< G4int, G4MaterialPropertyVector *, std::less< G4int > > * GetPropertyMap() const
G4double GetPressure() const
Definition: G4Material.hh:181
G4double GetDensity() const
Definition: G4Material.hh:178
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:254
G4State GetState() const
Definition: G4Material.hh:179
G4double GetTemperature() const
Definition: G4Material.hh:180
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
G4double GetZ() const
Definition: G4Material.cc:701
const G4double * GetFractionVector() const
Definition: G4Material.hh:192
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4double GetA() const
Definition: G4Material.cc:714
const G4String & GetName() const
Definition: G4Material.hh:175
G4double Energy(std::size_t index) const
std::size_t GetVectorLength() const