Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
27// $Id$
28//
29// class G4GDMLWriteMaterials Implementation
30//
31// Original author: Zoltan Torzsok, November 2007
32//
33// --------------------------------------------------------------------
34
35#include <sstream>
36
38
40#include "G4SystemOfUnits.hh"
41#include "G4Element.hh"
42#include "G4Isotope.hh"
43#include "G4Material.hh"
44
46 : G4GDMLWriteDefine(), materialsElement(0)
47{
48}
49
51{
52}
53
55AtomWrite(xercesc::DOMElement* element,const G4double& a)
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}
62
64DWrite(xercesc::DOMElement* element,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
73PWrite(xercesc::DOMElement* element,const G4double& P)
74{
75 xercesc::DOMElement* PElement = NewElement("P");
76 PElement->setAttributeNode(NewAttribute("unit","pascal"));
77 PElement->setAttributeNode(NewAttribute("value",P/pascal));
78 element->appendChild(PElement);
79}
80
82TWrite(xercesc::DOMElement* element,const G4double& T)
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}
89
91MEEWrite(xercesc::DOMElement* element,const G4double& MEE)
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}
98
100IsotopeWrite(const G4Isotope* const isotopePtr)
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}
111
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}
146
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}
209
211 const G4PhysicsOrderedFreeVector* const pvec)
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}
227
228void G4GDMLWriteMaterials::PropertyWrite(xercesc::DOMElement* matElement,
229 const G4Material* const mat)
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();
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}
273
274void G4GDMLWriteMaterials::MaterialsWrite(xercesc::DOMElement* element)
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}
285
286void G4GDMLWriteMaterials::AddIsotope(const G4Isotope* const isotopePtr)
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}
295
296void G4GDMLWriteMaterials::AddElement(const G4Element* const elementPtr)
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}
305
306void G4GDMLWriteMaterials::AddMaterial(const G4Material* const materialPtr)
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}
@ JustWarning
G4State
Definition: G4Material.hh:114
@ kStateSolid
Definition: G4Material.hh:114
@ kStateLiquid
Definition: G4Material.hh:114
@ kStateGas
Definition: G4Material.hh:114
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define pascal
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:131
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:127
xercesc::DOMElement * defineElement
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 PropertyWrite(xercesc::DOMElement *, const G4Material *const)
void IsotopeWrite(const G4Isotope *const)
std::vector< const G4Isotope * > isotopeList
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:127
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:90
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:103
G4double GetMeanExcitationEnergy() const
G4int GetZ() const
Definition: G4Isotope.hh:91
G4int GetN() const
Definition: G4Isotope.hh:94
const G4String & GetName() const
Definition: G4Isotope.hh:88
G4double GetA() const
Definition: G4Isotope.hh:97
const std::map< G4String, G4MaterialPropertyVector *, std::less< G4String > > * GetPropertiesMap() const
const std::map< G4String, G4double, std::less< G4String > > * GetPropertiesCMap() const
G4double GetPressure() const
Definition: G4Material.hh:182
G4double GetDensity() const
Definition: G4Material.hh:179
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
G4State GetState() const
Definition: G4Material.hh:180
G4double GetTemperature() const
Definition: G4Material.hh:181
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:201
G4double GetZ() const
Definition: G4Material.cc:604
const G4double * GetFractionVector() const
Definition: G4Material.hh:193
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:225
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
G4double GetA() const
Definition: G4Material.cc:617
const G4String & GetName() const
Definition: G4Material.hh:177
size_t GetVectorLength() const
G4double Energy(size_t index) const
G4int first(char) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41