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