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