Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/persistency/gdml/src/G4GDMLWriteMaterials.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 33 #include "G4GDMLWriteMaterials.hh"
 34 
 35 #include "G4PhysicalConstants.hh"
 36 #include "G4SystemOfUnits.hh"
 37 #include "G4Element.hh"
 38 #include "G4Isotope.hh"
 39 #include "G4Material.hh"
 40 
 41 // --------------------------------------------------------------------
 42 G4GDMLWriteMaterials::G4GDMLWriteMaterials()
 43   : G4GDMLWriteDefine()
 44 {
 45 }
 46 
 47 // --------------------------------------------------------------------
 48 G4GDMLWriteMaterials::~G4GDMLWriteMaterials()
 49 {
 50 }
 51 
 52 // --------------------------------------------------------------------
 53 void 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 // --------------------------------------------------------------------
 63 void 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 // --------------------------------------------------------------------
 73 void 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 // --------------------------------------------------------------------
 83 void 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 // --------------------------------------------------------------------
 93 void 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 // --------------------------------------------------------------------
103 void G4GDMLWriteMaterials::IsotopeWrite(const G4Isotope* const isotopePtr)
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 // --------------------------------------------------------------------
116 void G4GDMLWriteMaterials::ElementWrite(const G4Element* const elementPtr)
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 }
150 
151 // --------------------------------------------------------------------
152 void G4GDMLWriteMaterials::MaterialWrite(const G4Material* const materialPtr)
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 }
228 
229 // --------------------------------------------------------------------
230 void G4GDMLWriteMaterials::PropertyVectorWrite(
231   const G4String& key, const G4PhysicsFreeVector* 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 // --------------------------------------------------------------------
261 void G4GDMLWriteMaterials::PropertyConstWrite(
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 // --------------------------------------------------------------------
278 void G4GDMLWriteMaterials::PropertyWrite(xercesc::DOMElement* matElement,
279                                          const G4Material* const mat)
280 {
281   xercesc::DOMElement* propElement;
282   G4MaterialPropertiesTable* ptable = mat->GetMaterialPropertiesTable();
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])));
297       PropertyVectorWrite(ptable->GetMaterialPropertyNames()[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(
311         "ref", GenerateName(ptable->GetMaterialConstPropertyNames()[i],
312                             ptable)));
313       PropertyConstWrite(ptable->GetMaterialConstPropertyNames()[i],
314                          cvec[i].first, ptable);
315       matElement->appendChild(propElement);
316     }
317   }
318 }
319 
320 // --------------------------------------------------------------------
321 void G4GDMLWriteMaterials::MaterialsWrite(xercesc::DOMElement* element)
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 }
334 
335 // --------------------------------------------------------------------
336 void G4GDMLWriteMaterials::AddIsotope(const G4Isotope* const isotopePtr)
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 }
348 
349 // --------------------------------------------------------------------
350 void G4GDMLWriteMaterials::AddElement(const G4Element* const elementPtr)
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 }
362 
363 // --------------------------------------------------------------------
364 void G4GDMLWriteMaterials::AddMaterial(const G4Material* const materialPtr)
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 }
376