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