Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// G4GDMLReadMaterials implementation
27//
28// Author: Zoltan Torzsok, November 2007
29// --------------------------------------------------------------------
30
32
34#include "G4SystemOfUnits.hh"
35#include "G4UnitsTable.hh"
36#include "G4Element.hh"
37#include "G4Isotope.hh"
38#include "G4Material.hh"
39#include "G4NistManager.hh"
40
41// --------------------------------------------------------------------
44{
45}
46
47// --------------------------------------------------------------------
49{
50}
51
52// --------------------------------------------------------------------
54 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; attribute_index < attributeCount;
64 ++attribute_index)
65 {
66 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
67
68 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
69 {
70 continue;
71 }
72
73 const xercesc::DOMAttr* const attribute =
74 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
75 if(attribute == nullptr)
76 {
77 G4Exception("G4GDMLReadMaterials::AtomRead()", "InvalidRead",
78 FatalException, "No attribute found!");
79 return value;
80 }
81 const G4String attName = Transcode(attribute->getName());
82 const G4String attValue = Transcode(attribute->getValue());
83
84 if(attName == "value")
85 {
86 value = eval.Evaluate(attValue);
87 }
88 else if(attName == "unit")
89 {
90 unit = G4UnitDefinition::GetValueOf(attValue);
91 if(G4UnitDefinition::GetCategory(attValue) != "Molar mass")
92 {
93 G4Exception("G4GDMLReadMaterials::AtomRead()", "InvalidRead",
94 FatalException, "Invalid unit for atomic mass!");
95 }
96 }
97 }
98
99 return value * unit;
100}
101
102// --------------------------------------------------------------------
104 const xercesc::DOMElement* const compositeElement, G4String& ref)
105{
106 G4int n = 0;
107
108 const xercesc::DOMNamedNodeMap* const attributes =
109 compositeElement->getAttributes();
110 XMLSize_t attributeCount = attributes->getLength();
111
112 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
113 ++attribute_index)
114 {
115 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
116
117 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
118 {
119 continue;
120 }
121
122 const xercesc::DOMAttr* const attribute =
123 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
124 if(attribute == nullptr)
125 {
126 G4Exception("G4GDMLReadMaterials::CompositeRead()", "InvalidRead",
127 FatalException, "No attribute found!");
128 return n;
129 }
130 const G4String attName = Transcode(attribute->getName());
131 const G4String attValue = Transcode(attribute->getValue());
132
133 if(attName == "n")
134 {
135 n = eval.EvaluateInteger(attValue);
136 }
137 else if(attName == "ref")
138 {
139 ref = attValue;
140 }
141 }
142
143 return n;
144}
145
146// --------------------------------------------------------------------
147G4double G4GDMLReadMaterials::DRead(const xercesc::DOMElement* const DElement)
148{
149 G4double value = 0.0;
150 G4double unit = g / cm3;
151
152 const xercesc::DOMNamedNodeMap* const attributes = DElement->getAttributes();
153 XMLSize_t attributeCount = attributes->getLength();
154
155 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
156 ++attribute_index)
157 {
158 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
159
160 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
161 {
162 continue;
163 }
164
165 const xercesc::DOMAttr* const attribute =
166 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
167 if(attribute == nullptr)
168 {
169 G4Exception("G4GDMLReadMaterials::DRead()", "InvalidRead", FatalException,
170 "No attribute found!");
171 return value;
172 }
173 const G4String attName = Transcode(attribute->getName());
174 const G4String attValue = Transcode(attribute->getValue());
175
176 if(attName == "value")
177 {
178 value = eval.Evaluate(attValue);
179 }
180 else if(attName == "unit")
181 {
182 unit = G4UnitDefinition::GetValueOf(attValue);
183 if(G4UnitDefinition::GetCategory(attValue) != "Volumic Mass")
184 {
185 G4Exception("G4GDMLReadMaterials::DRead()", "InvalidRead",
186 FatalException, "Invalid unit for density!");
187 }
188 }
189 }
190
191 return value * unit;
192}
193
194// --------------------------------------------------------------------
195G4double G4GDMLReadMaterials::PRead(const xercesc::DOMElement* const PElement)
196{
197 G4double value = STP_Pressure;
198 G4double unit = hep_pascal;
199
200 const xercesc::DOMNamedNodeMap* const attributes = PElement->getAttributes();
201 XMLSize_t attributeCount = attributes->getLength();
202
203 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
204 ++attribute_index)
205 {
206 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
207
208 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
209 {
210 continue;
211 }
212
213 const xercesc::DOMAttr* const attribute =
214 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
215 if(attribute == nullptr)
216 {
217 G4Exception("G4GDMLReadMaterials::PRead()", "InvalidRead", FatalException,
218 "No attribute found!");
219 return value;
220 }
221 const G4String attName = Transcode(attribute->getName());
222 const G4String attValue = Transcode(attribute->getValue());
223
224 if(attName == "value")
225 {
226 value = eval.Evaluate(attValue);
227 }
228 else if(attName == "unit")
229 {
230 unit = G4UnitDefinition::GetValueOf(attValue);
231 if(G4UnitDefinition::GetCategory(attValue) != "Pressure")
232 {
233 G4Exception("G4GDMLReadMaterials::PRead()", "InvalidRead",
234 FatalException, "Invalid unit for pressure!");
235 }
236 }
237 }
238
239 return value * unit;
240}
241
242// --------------------------------------------------------------------
243G4double G4GDMLReadMaterials::TRead(const xercesc::DOMElement* const TElement)
244{
245 G4double value = NTP_Temperature;
246 G4double unit = kelvin;
247
248 const xercesc::DOMNamedNodeMap* const attributes = TElement->getAttributes();
249 XMLSize_t attributeCount = attributes->getLength();
250
251 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
252 ++attribute_index)
253 {
254 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
255
256 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
257 {
258 continue;
259 }
260
261 const xercesc::DOMAttr* const attribute =
262 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
263 if(attribute == nullptr)
264 {
265 G4Exception("G4GDMLReadMaterials::TRead()", "InvalidRead", FatalException,
266 "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")
273 {
274 value = eval.Evaluate(attValue);
275 }
276 else if(attName == "unit")
277 {
278 unit = G4UnitDefinition::GetValueOf(attValue);
279 if(G4UnitDefinition::GetCategory(attValue) != "Temperature")
280 {
281 G4Exception("G4GDMLReadMaterials::TRead()", "InvalidRead",
282 FatalException, "Invalid unit for temperature!");
283 }
284 }
285 }
286
287 return value * unit;
288}
289
290// --------------------------------------------------------------------
291G4double G4GDMLReadMaterials::MEERead(const xercesc::DOMElement* const PElement)
292{
293 G4double value = -1;
294 G4double unit = eV;
295
296 const xercesc::DOMNamedNodeMap* const attributes = PElement->getAttributes();
297 XMLSize_t attributeCount = attributes->getLength();
298
299 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
300 ++attribute_index)
301 {
302 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
303
304 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
305 {
306 continue;
307 }
308
309 const xercesc::DOMAttr* const attribute =
310 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
311 if(attribute == nullptr)
312 {
313 G4Exception("G4GDMLReadMaterials::MEERead()", "InvalidRead",
314 FatalException, "No attribute found!");
315 return value;
316 }
317 const G4String attName = Transcode(attribute->getName());
318 const G4String attValue = Transcode(attribute->getValue());
319
320 if(attName == "value")
321 {
322 value = eval.Evaluate(attValue);
323 }
324 else if(attName == "unit")
325 {
326 unit = G4UnitDefinition::GetValueOf(attValue);
327 if(G4UnitDefinition::GetCategory(attValue) != "Energy")
328 {
329 G4Exception("G4GDMLReadMaterials::MEERead()", "InvalidRead",
330 FatalException, "Invalid unit for energy!");
331 }
332 }
333 }
334
335 return value * unit;
336}
337
338// --------------------------------------------------------------------
340 const xercesc::DOMElement* const elementElement)
341{
342 G4String name;
343 G4String formula;
344 G4double a = 0.0;
345 G4double Z = 0.0;
346
347 const xercesc::DOMNamedNodeMap* const attributes =
348 elementElement->getAttributes();
349 XMLSize_t attributeCount = attributes->getLength();
350
351 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
352 ++attribute_index)
353 {
354 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
355
356 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
357 {
358 continue;
359 }
360
361 const xercesc::DOMAttr* const attribute =
362 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
363 if(attribute == nullptr)
364 {
365 G4Exception("G4GDMLReadMaterials::ElementRead()", "InvalidRead",
366 FatalException, "No attribute found!");
367 return;
368 }
369 const G4String attName = Transcode(attribute->getName());
370 const G4String attValue = Transcode(attribute->getValue());
371
372 if(attName == "name")
373 {
374 name = GenerateName(attValue);
375 }
376 else if(attName == "formula")
377 {
378 formula = attValue;
379 }
380 else if(attName == "Z")
381 {
382 Z = eval.Evaluate(attValue);
383 }
384 }
385
386 G4int nComponents = 0;
387
388 for(xercesc::DOMNode* iter = elementElement->getFirstChild(); iter != nullptr;
389 iter = iter->getNextSibling())
390 {
391 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
392 {
393 continue;
394 }
395
396 const xercesc::DOMElement* const child =
397 dynamic_cast<xercesc::DOMElement*>(iter);
398 if(child == nullptr)
399 {
400 G4Exception("G4GDMLReadMaterials::ElementRead()", "InvalidRead",
401 FatalException, "No child found!");
402 return;
403 }
404 const G4String tag = Transcode(child->getTagName());
405
406 if(tag == "atom")
407 {
408 a = AtomRead(child);
409 }
410 else if(tag == "fraction")
411 {
412 nComponents++;
413 }
414 }
415
416 if(nComponents > 0)
417 {
418 MixtureRead(elementElement,
419 new G4Element(Strip(name), formula, nComponents));
420 }
421 else
422 {
423 new G4Element(Strip(name), formula, Z, a);
424 }
425}
426
427// --------------------------------------------------------------------
429 const xercesc::DOMElement* const fractionElement, G4String& ref)
430{
431 G4double n = 0.0;
432
433 const xercesc::DOMNamedNodeMap* const attributes =
434 fractionElement->getAttributes();
435 XMLSize_t attributeCount = attributes->getLength();
436
437 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
438 ++attribute_index)
439 {
440 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
441
442 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
443 {
444 continue;
445 }
446
447 const xercesc::DOMAttr* const attribute =
448 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
449 if(attribute == nullptr)
450 {
451 G4Exception("G4GDMLReadMaterials::FractionRead()", "InvalidRead",
452 FatalException, "No attribute found!");
453 return n;
454 }
455 const G4String attName = Transcode(attribute->getName());
456 const G4String attValue = Transcode(attribute->getValue());
457
458 if(attName == "n")
459 {
460 n = eval.Evaluate(attValue);
461 }
462 else if(attName == "ref")
463 {
464 ref = attValue;
465 }
466 }
467
468 return n;
469}
470
471// --------------------------------------------------------------------
473 const xercesc::DOMElement* const isotopeElement)
474{
475 G4String name;
476 G4int Z = 0;
477 G4int N = 0;
478 G4double a = 0.0;
479
480 const xercesc::DOMNamedNodeMap* const attributes =
481 isotopeElement->getAttributes();
482 XMLSize_t attributeCount = attributes->getLength();
483
484 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
485 ++attribute_index)
486 {
487 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
488
489 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
490 {
491 continue;
492 }
493
494 const xercesc::DOMAttr* const attribute =
495 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
496 if(attribute == nullptr)
497 {
498 G4Exception("G4GDMLReadMaterials::IsotopeRead()", "InvalidRead",
499 FatalException, "No attribute found!");
500 return;
501 }
502 const G4String attName = Transcode(attribute->getName());
503 const G4String attValue = Transcode(attribute->getValue());
504
505 if(attName == "name")
506 {
507 name = GenerateName(attValue);
508 }
509 else if(attName == "Z")
510 {
511 Z = eval.EvaluateInteger(attValue);
512 }
513 else if(attName == "N")
514 {
515 N = eval.EvaluateInteger(attValue);
516 }
517 }
518
519 for(xercesc::DOMNode* iter = isotopeElement->getFirstChild(); iter != nullptr;
520 iter = iter->getNextSibling())
521 {
522 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
523 {
524 continue;
525 }
526
527 const xercesc::DOMElement* const child =
528 dynamic_cast<xercesc::DOMElement*>(iter);
529 if(child == nullptr)
530 {
531 G4Exception("G4GDMLReadMaterials::IsotopeRead()", "InvalidRead",
532 FatalException, "No child found!");
533 return;
534 }
535 const G4String tag = Transcode(child->getTagName());
536
537 if(tag == "atom")
538 {
539 a = AtomRead(child);
540 }
541 }
542
543 new G4Isotope(Strip(name), Z, N, a);
544}
545
546// --------------------------------------------------------------------
548 const xercesc::DOMElement* const materialElement)
549{
550 G4String name;
551 G4double Z = 0.0;
552 G4double a = 0.0;
553 G4double D = 0.0;
554 G4State state = kStateUndefined;
555 G4double T = NTP_Temperature;
556 G4double P = STP_Pressure;
557 G4double MEE = -1.0;
558
559 const xercesc::DOMNamedNodeMap* const attributes =
560 materialElement->getAttributes();
561 XMLSize_t attributeCount = attributes->getLength();
562
563 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
564 ++attribute_index)
565 {
566 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
567
568 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
569 {
570 continue;
571 }
572
573 const xercesc::DOMAttr* const attribute =
574 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
575 if(attribute == nullptr)
576 {
577 G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
578 FatalException, "No attribute found!");
579 return;
580 }
581 const G4String attName = Transcode(attribute->getName());
582 const G4String attValue = Transcode(attribute->getValue());
583
584 if(attName == "name")
585 {
586 name = GenerateName(attValue);
587 }
588 else if(attName == "Z")
589 {
590 Z = eval.Evaluate(attValue);
591 }
592 else if(attName == "state")
593 {
594 if(attValue == "solid")
595 {
596 state = kStateSolid;
597 }
598 else if(attValue == "liquid")
599 {
600 state = kStateLiquid;
601 }
602 else if(attValue == "gas")
603 {
604 state = kStateGas;
605 }
606 }
607 }
608
609 std::size_t nComponents = 0;
610
611 for(xercesc::DOMNode* iter = materialElement->getFirstChild();
612 iter != nullptr; iter = iter->getNextSibling())
613 {
614 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
615 {
616 continue;
617 }
618
619 const xercesc::DOMElement* const child =
620 dynamic_cast<xercesc::DOMElement*>(iter);
621 if(child == nullptr)
622 {
623 G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
624 FatalException, "No child found!");
625 return;
626 }
627 const G4String tag = Transcode(child->getTagName());
628
629 if(tag == "atom")
630 {
631 a = AtomRead(child);
632 }
633 else if(tag == "Dref")
634 {
636 }
637 else if(tag == "Pref")
638 {
639 P = GetQuantity(GenerateName(RefRead(child)));
640 }
641 else if(tag == "Tref")
642 {
643 T = GetQuantity(GenerateName(RefRead(child)));
644 }
645 else if(tag == "MEEref")
646 {
647 MEE = GetQuantity(GenerateName(RefRead(child)));
648 }
649 else if(tag == "D")
650 {
651 D = DRead(child);
652 }
653 else if(tag == "P")
654 {
655 P = PRead(child);
656 }
657 else if(tag == "T")
658 {
659 T = TRead(child);
660 }
661 else if(tag == "MEE")
662 {
663 MEE = MEERead(child);
664 }
665 else if(tag == "fraction" || tag == "composite")
666 {
667 nComponents++;
668 }
669 }
670
671 G4Material* material = nullptr;
672
673 if(nComponents == 0)
674 {
675 material = new G4Material(Strip(name), Z, a, D, state, T, P);
676 }
677 else
678 {
679 material = new G4Material(Strip(name), D, (G4int)nComponents, state, T, P);
680 MixtureRead(materialElement, material);
681 }
682 if(MEE != -1) // ionisation potential (mean excitation energy)
683 {
684 material->GetIonisation()->SetMeanExcitationEnergy(MEE);
685 }
686
687 for(xercesc::DOMNode* iter = materialElement->getFirstChild();
688 iter != nullptr; iter = iter->getNextSibling())
689 {
690 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
691 {
692 continue;
693 }
694
695 const xercesc::DOMElement* const child =
696 dynamic_cast<xercesc::DOMElement*>(iter);
697 if(child == nullptr)
698 {
699 G4Exception("G4GDMLReadMaterials::MaterialRead()", "InvalidRead",
700 FatalException, "No child found!");
701 return;
702 }
703 const G4String tag = Transcode(child->getTagName());
704
705 if(tag == "property")
706 {
707 PropertyRead(child, material);
708 }
709 }
710}
711
712// --------------------------------------------------------------------
714 const xercesc::DOMElement* const mixtureElement, G4Element* element)
715{
716 for(xercesc::DOMNode* iter = mixtureElement->getFirstChild(); iter != nullptr;
717 iter = iter->getNextSibling())
718 {
719 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
720 {
721 continue;
722 }
723
724 const xercesc::DOMElement* const child =
725 dynamic_cast<xercesc::DOMElement*>(iter);
726 if(child == nullptr)
727 {
728 G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidRead",
729 FatalException, "No child found!");
730 return;
731 }
732 const G4String tag = Transcode(child->getTagName());
733
734 if(tag == "fraction")
735 {
736 G4String ref;
737 G4double n = FractionRead(child, ref);
738 element->AddIsotope(GetIsotope(GenerateName(ref, true)), n);
739 }
740 }
741}
742
743// --------------------------------------------------------------------
745 const xercesc::DOMElement* const mixtureElement, G4Material* material)
746{
747 for(xercesc::DOMNode* iter = mixtureElement->getFirstChild(); iter != nullptr;
748 iter = iter->getNextSibling())
749 {
750 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
751 {
752 continue;
753 }
754
755 const xercesc::DOMElement* const child =
756 dynamic_cast<xercesc::DOMElement*>(iter);
757 if(child == nullptr)
758 {
759 G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidRead",
760 FatalException, "No child found!");
761 return;
762 }
763 const G4String tag = Transcode(child->getTagName());
764
765 if(tag == "fraction")
766 {
767 G4String ref;
768 G4double n = FractionRead(child, ref);
769
770 G4Material* materialPtr = GetMaterial(GenerateName(ref, true), false);
771 G4Element* elementPtr = GetElement(GenerateName(ref, true), false);
772
773 if(elementPtr != nullptr)
774 {
775 material->AddElement(elementPtr, n);
776 }
777 else if(materialPtr != nullptr)
778 {
779 material->AddMaterial(materialPtr, n);
780 }
781
782 if((materialPtr == nullptr) && (elementPtr == nullptr))
783 {
784 G4String error_msg = "Referenced material/element '" +
785 GenerateName(ref, true) + "' was not found!";
786 G4Exception("G4GDMLReadMaterials::MixtureRead()", "InvalidSetup",
787 FatalException, error_msg);
788 }
789 }
790 else if(tag == "composite")
791 {
792 G4String ref;
793 G4int n = CompositeRead(child, ref);
794
795 G4Element* elementPtr = GetElement(GenerateName(ref, true));
796 material->AddElement(elementPtr, n);
797 }
798 }
799}
800
801// --------------------------------------------------------------------
803 const xercesc::DOMElement* const propertyElement, G4Material* material)
804{
805 G4String name;
806 G4String ref;
807 G4GDMLMatrix matrix;
808
809 const xercesc::DOMNamedNodeMap* const attributes =
810 propertyElement->getAttributes();
811 XMLSize_t attributeCount = attributes->getLength();
812
813 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
814 ++attribute_index)
815 {
816 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
817
818 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
819 {
820 continue;
821 }
822
823 const xercesc::DOMAttr* const attribute =
824 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
825 if(attribute == nullptr)
826 {
827 G4Exception("G4GDMLReadMaterials::PropertyRead()", "InvalidRead",
828 FatalException, "No attribute found!");
829 return;
830 }
831 const G4String attName = Transcode(attribute->getName());
832 const G4String attValue = Transcode(attribute->getValue());
833
834 if(attName == "name")
835 {
836 name = GenerateName(attValue);
837 }
838 else if(attName == "ref")
839 {
840 matrix = GetMatrix(ref = attValue);
841 }
842 }
843
844 /*
845 if (matrix.GetCols() != 2)
846 {
847 G4String error_msg = "Referenced matrix '" + ref
848 + "' should have \n two columns as a property table for material: "
849 + material->GetName();
850 G4Exception("G4GDMLReadMaterials::PropertyRead()", "InvalidRead",
851 FatalException, error_msg);
852 }
853 */
854
855 if(matrix.GetRows() == 0)
856 {
857 return;
858 }
859
861 if(matprop == nullptr)
862 {
863 matprop = new G4MaterialPropertiesTable();
864 material->SetMaterialPropertiesTable(matprop);
865 }
866 if(matrix.GetCols() == 1) // constant property assumed
867 {
868 matprop->AddConstProperty(Strip(name), matrix.Get(0, 0), true);
869 }
870 else // build the material properties vector
871 {
873 for(std::size_t i = 0; i < matrix.GetRows(); ++i)
874 {
875 propvect->InsertValues(matrix.Get(i, 0), matrix.Get(i, 1));
876 }
877 matprop->AddProperty(Strip(name), propvect, true);
878 }
879}
880
881// --------------------------------------------------------------------
883 const xercesc::DOMElement* const materialsElement)
884{
885#ifdef G4VERBOSE
886 G4cout << "G4GDML: Reading materials..." << G4endl;
887#endif
888 for(xercesc::DOMNode* iter = materialsElement->getFirstChild();
889 iter != nullptr; iter = iter->getNextSibling())
890 {
891 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
892 {
893 continue;
894 }
895
896 const xercesc::DOMElement* const child =
897 dynamic_cast<xercesc::DOMElement*>(iter);
898 if(child == nullptr)
899 {
900 G4Exception("G4GDMLReadMaterials::MaterialsRead()", "InvalidRead",
901 FatalException, "No child found!");
902 return;
903 }
904 const G4String tag = Transcode(child->getTagName());
905
906 if(tag == "define")
907 {
908 DefineRead(child);
909 }
910 else if(tag == "element")
911 {
912 ElementRead(child);
913 }
914 else if(tag == "isotope")
915 {
916 IsotopeRead(child);
917 }
918 else if(tag == "material")
919 {
920 MaterialRead(child);
921 }
922 else
923 {
924 G4String error_msg = "Unknown tag in materials: " + tag;
925 G4Exception("G4GDMLReadMaterials::MaterialsRead()", "InvalidSetup",
926 FatalException, error_msg);
927 }
928 }
929}
930
931// --------------------------------------------------------------------
933 G4bool verbose) const
934{
935 G4Element* elementPtr = G4Element::GetElement(ref, false);
936
937 if(elementPtr == nullptr)
938 {
939 elementPtr = G4NistManager::Instance()->FindOrBuildElement(ref);
940 }
941
942 if(verbose && elementPtr == nullptr)
943 {
944 G4String error_msg = "Referenced element '" + ref + "' was not found!";
945 G4Exception("G4GDMLReadMaterials::GetElement()", "InvalidRead",
946 FatalException, error_msg);
947 }
948
949 return elementPtr;
950}
951
952// --------------------------------------------------------------------
954 G4bool verbose) const
955{
956 G4Isotope* isotopePtr = G4Isotope::GetIsotope(ref, false);
957
958 if(verbose && isotopePtr == nullptr)
959 {
960 G4String error_msg = "Referenced isotope '" + ref + "' was not found!";
961 G4Exception("G4GDMLReadMaterials::GetIsotope()", "InvalidRead",
962 FatalException, error_msg);
963 }
964
965 return isotopePtr;
966}
967
968// --------------------------------------------------------------------
970 G4bool verbose) const
971{
972 G4Material* materialPtr = G4Material::GetMaterial(ref, false);
973
974 if(materialPtr == nullptr)
975 {
976 materialPtr = G4NistManager::Instance()->FindOrBuildMaterial(ref);
977 }
978
979 if(verbose && materialPtr == nullptr)
980 {
981 G4String error_msg = "Referenced material '" + ref + "' was not found!";
982 G4Exception("G4GDMLReadMaterials::GetMaterial()", "InvalidRead",
983 FatalException, error_msg);
984 }
985
986 return materialPtr;
987}
G4double D(G4double temp)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4PhysicsFreeVector G4MaterialPropertyVector
G4State
Definition: G4Material.hh:110
@ kStateSolid
Definition: G4Material.hh:110
@ kStateLiquid
Definition: G4Material.hh:110
@ kStateGas
Definition: G4Material.hh:110
@ kStateUndefined
Definition: G4Material.hh:110
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void AddIsotope(G4Isotope *isotope, G4double RelativeAbundance)
Definition: G4Element.cc:154
static G4Element * GetElement(const G4String &name, G4bool warning=true)
Definition: G4Element.cc:417
G4double Evaluate(const G4String &)
G4int EvaluateInteger(const G4String &)
std::size_t GetCols() const
std::size_t GetRows() const
G4double Get(std::size_t r, std::size_t c) const
G4double GetQuantity(const G4String &)
G4GDMLMatrix GetMatrix(const G4String &)
G4String RefRead(const xercesc::DOMElement *const)
virtual void DefineRead(const xercesc::DOMElement *const)
void PropertyRead(const xercesc::DOMElement *const, G4Material *)
G4double MEERead(const xercesc::DOMElement *const)
void MixtureRead(const xercesc::DOMElement *const, G4Element *)
G4double DRead(const xercesc::DOMElement *const)
G4double AtomRead(const xercesc::DOMElement *const)
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
G4Element * GetElement(const G4String &, G4bool verbose=true) const
G4double PRead(const xercesc::DOMElement *const)
G4Isotope * GetIsotope(const G4String &, G4bool verbose=true) const
G4int CompositeRead(const xercesc::DOMElement *const, G4String &)
virtual void MaterialsRead(const xercesc::DOMElement *const)
G4double TRead(const xercesc::DOMElement *const)
G4double FractionRead(const xercesc::DOMElement *const, G4String &)
void ElementRead(const xercesc::DOMElement *const)
void IsotopeRead(const xercesc::DOMElement *const)
void MaterialRead(const xercesc::DOMElement *const)
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:104
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:70
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:156
void SetMeanExcitationEnergy(G4double value)
static G4Isotope * GetIsotope(const G4String &name, G4bool warning=false)
Definition: G4Isotope.cc:197
void AddConstProperty(const G4String &key, G4double propertyValue, G4bool createNewKey=false)
G4MaterialPropertyVector * AddProperty(const G4String &key, const std::vector< G4double > &photonEnergies, const std::vector< G4double > &propertyValues, G4bool createNewKey=false, G4bool spline=false)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:221
void AddElement(G4Element *elm, G4int nAtoms)
Definition: G4Material.hh:156
void AddMaterial(G4Material *material, G4double fraction)
Definition: G4Material.cc:515
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
Definition: G4Material.cc:849
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
G4Element * FindOrBuildElement(G4int Z, G4bool isotopes=true)
void InsertValues(const G4double energy, const G4double value)
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)
#define N
Definition: crc32.c:56
Definition: xmlparse.c:284