Geant4  10.02
G4GDMLReadMaterials.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: G4GDMLReadMaterials.cc 89243 2015-03-27 16:24:39Z gcosmo $
28 // GEANT4 tag $ Name:$
29 //
30 // class G4GDMLReadMaterials Implementation
31 //
32 // Original author: Zoltan Torzsok, November 2007
33 //
34 // --------------------------------------------------------------------
35 
36 #include "G4GDMLReadMaterials.hh"
37 
38 #include "G4PhysicalConstants.hh"
39 #include "G4SystemOfUnits.hh"
40 #include "G4UnitsTable.hh"
41 #include "G4Element.hh"
42 #include "G4Isotope.hh"
43 #include "G4Material.hh"
44 #include "G4NistManager.hh"
45 
47 {
48 }
49 
51 {
52 }
53 
55 G4GDMLReadMaterials::AtomRead(const xercesc::DOMElement* const atomElement)
56 {
57  G4double value = 0.0;
58  G4double unit = g/mole;
59 
60  const xercesc::DOMNamedNodeMap* const attributes
61  = atomElement->getAttributes();
62  XMLSize_t attributeCount = attributes->getLength();
63 
64  for (XMLSize_t attribute_index=0;
65  attribute_index<attributeCount; attribute_index++)
66  {
67  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
68 
69  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
70  { continue; }
71 
72  const xercesc::DOMAttr* const attribute
73  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
74  if (!attribute)
75  {
76  G4Exception("G4GDMLReadMaterials::AtomRead()", "InvalidRead",
77  FatalException, "No attribute found!");
78  return value;
79  }
80  const G4String attName = Transcode(attribute->getName());
81  const G4String attValue = Transcode(attribute->getValue());
82 
83  if (attName=="value") { value = eval.Evaluate(attValue); } else
84  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); }
85  }
86 
87  return value*unit;
88 }
89 
91 CompositeRead(const xercesc::DOMElement* const compositeElement,G4String& ref)
92 {
93  G4int n = 0;
94 
95  const xercesc::DOMNamedNodeMap* const attributes
96  = compositeElement->getAttributes();
97  XMLSize_t attributeCount = attributes->getLength();
98 
99  for (XMLSize_t attribute_index=0;
100  attribute_index<attributeCount; attribute_index++)
101  {
102  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
103 
104  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
105  { continue; }
106 
107  const xercesc::DOMAttr* const attribute
108  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
109  if (!attribute)
110  {
111  G4Exception("G4GDMLReadMaterials::CompositeRead()", "InvalidRead",
112  FatalException, "No attribute found!");
113  return n;
114  }
115  const G4String attName = Transcode(attribute->getName());
116  const G4String attValue = Transcode(attribute->getValue());
117 
118  if (attName=="n") { n = eval.EvaluateInteger(attValue); } else
119  if (attName=="ref") { ref = attValue; }
120  }
121 
122  return n;
123 }
124 
125 G4double G4GDMLReadMaterials::DRead(const xercesc::DOMElement* const DElement)
126 {
127  G4double value = 0.0;
128  G4double unit = g/cm3;
129 
130  const xercesc::DOMNamedNodeMap* const attributes
131  = DElement->getAttributes();
132  XMLSize_t attributeCount = attributes->getLength();
133 
134  for (XMLSize_t attribute_index=0;
135  attribute_index<attributeCount; attribute_index++)
136  {
137  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
138 
139  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
140  { continue; }
141 
142  const xercesc::DOMAttr* const attribute
143  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
144  if (!attribute)
145  {
146  G4Exception("G4GDMLReadMaterials::DRead()", "InvalidRead",
147  FatalException, "No attribute found!");
148  return value;
149  }
150  const G4String attName = Transcode(attribute->getName());
151  const G4String attValue = Transcode(attribute->getValue());
152 
153  if (attName=="value") { value = eval.Evaluate(attValue); } else
154  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); }
155  }
156 
157  return value*unit;
158 }
159 
160 G4double G4GDMLReadMaterials::PRead(const xercesc::DOMElement* const PElement)
161 {
162  G4double value = STP_Pressure;
163  G4double unit = hep_pascal;
164 
165  const xercesc::DOMNamedNodeMap* const attributes = PElement->getAttributes();
166  XMLSize_t attributeCount = attributes->getLength();
167 
168  for (XMLSize_t attribute_index=0;
169  attribute_index<attributeCount; attribute_index++)
170  {
171  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
172 
173  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
174  { continue; }
175 
176  const xercesc::DOMAttr* const attribute
177  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
178  if (!attribute)
179  {
180  G4Exception("G4GDMLReadMaterials::PRead()", "InvalidRead",
181  FatalException, "No attribute found!");
182  return value;
183  }
184  const G4String attName = Transcode(attribute->getName());
185  const G4String attValue = Transcode(attribute->getValue());
186 
187  if (attName=="value") { value = eval.Evaluate(attValue); } else
188  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); }
189  }
190 
191  return value*unit;
192 }
193 
194 G4double G4GDMLReadMaterials::TRead(const xercesc::DOMElement* const TElement)
195 {
196  G4double value = STP_Temperature;
197  G4double unit = kelvin;
198 
199  const xercesc::DOMNamedNodeMap* const attributes = TElement->getAttributes();
200  XMLSize_t attributeCount = attributes->getLength();
201 
202  for (XMLSize_t attribute_index=0;
203  attribute_index<attributeCount; attribute_index++)
204  {
205  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
206 
207  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
208  { continue; }
209 
210  const xercesc::DOMAttr* const attribute
211  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
212  if (!attribute)
213  {
214  G4Exception("G4GDMLReadMaterials::TRead()", "InvalidRead",
215  FatalException, "No attribute found!");
216  return value;
217  }
218  const G4String attName = Transcode(attribute->getName());
219  const G4String attValue = Transcode(attribute->getValue());
220 
221  if (attName=="value") { value = eval.Evaluate(attValue); } else
222  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); }
223  }
224 
225  return value*unit;
226 }
227 
228 G4double G4GDMLReadMaterials::MEERead(const xercesc::DOMElement* const PElement)
229 {
230  G4double value = -1;
231  G4double unit = eV;
232 
233  const xercesc::DOMNamedNodeMap* const attributes = PElement->getAttributes();
234  XMLSize_t attributeCount = attributes->getLength();
235 
236  for (XMLSize_t attribute_index=0;
237  attribute_index<attributeCount; attribute_index++)
238  {
239  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
240 
241  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
242  { continue; }
243 
244  const xercesc::DOMAttr* const attribute
245  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
246  if (!attribute)
247  {
248  G4Exception("G4GDMLReadMaterials::MEERead()", "InvalidRead",
249  FatalException, "No attribute found!");
250  return value;
251  }
252  const G4String attName = Transcode(attribute->getName());
253  const G4String attValue = Transcode(attribute->getValue());
254 
255  if (attName=="value") { value = eval.Evaluate(attValue); } else
256  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue); }
257  }
258 
259  return value*unit;
260 }
261 
263 ElementRead(const xercesc::DOMElement* const elementElement)
264 {
265  G4String name;
266  G4String formula;
267  G4double a = 0.0;
268  G4double Z = 0.0;
269 
270  const xercesc::DOMNamedNodeMap* const attributes
271  = elementElement->getAttributes();
272  XMLSize_t attributeCount = attributes->getLength();
273 
274  for (XMLSize_t attribute_index=0;
275  attribute_index<attributeCount; attribute_index++)
276  {
277  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
278 
279  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
280  { continue; }
281 
282  const xercesc::DOMAttr* const attribute
283  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
284  if (!attribute)
285  {
286  G4Exception("G4GDMLReadMaterials::ElementRead()", "InvalidRead",
287  FatalException, "No attribute found!");
288  return;
289  }
290  const G4String attName = Transcode(attribute->getName());
291  const G4String attValue = Transcode(attribute->getValue());
292 
293  if (attName=="name") { name = GenerateName(attValue); } else
294  if (attName=="formula") { formula = attValue; } else
295  if (attName=="Z") { Z = eval.Evaluate(attValue); }
296  }
297 
298  G4int nComponents = 0;
299 
300  for (xercesc::DOMNode* iter = elementElement->getFirstChild();
301  iter != 0; iter = iter->getNextSibling())
302  {
303  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
304 
305  const xercesc::DOMElement* const child
306  = dynamic_cast<xercesc::DOMElement*>(iter);
307  if (!child)
308  {
309  G4Exception("G4GDMLReadMaterials::ElementRead()", "InvalidRead",
310  FatalException, "No child found!");
311  return;
312  }
313  const G4String tag = Transcode(child->getTagName());
314 
315  if (tag=="atom") { a = AtomRead(child); } else
316  if (tag=="fraction") { nComponents++; }
317  }
318 
319  if (nComponents>0)
320  {
321  MixtureRead(elementElement,
322  new G4Element(Strip(name),formula,nComponents));
323  }
324  else
325  {
326  new G4Element(Strip(name),formula,Z,a);
327  }
328 }
329 
331 FractionRead(const xercesc::DOMElement* const fractionElement, G4String& ref)
332 {
333  G4double n = 0.0;
334 
335  const xercesc::DOMNamedNodeMap* const attributes
336  = fractionElement->getAttributes();
337  XMLSize_t attributeCount = attributes->getLength();
338 
339  for (XMLSize_t attribute_index=0;
340  attribute_index<attributeCount; attribute_index++)
341  {
342  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
343 
344  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
345  { continue; }
346 
347  const xercesc::DOMAttr* const attribute
348  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
349  if (!attribute)
350  {
351  G4Exception("G4GDMLReadMaterials::FractionRead()", "InvalidRead",
352  FatalException, "No attribute found!");
353  return n;
354  }
355  const G4String attName = Transcode(attribute->getName());
356  const G4String attValue = Transcode(attribute->getValue());
357 
358  if (attName=="n") { n = eval.Evaluate(attValue); } else
359  if (attName=="ref") { ref = attValue; }
360  }
361 
362  return n;
363 }
364 
366 IsotopeRead(const xercesc::DOMElement* const isotopeElement)
367 {
368  G4String name;
369  G4int Z = 0;
370  G4int N = 0;
371  G4double a = 0.0;
372 
373  const xercesc::DOMNamedNodeMap* const attributes
374  = isotopeElement->getAttributes();
375  XMLSize_t attributeCount = attributes->getLength();
376 
377  for (XMLSize_t attribute_index=0;
378  attribute_index<attributeCount;attribute_index++)
379  {
380  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
381 
382  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
383  { continue; }
384 
385  const xercesc::DOMAttr* const attribute
386  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
387  if (!attribute)
388  {
389  G4Exception("G4GDMLReadMaterials::IsotopeRead()", "InvalidRead",
390  FatalException, "No attribute found!");
391  return;
392  }
393  const G4String attName = Transcode(attribute->getName());
394  const G4String attValue = Transcode(attribute->getValue());
395 
396  if (attName=="name") { name = GenerateName(attValue); } else
397  if (attName=="Z") { Z = eval.EvaluateInteger(attValue); } else
398  if (attName=="N") { N = eval.EvaluateInteger(attValue); }
399  }
400 
401  for (xercesc::DOMNode* iter = isotopeElement->getFirstChild();
402  iter != 0; iter = iter->getNextSibling())
403  {
404  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
405 
406  const xercesc::DOMElement* const child
407  = dynamic_cast<xercesc::DOMElement*>(iter);
408  if (!child)
409  {
410  G4Exception("G4GDMLReadMaterials::IsotopeRead()", "InvalidRead",
411  FatalException, "No child found!");
412  return;
413  }
414  const G4String tag = Transcode(child->getTagName());
415 
416  if (tag=="atom") { a = AtomRead(child); }
417  }
418 
419  new G4Isotope(Strip(name),Z,N,a);
420 }
421 
423 MaterialRead(const xercesc::DOMElement* const materialElement)
424 {
425  G4String name;
426  G4double Z = 0.0;
427  G4double a = 0.0;
428  G4double D = 0.0;
429  G4State state = kStateUndefined;
430  G4double T = STP_Temperature;
431  G4double P = STP_Pressure;
432  G4double MEE = -1.0;
433 
434  const xercesc::DOMNamedNodeMap* const attributes
435  = materialElement->getAttributes();
436  XMLSize_t attributeCount = attributes->getLength();
437 
438  for (XMLSize_t attribute_index=0;
439  attribute_index<attributeCount; attribute_index++)
440  {
441  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
442 
443  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
444  { continue; }
445 
446  const xercesc::DOMAttr* const attribute
447  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
448  if (!attribute)
449  {
450  G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
451  FatalException, "No attribute found!");
452  return;
453  }
454  const G4String attName = Transcode(attribute->getName());
455  const G4String attValue = Transcode(attribute->getValue());
456 
457  if (attName=="name") { name = GenerateName(attValue); } else
458  if (attName=="Z") { Z = eval.Evaluate(attValue); } else
459  if (attName=="state")
460  {
461  if (attValue=="solid") { state = kStateSolid; } else
462  if (attValue=="liquid") { state = kStateLiquid; } else
463  if (attValue=="gas") { state = kStateGas; }
464  }
465  }
466 
467  size_t nComponents = 0;
468 
469  for (xercesc::DOMNode* iter = materialElement->getFirstChild();
470  iter != 0; iter = iter->getNextSibling())
471  {
472  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
473 
474  const xercesc::DOMElement* const child
475  = dynamic_cast<xercesc::DOMElement*>(iter);
476  if (!child)
477  {
478  G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
479  FatalException, "No child found!");
480  return;
481  }
482  const G4String tag = Transcode(child->getTagName());
483 
484  if (tag=="atom") { a = AtomRead(child); } else
485  if (tag=="Dref") { D = GetQuantity(GenerateName(RefRead(child))); } else
486  if (tag=="Pref") { P = GetQuantity(GenerateName(RefRead(child))); } else
487  if (tag=="Tref") { T = GetQuantity(GenerateName(RefRead(child))); } else
488  if (tag=="MEEref") { MEE = GetQuantity(GenerateName(RefRead(child))); } else
489  if (tag=="D") { D = DRead(child); } else
490  if (tag=="P") { P = PRead(child); } else
491  if (tag=="T") { T = TRead(child); } else
492  if (tag=="MEE") { MEE = MEERead(child); } else
493  if (tag=="fraction" || tag=="composite") { nComponents++; }
494  }
495 
496  G4Material* material = 0;
497 
498  if (nComponents==0)
499  {
500  material = new G4Material(Strip(name),Z,a,D,state,T,P);
501  }
502  else
503  {
504  material = new G4Material(Strip(name),D,nComponents,state,T,P);
505  MixtureRead(materialElement, material);
506  }
507  if (MEE != -1) // ionisation potential (mean excitation energy)
508  {
509  material->GetIonisation()->SetMeanExcitationEnergy(MEE);
510  }
511 
512  for (xercesc::DOMNode* iter = materialElement->getFirstChild();
513  iter != 0; iter = iter->getNextSibling())
514  {
515  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
516 
517  const xercesc::DOMElement* const child
518  = dynamic_cast<xercesc::DOMElement*>(iter);
519  if (!child)
520  {
521  G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
522  FatalException, "No child found!");
523  return;
524  }
525  const G4String tag = Transcode(child->getTagName());
526 
527  if (tag=="property") { PropertyRead(child,material); }
528  }
529 }
530 
532 MixtureRead(const xercesc::DOMElement *const mixtureElement, G4Element *element)
533 {
534  for (xercesc::DOMNode* iter = mixtureElement->getFirstChild();
535  iter != 0; iter = iter->getNextSibling())
536  {
537  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
538 
539  const xercesc::DOMElement* const child
540  = dynamic_cast<xercesc::DOMElement*>(iter);
541  if (!child)
542  {
543  G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidRead",
544  FatalException, "No child found!");
545  return;
546  }
547  const G4String tag = Transcode(child->getTagName());
548 
549  if (tag=="fraction")
550  {
551  G4String ref;
552  G4double n = FractionRead(child,ref);
553  element->AddIsotope(GetIsotope(GenerateName(ref,true)),n);
554  }
555  }
556 }
557 
559 MixtureRead(const xercesc::DOMElement *const mixtureElement,
560  G4Material *material)
561 {
562  for (xercesc::DOMNode* iter = mixtureElement->getFirstChild();
563  iter != 0; iter = iter->getNextSibling())
564  {
565  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
566 
567  const xercesc::DOMElement* const child
568  = dynamic_cast<xercesc::DOMElement*>(iter);
569  if (!child)
570  {
571  G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidRead",
572  FatalException, "No child found!");
573  return;
574  }
575  const G4String tag = Transcode(child->getTagName());
576 
577  if (tag=="fraction")
578  {
579  G4String ref;
580  G4double n = FractionRead(child,ref);
581 
582  G4Material *materialPtr = GetMaterial(GenerateName(ref,true), false);
583  G4Element *elementPtr = GetElement(GenerateName(ref,true), false);
584 
585  if (elementPtr != 0) { material->AddElement(elementPtr,n); } else
586  if (materialPtr != 0) { material->AddMaterial(materialPtr,n); }
587 
588  if ((materialPtr == 0) && (elementPtr == 0))
589  {
590  G4String error_msg = "Referenced material/element '"
591  + GenerateName(ref,true) + "' was not found!";
592  G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidSetup",
593  FatalException, error_msg);
594  }
595  }
596  else if (tag=="composite")
597  {
598  G4String ref;
599  G4int n = CompositeRead(child,ref);
600 
601  G4Element *elementPtr = GetElement(GenerateName(ref,true));
602  material->AddElement(elementPtr,n);
603  }
604  }
605 }
606 
608 PropertyRead(const xercesc::DOMElement* const propertyElement,
609  G4Material* material)
610 {
611  G4String name;
612  G4String ref;
613  G4GDMLMatrix matrix;
614 
615  const xercesc::DOMNamedNodeMap* const attributes
616  = propertyElement->getAttributes();
617  XMLSize_t attributeCount = attributes->getLength();
618 
619  for (XMLSize_t attribute_index=0;
620  attribute_index<attributeCount; attribute_index++)
621  {
622  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
623 
624  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
625  { continue; }
626 
627  const xercesc::DOMAttr* const attribute
628  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
629  if (!attribute)
630  {
631  G4Exception("G4GDMLReadMaterials::PropertyRead()", "InvalidRead",
632  FatalException, "No attribute found!");
633  return;
634  }
635  const G4String attName = Transcode(attribute->getName());
636  const G4String attValue = Transcode(attribute->getValue());
637 
638  if (attName=="name") { name = GenerateName(attValue); } else
639  if (attName=="ref") { matrix = GetMatrix(ref=attValue); }
640  }
641 
642  /*
643  if (matrix.GetCols() != 2)
644  {
645  G4String error_msg = "Referenced matrix '" + ref
646  + "' should have \n two columns as a property table for material: "
647  + material->GetName();
648  G4Exception("G4GDMLReadMaterials::PropertyRead()", "InvalidRead",
649  FatalException, error_msg);
650  }
651  */
652 
653  if (matrix.GetRows() == 0) { return; }
654 
656  if (!matprop)
657  {
658  matprop = new G4MaterialPropertiesTable();
659  material->SetMaterialPropertiesTable(matprop);
660  }
661  if (matrix.GetCols() == 1) // constant property assumed
662  {
663  matprop->AddConstProperty(Strip(name), matrix.Get(0,0));
664  }
665  else // build the material properties vector
666  {
668  for (size_t i=0; i<matrix.GetRows(); i++)
669  {
670  propvect->InsertValues(matrix.Get(i,0),matrix.Get(i,1));
671  }
672  matprop->AddProperty(Strip(name),propvect);
673  }
674 }
675 
677 MaterialsRead(const xercesc::DOMElement* const materialsElement)
678 {
679  G4cout << "G4GDML: Reading materials..." << G4endl;
680 
681  for (xercesc::DOMNode* iter = materialsElement->getFirstChild();
682  iter != 0; iter = iter->getNextSibling())
683  {
684  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
685 
686  const xercesc::DOMElement* const child
687  = dynamic_cast<xercesc::DOMElement*>(iter);
688  if (!child)
689  {
690  G4Exception("G4GDMLReadMaterials::MaterialsRead()", "InvalidRead",
691  FatalException, "No child found!");
692  return;
693  }
694  const G4String tag = Transcode(child->getTagName());
695 
696  if (tag=="define") { DefineRead(child); } else
697  if (tag=="element") { ElementRead(child); } else
698  if (tag=="isotope") { IsotopeRead(child); } else
699  if (tag=="material") { MaterialRead(child); }
700  else
701  {
702  G4String error_msg = "Unknown tag in materials: " + tag;
703  G4Exception("G4GDMLReadMaterials::MaterialsRead()", "InvalidSetup",
704  FatalException, error_msg);
705  }
706  }
707 }
708 
710 GetElement(const G4String& ref, G4bool verbose) const
711 {
712  G4Element* elementPtr = G4Element::GetElement(ref,false);
713 
714  if (!elementPtr)
715  {
716  elementPtr = G4NistManager::Instance()->FindOrBuildElement(ref);
717  }
718 
719  if (verbose && !elementPtr)
720  {
721  G4String error_msg = "Referenced element '" + ref + "' was not found!";
722  G4Exception("G4GDMLReadMaterials::GetElement()", "InvalidRead",
723  FatalException, error_msg);
724  }
725 
726  return elementPtr;
727 }
728 
730  G4bool verbose) const
731 {
732  G4Isotope* isotopePtr = G4Isotope::GetIsotope(ref,false);
733 
734  if (verbose && !isotopePtr)
735  {
736  G4String error_msg = "Referenced isotope '" + ref + "' was not found!";
737  G4Exception("G4GDMLReadMaterials::GetIsotope()", "InvalidRead",
738  FatalException, error_msg);
739  }
740 
741  return isotopePtr;
742 }
743 
745  G4bool verbose) const
746 {
747  G4Material *materialPtr = G4Material::GetMaterial(ref,false);
748 
749  if (!materialPtr)
750  {
751  materialPtr = G4NistManager::Instance()->FindOrBuildMaterial(ref);
752  }
753 
754  if (verbose && !materialPtr)
755  {
756  G4String error_msg = "Referenced material '" + ref + "' was not found!";
757  G4Exception("G4GDMLReadMaterials::GetMaterial()", "InvalidRead",
758  FatalException, error_msg);
759  }
760 
761  return materialPtr;
762 }
G4double TRead(const xercesc::DOMElement *const)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
G4int EvaluateInteger(const G4String &)
G4double GetQuantity(const G4String &)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
G4double DRead(const xercesc::DOMElement *const)
static G4Element * GetElement(G4String name, G4bool warning=true)
Definition: G4Element.cc:409
void AddMaterial(G4Material *material, G4double fraction)
Definition: G4Material.cc:469
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
G4State
Definition: G4Material.hh:114
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
virtual void DefineRead(const xercesc::DOMElement *const)
void SetMeanExcitationEnergy(G4double value)
G4String name
Definition: TRTMaterials.hh:40
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
Definition: G4Material.hh:249
Definition: xmlparse.cc:187
static const double hep_pascal
Definition: G4SIunits.hh:232
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
void InsertValues(G4double energy, G4double value)
G4double a
Definition: TRTMaterials.hh:39
G4int nComponents
Definition: TRTMaterials.hh:41
G4double Get(size_t r, size_t c) const
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
static G4Isotope * GetIsotope(const G4String &name, G4bool warning=false)
Definition: G4Isotope.cc:196
size_t GetCols() const
static double P[]
G4MaterialPropertyVector * AddProperty(const char *key, G4double *PhotonEnergies, G4double *PropertyValues, G4int NumEntries)
G4String RefRead(const xercesc::DOMElement *const)
G4double FractionRead(const xercesc::DOMElement *const, G4String &)
void MixtureRead(const xercesc::DOMElement *const, G4Element *)
G4Isotope * GetIsotope(const G4String &, G4bool verbose=true) const
static G4double GetValueOf(const G4String &)
G4int CompositeRead(const xercesc::DOMElement *const, G4String &)
G4GLOB_DLL std::ostream G4cout
void IsotopeRead(const xercesc::DOMElement *const)
G4PhysicsOrderedFreeVector G4MaterialPropertyVector
bool G4bool
Definition: G4Types.hh:79
G4Element * GetElement(const G4String &, G4bool verbose=true) const
void AddIsotope(G4Isotope *isotope, G4double RelativeAbundance)
Definition: G4Element.cc:151
static const double cm3
Definition: G4SIunits.hh:120
const G4int n
void AddConstProperty(const char *key, G4double PropertyValue)
virtual void MaterialsRead(const xercesc::DOMElement *const)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double kelvin
Definition: G4SIunits.hh:278
static const double eV
Definition: G4SIunits.hh:212
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
G4double PRead(const xercesc::DOMElement *const)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
static const double g
Definition: G4SIunits.hh:180
G4double MEERead(const xercesc::DOMElement *const)
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:97
double D(double temp)
size_t GetRows() const
void PropertyRead(const xercesc::DOMElement *const, G4Material *)
static const double mole
Definition: G4SIunits.hh:283
#define G4endl
Definition: G4ios.hh:61
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:364
double G4double
Definition: G4Types.hh:76
void ElementRead(const xercesc::DOMElement *const)
G4Element * FindOrBuildElement(G4int Z, G4bool isotopes=true)
void MaterialRead(const xercesc::DOMElement *const)
G4double Evaluate(const G4String &)
G4GDMLMatrix GetMatrix(const G4String &)
G4double AtomRead(const xercesc::DOMElement *const)