Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLReadStructure.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// G4GDMLReadStructure implementation
27//
28// Author: Zoltan Torzsok, November 2007
29// --------------------------------------------------------------------
30
32
33#include "G4UnitsTable.hh"
34#include "G4LogicalVolume.hh"
35#include "G4VPhysicalVolume.hh"
36#include "G4PVPlacement.hh"
39#include "G4AssemblyVolume.hh"
44#include "G4VisAttributes.hh"
45
46// --------------------------------------------------------------------
49{
50}
51
52// --------------------------------------------------------------------
54{
55}
56
57// --------------------------------------------------------------------
59 const xercesc::DOMElement* const bordersurfaceElement)
60{
61 G4String name;
62 G4VPhysicalVolume* pv1 = nullptr;
63 G4VPhysicalVolume* pv2 = nullptr;
64 G4SurfaceProperty* prop = nullptr;
65 G4int index = 0;
66
67 const xercesc::DOMNamedNodeMap* const attributes =
68 bordersurfaceElement->getAttributes();
69 XMLSize_t attributeCount = attributes->getLength();
70
71 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
72 ++attribute_index)
73 {
74 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
75
76 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
77 {
78 continue;
79 }
80
81 const xercesc::DOMAttr* const attribute =
82 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
83 if(attribute == nullptr)
84 {
85 G4Exception("G4GDMLReadStructure::BorderSurfaceRead()", "InvalidRead",
86 FatalException, "No attribute found!");
87 return;
88 }
89 const G4String attName = Transcode(attribute->getName());
90 const G4String attValue = Transcode(attribute->getValue());
91
92 if(attName == "name")
93 {
94 name = GenerateName(attValue);
95 }
96 else if(attName == "surfaceproperty")
97 {
98 prop = GetSurfaceProperty(GenerateName(attValue));
99 }
100 }
101
102 for(xercesc::DOMNode* iter = bordersurfaceElement->getFirstChild();
103 iter != nullptr; iter = iter->getNextSibling())
104 {
105 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
106 {
107 continue;
108 }
109
110 const xercesc::DOMElement* const child =
111 dynamic_cast<xercesc::DOMElement*>(iter);
112 if(child == nullptr)
113 {
114 G4Exception("G4GDMLReadStructure::BorderSurfaceRead()", "InvalidRead",
115 FatalException, "No child found!");
116 return;
117 }
118 const G4String tag = Transcode(child->getTagName());
119
120 if(tag != "physvolref")
121 {
122 continue;
123 }
124
125 if(index == 0)
126 {
127 pv1 = GetPhysvol(GenerateName(RefRead(child)));
128 ++index;
129 }
130 else if(index == 1)
131 {
132 pv2 = GetPhysvol(GenerateName(RefRead(child)));
133 ++index;
134 }
135 else
136 break;
137 }
138
139 new G4LogicalBorderSurface(Strip(name), pv1, pv2, prop);
140}
141
142// --------------------------------------------------------------------
144 const xercesc::DOMElement* const divisionvolElement)
145{
146 G4String name;
147 G4double unit = 1.0;
148 G4double width = 0.0;
149 G4double offset = 0.0;
150 G4int number = 0;
151 EAxis axis = kUndefined;
152 G4LogicalVolume* logvol = nullptr;
153
154 const xercesc::DOMNamedNodeMap* const attributes =
155 divisionvolElement->getAttributes();
156 XMLSize_t attributeCount = attributes->getLength();
157 G4String unitname;
158
159 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
160 ++attribute_index)
161 {
162 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
163
164 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
165 {
166 continue;
167 }
168
169 const xercesc::DOMAttr* const attribute =
170 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
171 if(attribute == nullptr)
172 {
173 G4Exception("G4GDMLReadStructure::DivisionvolRead()", "InvalidRead",
174 FatalException, "No attribute found!");
175 return;
176 }
177 const G4String attName = Transcode(attribute->getName());
178 const G4String attValue = Transcode(attribute->getValue());
179
180 if(attName == "name")
181 {
182 name = attValue;
183 }
184 else if(attName == "unit")
185 {
186 unit = G4UnitDefinition::GetValueOf(attValue);
187 unitname = G4UnitDefinition::GetCategory(attValue);
188 }
189 else if(attName == "width")
190 {
191 width = eval.Evaluate(attValue);
192 }
193 else if(attName == "offset")
194 {
195 offset = eval.Evaluate(attValue);
196 }
197 else if(attName == "number")
198 {
199 number = eval.EvaluateInteger(attValue);
200 }
201 else if(attName == "axis")
202 {
203 if(attValue == "kXAxis")
204 {
205 axis = kXAxis;
206 }
207 else if(attValue == "kYAxis")
208 {
209 axis = kYAxis;
210 }
211 else if(attValue == "kZAxis")
212 {
213 axis = kZAxis;
214 }
215 else if(attValue == "kRho")
216 {
217 axis = kRho;
218 }
219 else if(attValue == "kPhi")
220 {
221 axis = kPhi;
222 }
223 }
224 }
225
226 if(((axis == kXAxis || axis == kYAxis || axis == kZAxis) &&
227 unitname != "Length") ||
228 ((axis == kRho || axis == kPhi) && unitname != "Angle"))
229 {
230 G4Exception("G4GDMLReadStructure::DivisionvolRead()", "InvalidRead",
231 FatalException, "Invalid unit!");
232 }
233
234 width *= unit;
235 offset *= unit;
236
237 for(xercesc::DOMNode* iter = divisionvolElement->getFirstChild();
238 iter != nullptr; iter = iter->getNextSibling())
239 {
240 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
241 {
242 continue;
243 }
244
245 const xercesc::DOMElement* const child =
246 dynamic_cast<xercesc::DOMElement*>(iter);
247 if(child == nullptr)
248 {
249 G4Exception("G4GDMLReadStructure::DivisionvolRead()", "InvalidRead",
250 FatalException, "No child found!");
251 return;
252 }
253 const G4String tag = Transcode(child->getTagName());
254
255 if(tag == "volumeref")
256 {
257 logvol = GetVolume(GenerateName(RefRead(child)));
258 }
259 }
260
261 if(logvol == nullptr)
262 {
263 return;
264 }
265
268
269 G4String pv_name = logvol->GetName() + "_div";
270 if((number != 0) && (width == 0.0))
271 {
273 pv_name, logvol, pMotherLogical, axis, number, offset);
274 }
275 else if((number == 0) && (width != 0.0))
276 {
278 pv_name, logvol, pMotherLogical, axis, width, offset);
279 }
280 else
281 {
283 pv_name, logvol, pMotherLogical, axis, number, width, offset);
284 }
285
286 if(pair.first != nullptr)
287 {
288 GeneratePhysvolName(name, pair.first);
289 }
290 if(pair.second != nullptr)
291 {
292 GeneratePhysvolName(name, pair.second);
293 }
294}
295
296// --------------------------------------------------------------------
298 const xercesc::DOMElement* const fileElement)
299{
300 G4String name;
301 G4String volname;
302
303 const xercesc::DOMNamedNodeMap* const attributes =
304 fileElement->getAttributes();
305 XMLSize_t attributeCount = attributes->getLength();
306
307 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
308 ++attribute_index)
309 {
310 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
311
312 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
313 {
314 continue;
315 }
316
317 const xercesc::DOMAttr* const attribute =
318 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
319 if(attribute == nullptr)
320 {
321 G4Exception("G4GDMLReadStructure::FileRead()", "InvalidRead",
322 FatalException, "No attribute found!");
323 return nullptr;
324 }
325 const G4String attName = Transcode(attribute->getName());
326 const G4String attValue = Transcode(attribute->getValue());
327
328 if(attName == "name")
329 {
330 name = attValue;
331 }
332 else if(attName == "volname")
333 {
334 volname = attValue;
335 }
336 }
337
338 const G4bool isModule = true;
339 G4GDMLReadStructure structure;
340 structure.Read(name, validate, isModule);
341
342 // Register existing auxiliar information defined in child module
343 //
344 const G4GDMLAuxMapType* aux = structure.GetAuxMap();
345 if(!aux->empty())
346 {
347 for(auto pos = aux->cbegin(); pos != aux->cend(); ++pos)
348 {
349 auxMap.insert(std::make_pair(pos->first, pos->second));
350 }
351 }
352
353 // Return volume structure from child module
354 //
355 if(volname.empty())
356 {
357 return structure.GetVolume(structure.GetSetup("Default"));
358 }
359 else
360 {
361 return structure.GetVolume(structure.GenerateName(volname));
362 }
363}
364
365// --------------------------------------------------------------------
367 const xercesc::DOMElement* const physvolElement, G4AssemblyVolume* pAssembly)
368{
369 G4String name;
370 G4LogicalVolume* logvol = nullptr;
371 G4AssemblyVolume* assembly = nullptr;
372 G4ThreeVector position(0.0, 0.0, 0.0);
373 G4ThreeVector rotation(0.0, 0.0, 0.0);
374 G4ThreeVector scale(1.0, 1.0, 1.0);
375 G4int copynumber = 0;
376
377 const xercesc::DOMNamedNodeMap* const attributes =
378 physvolElement->getAttributes();
379 XMLSize_t attributeCount = attributes->getLength();
380
381 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
382 ++attribute_index)
383 {
384 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
385
386 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
387 {
388 continue;
389 }
390
391 const xercesc::DOMAttr* const attribute =
392 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
393 if(attribute == nullptr)
394 {
395 G4Exception("G4GDMLReadStructure::PhysvolRead()", "InvalidRead",
396 FatalException, "No attribute found!");
397 return;
398 }
399 const G4String attName = Transcode(attribute->getName());
400 const G4String attValue = Transcode(attribute->getValue());
401
402 if(attName == "name")
403 {
404 name = attValue;
405 }
406 if(attName == "copynumber")
407 {
408 copynumber = eval.EvaluateInteger(attValue);
409 }
410 }
411
412 for(xercesc::DOMNode* iter = physvolElement->getFirstChild(); iter != nullptr;
413 iter = iter->getNextSibling())
414 {
415 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
416 {
417 continue;
418 }
419
420 const xercesc::DOMElement* const child =
421 dynamic_cast<xercesc::DOMElement*>(iter);
422 if(child == nullptr)
423 {
424 G4Exception("G4GDMLReadStructure::PhysvolRead()", "InvalidRead",
425 FatalException, "No child found!");
426 return;
427 }
428 const G4String tag = Transcode(child->getTagName());
429
430 if(tag == "volumeref")
431 {
432 const G4String& child_name = GenerateName(RefRead(child));
433 assembly = GetAssembly(child_name);
434 if(assembly == nullptr)
435 {
436 logvol = GetVolume(child_name);
437 }
438 }
439 else if(tag == "file")
440 {
441 logvol = FileRead(child);
442 }
443 else if(tag == "position")
444 {
445 VectorRead(child, position);
446 }
447 else if(tag == "rotation")
448 {
449 VectorRead(child, rotation);
450 }
451 else if(tag == "scale")
452 {
453 VectorRead(child, scale);
454 }
455 else if(tag == "positionref")
456 {
458 }
459 else if(tag == "rotationref")
460 {
461 rotation = GetRotation(GenerateName(RefRead(child)));
462 }
463 else if(tag == "scaleref")
464 {
465 scale = GetScale(GenerateName(RefRead(child)));
466 }
467 else
468 {
469 G4String error_msg = "Unknown tag in physvol: " + tag;
470 G4Exception("G4GDMLReadStructure::PhysvolRead()", "ReadError",
471 FatalException, error_msg);
472 return;
473 }
474 }
475
476 G4Transform3D transform(GetRotationMatrix(rotation).inverse(), position);
477 transform = transform * G4Scale3D(scale.x(), scale.y(), scale.z());
478
479 if(pAssembly != nullptr) // Fill assembly structure
480 {
481 if(assembly != nullptr) // Case of recursive assemblies
482 {
483 pAssembly->AddPlacedAssembly(assembly, transform);
484 }
485 if(logvol == nullptr)
486 {
487 return;
488 }
489 pAssembly->AddPlacedVolume(logvol, transform);
490 }
491 else // Generate physical volume tree or do assembly imprint
492 {
493 if(assembly != nullptr)
494 {
495 assembly->MakeImprint(pMotherLogical, transform, 0, check);
496 }
497 else
498 {
499 if(logvol == nullptr)
500 {
501 return;
502 }
503 G4String pv_name = logvol->GetName() + "_PV";
505 transform, pv_name, logvol, pMotherLogical, false, copynumber, check);
506
507 if(pair.first != nullptr)
508 {
509 GeneratePhysvolName(name, pair.first);
510 }
511 if(pair.second != nullptr)
512 {
513 GeneratePhysvolName(name, pair.second);
514 }
515 }
516 }
517}
518
519// --------------------------------------------------------------------
521 const xercesc::DOMElement* const replicavolElement, G4int number)
522{
523 G4LogicalVolume* logvol = nullptr;
524 for(xercesc::DOMNode* iter = replicavolElement->getFirstChild();
525 iter != nullptr; iter = iter->getNextSibling())
526 {
527 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
528 {
529 continue;
530 }
531
532 const xercesc::DOMElement* const child =
533 dynamic_cast<xercesc::DOMElement*>(iter);
534 if(child == nullptr)
535 {
536 G4Exception("G4GDMLReadStructure::ReplicavolRead()", "InvalidRead",
537 FatalException, "No child found!");
538 return;
539 }
540 const G4String tag = Transcode(child->getTagName());
541
542 if(tag == "volumeref")
543 {
544 logvol = GetVolume(GenerateName(RefRead(child)));
545 }
546 else if(tag == "replicate_along_axis")
547 {
548 if(logvol)
549 {
550 ReplicaRead(child, logvol, number);
551 }
552 }
553 else
554 {
555 G4String error_msg = "Unknown tag in ReplicavolRead: " + tag;
556 G4Exception("G4GDMLReadStructure::ReplicavolRead()", "ReadError",
557 FatalException, error_msg);
558 }
559 }
560}
561
562// --------------------------------------------------------------------
564 const xercesc::DOMElement* const replicaElement, G4LogicalVolume* logvol,
565 G4int number)
566{
567 G4double width = 0.0;
568 G4double offset = 0.0;
569 G4ThreeVector position(0.0, 0.0, 0.0);
570 G4ThreeVector rotation(0.0, 0.0, 0.0);
571 EAxis axis = kUndefined;
572 G4String name;
573
574 for(xercesc::DOMNode* iter = replicaElement->getFirstChild(); iter != nullptr;
575 iter = iter->getNextSibling())
576 {
577 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
578 {
579 continue;
580 }
581
582 const xercesc::DOMElement* const child =
583 dynamic_cast<xercesc::DOMElement*>(iter);
584 if(child == nullptr)
585 {
586 G4Exception("G4GDMLReadStructure::ReplicaRead()", "InvalidRead",
587 FatalException, "No child found!");
588 return;
589 }
590 const G4String tag = Transcode(child->getTagName());
591
592 if(tag == "position")
593 {
594 VectorRead(child, position);
595 }
596 else if(tag == "rotation")
597 {
598 VectorRead(child, rotation);
599 }
600 else if(tag == "positionref")
601 {
603 }
604 else if(tag == "rotationref")
605 {
606 rotation = GetRotation(GenerateName(RefRead(child)));
607 }
608 else if(tag == "direction")
609 {
610 axis = AxisRead(child);
611 }
612 else if(tag == "width")
613 {
614 width = QuantityRead(child);
615 }
616 else if(tag == "offset")
617 {
618 offset = QuantityRead(child);
619 }
620 else
621 {
622 G4String error_msg = "Unknown tag in ReplicaRead: " + tag;
623 G4Exception("G4GDMLReadStructure::ReplicaRead()", "ReadError",
624 FatalException, error_msg);
625 }
626 }
627
628 G4String pv_name = logvol->GetName() + "_PV";
630 pv_name, logvol, pMotherLogical, axis, number, width, offset);
631
632 if(pair.first != nullptr)
633 {
634 GeneratePhysvolName(name, pair.first);
635 }
636 if(pair.second != nullptr)
637 {
638 GeneratePhysvolName(name, pair.second);
639 }
640}
641
642// --------------------------------------------------------------------
644 const xercesc::DOMElement* const axisElement)
645{
646 EAxis axis = kUndefined;
647
648 const xercesc::DOMNamedNodeMap* const attributes =
649 axisElement->getAttributes();
650 XMLSize_t attributeCount = attributes->getLength();
651
652 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
653 ++attribute_index)
654 {
655 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
656
657 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
658 {
659 continue;
660 }
661
662 const xercesc::DOMAttr* const attribute =
663 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
664 if(attribute == nullptr)
665 {
666 G4Exception("G4GDMLReadStructure::AxisRead()", "InvalidRead",
667 FatalException, "No attribute found!");
668 return axis;
669 }
670 const G4String attName = Transcode(attribute->getName());
671 const G4String attValue = Transcode(attribute->getValue());
672 if(attName == "x")
673 {
674 if(eval.Evaluate(attValue) == 1.)
675 {
676 axis = kXAxis;
677 }
678 }
679 else if(attName == "y")
680 {
681 if(eval.Evaluate(attValue) == 1.)
682 {
683 axis = kYAxis;
684 }
685 }
686 else if(attName == "z")
687 {
688 if(eval.Evaluate(attValue) == 1.)
689 {
690 axis = kZAxis;
691 }
692 }
693 else if(attName == "rho")
694 {
695 if(eval.Evaluate(attValue) == 1.)
696 {
697 axis = kRho;
698 }
699 }
700 else if(attName == "phi")
701 {
702 if(eval.Evaluate(attValue) == 1.)
703 {
704 axis = kPhi;
705 }
706 }
707 }
708
709 return axis;
710}
711
712// --------------------------------------------------------------------
714 const xercesc::DOMElement* const readElement)
715{
716 G4double value = 0.0;
717 G4double unit = 0.0;
718 const xercesc::DOMNamedNodeMap* const attributes =
719 readElement->getAttributes();
720 XMLSize_t attributeCount = attributes->getLength();
721
722 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
723 ++attribute_index)
724 {
725 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
726
727 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
728 {
729 continue;
730 }
731 const xercesc::DOMAttr* const attribute =
732 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
733 if(attribute == nullptr)
734 {
735 G4Exception("G4GDMLReadStructure::QuantityRead()", "InvalidRead",
736 FatalException, "No attribute found!");
737 return value;
738 }
739 const G4String attName = Transcode(attribute->getName());
740 const G4String attValue = Transcode(attribute->getValue());
741
742 if(attName == "unit")
743 {
744 unit = G4UnitDefinition::GetValueOf(attValue);
745 if(G4UnitDefinition::GetCategory(attValue) != "Length" &&
746 G4UnitDefinition::GetCategory(attValue) != "Angle")
747 {
748 G4Exception("G4GDMLReadStructure::QuantityRead()", "InvalidRead",
750 "Invalid unit for length or angle (width, offset)!");
751 }
752 }
753 else if(attName == "value")
754 {
755 value = eval.Evaluate(attValue);
756 }
757 }
758
759 return value * unit;
760}
761
762// --------------------------------------------------------------------
764 const xercesc::DOMElement* const volumeElement)
765{
766 G4VSolid* solidPtr = nullptr;
767 G4Material* materialPtr = nullptr;
768 G4GDMLAuxListType auxList;
769
770 XMLCh* name_attr = xercesc::XMLString::transcode("name");
771 const G4String name = Transcode(volumeElement->getAttribute(name_attr));
772 xercesc::XMLString::release(&name_attr);
773
774 for(xercesc::DOMNode* iter = volumeElement->getFirstChild(); iter != nullptr;
775 iter = iter->getNextSibling())
776 {
777 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
778 {
779 continue;
780 }
781
782 const xercesc::DOMElement* const child =
783 dynamic_cast<xercesc::DOMElement*>(iter);
784 if(child == nullptr)
785 {
786 G4Exception("G4GDMLReadStructure::VolumeRead()", "InvalidRead",
787 FatalException, "No child found!");
788 return;
789 }
790 const G4String tag = Transcode(child->getTagName());
791
792 if(tag == "auxiliary")
793 {
794 auxList.push_back(AuxiliaryRead(child));
795 }
796 else if(tag == "materialref")
797 {
798 materialPtr = GetMaterial(GenerateName(RefRead(child), true));
799 }
800 else if(tag == "solidref")
801 {
802 solidPtr = GetSolid(GenerateName(RefRead(child)));
803 }
804 }
805
807 new G4LogicalVolume(solidPtr, materialPtr, GenerateName(name), 0, 0, 0);
808
809 if(!auxList.empty())
810 {
811 auxMap[pMotherLogical] = auxList;
812 }
813
814 Volume_contentRead(volumeElement);
815}
816
817// --------------------------------------------------------------------
819 const xercesc::DOMElement* const assemblyElement)
820{
821 XMLCh* name_attr = xercesc::XMLString::transcode("name");
822 const G4String name = Transcode(assemblyElement->getAttribute(name_attr));
823 xercesc::XMLString::release(&name_attr);
824
825 G4AssemblyVolume* pAssembly = new G4AssemblyVolume();
826 assemblyMap.insert(std::make_pair(GenerateName(name), pAssembly));
827
828 for(xercesc::DOMNode* iter = assemblyElement->getFirstChild();
829 iter != nullptr; iter = iter->getNextSibling())
830 {
831 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
832 {
833 continue;
834 }
835 const xercesc::DOMElement* const child =
836 dynamic_cast<xercesc::DOMElement*>(iter);
837 if(child == nullptr)
838 {
839 G4Exception("G4GDMLReadStructure::AssemblyRead()", "InvalidRead",
840 FatalException, "No child found!");
841 return;
842 }
843 const G4String tag = Transcode(child->getTagName());
844
845 if(tag == "physvol")
846 {
847 PhysvolRead(child, pAssembly);
848 }
849 else
850 {
851 G4cout << "Unsupported GDML tag '" << tag
852 << "' for Geant4 assembly structure !" << G4endl;
853 }
854 }
855}
856
857// --------------------------------------------------------------------
859 const xercesc::DOMElement* const skinsurfaceElement)
860{
861 G4String name;
862 G4LogicalVolume* logvol = nullptr;
863 G4SurfaceProperty* prop = nullptr;
864
865 const xercesc::DOMNamedNodeMap* const attributes =
866 skinsurfaceElement->getAttributes();
867 XMLSize_t attributeCount = attributes->getLength();
868
869 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
870 ++attribute_index)
871 {
872 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
873
874 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
875 {
876 continue;
877 }
878
879 const xercesc::DOMAttr* const attribute =
880 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
881 if(attribute == nullptr)
882 {
883 G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "InvalidRead",
884 FatalException, "No attribute found!");
885 return;
886 }
887 const G4String attName = Transcode(attribute->getName());
888 const G4String attValue = Transcode(attribute->getValue());
889
890 if(attName == "name")
891 {
892 name = GenerateName(attValue);
893 }
894 else if(attName == "surfaceproperty")
895 {
896 prop = GetSurfaceProperty(GenerateName(attValue));
897 }
898 }
899
900 for(xercesc::DOMNode* iter = skinsurfaceElement->getFirstChild();
901 iter != nullptr; iter = iter->getNextSibling())
902 {
903 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
904 {
905 continue;
906 }
907
908 const xercesc::DOMElement* const child =
909 dynamic_cast<xercesc::DOMElement*>(iter);
910 if(child == nullptr)
911 {
912 G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "InvalidRead",
913 FatalException, "No child found!");
914 return;
915 }
916 const G4String tag = Transcode(child->getTagName());
917
918 if(tag == "volumeref")
919 {
920 logvol = GetVolume(GenerateName(RefRead(child)));
921 }
922 else
923 {
924 G4String error_msg = "Unknown tag in skinsurface: " + tag;
925 G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "ReadError",
926 FatalException, error_msg);
927 }
928 }
929
930 new G4LogicalSkinSurface(Strip(name), logvol, prop);
931}
932
933// --------------------------------------------------------------------
935 const xercesc::DOMElement* const volumeElement)
936{
937 for(xercesc::DOMNode* iter = volumeElement->getFirstChild(); iter != nullptr;
938 iter = iter->getNextSibling())
939 {
940 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
941 {
942 continue;
943 }
944
945 const xercesc::DOMElement* const child =
946 dynamic_cast<xercesc::DOMElement*>(iter);
947 if(child == nullptr)
948 {
949 G4Exception("G4GDMLReadStructure::Volume_contentRead()", "InvalidRead",
950 FatalException, "No child found!");
951 return;
952 }
953 const G4String tag = Transcode(child->getTagName());
954
955 if((tag == "auxiliary") || (tag == "materialref") || (tag == "solidref"))
956 {
957 // These are already processed in VolumeRead()
958 }
959 else if(tag == "paramvol")
960 {
962 }
963 else if(tag == "physvol")
964 {
965 PhysvolRead(child);
966 }
967 else if(tag == "replicavol")
968 {
969 G4int number = 1;
970 const xercesc::DOMNamedNodeMap* const attributes = child->getAttributes();
971 XMLSize_t attributeCount = attributes->getLength();
972 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
973 ++attribute_index)
974 {
975 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
976 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
977 {
978 continue;
979 }
980 const xercesc::DOMAttr* const attribute =
981 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
982 if(attribute == nullptr)
983 {
984 G4Exception("G4GDMLReadStructure::Volume_contentRead()",
985 "InvalidRead", FatalException, "No attribute found!");
986 return;
987 }
988 const G4String attName = Transcode(attribute->getName());
989 const G4String attValue = Transcode(attribute->getValue());
990 if(attName == "number")
991 {
992 number = eval.EvaluateInteger(attValue);
993 }
994 }
995 ReplicavolRead(child, number);
996 }
997 else if(tag == "divisionvol")
998 {
999 DivisionvolRead(child);
1000 }
1001 else if(tag == "loop")
1002 {
1004 }
1005 else
1006 {
1007 G4cout << "Treating unknown GDML tag in volume '" << tag
1008 << "' as GDML extension..." << G4endl;
1009 }
1010 }
1011}
1012
1013// --------------------------------------------------------------------
1015 const xercesc::DOMElement* const structureElement)
1016{
1017#ifdef G4VERBOSE
1018 G4cout << "G4GDML: Reading structure..." << G4endl;
1019#endif
1020 for(xercesc::DOMNode* iter = structureElement->getFirstChild();
1021 iter != nullptr; iter = iter->getNextSibling())
1022 {
1023 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
1024 {
1025 continue;
1026 }
1027
1028 const xercesc::DOMElement* const child =
1029 dynamic_cast<xercesc::DOMElement*>(iter);
1030 if(child == nullptr)
1031 {
1032 G4Exception("G4GDMLReadStructure::StructureRead()", "InvalidRead",
1033 FatalException, "No child found!");
1034 return;
1035 }
1036 const G4String tag = Transcode(child->getTagName());
1037
1038 if(tag == "bordersurface")
1039 {
1040 BorderSurfaceRead(child);
1041 }
1042 else if(tag == "skinsurface")
1043 {
1044 SkinSurfaceRead(child);
1045 }
1046 else if(tag == "volume")
1047 {
1048 VolumeRead(child);
1049 }
1050 else if(tag == "assembly")
1051 {
1052 AssemblyRead(child);
1053 }
1054 else if(tag == "loop")
1055 {
1057 }
1058 else
1059 {
1060 G4String error_msg = "Unknown tag in structure: " + tag;
1061 G4Exception("G4GDMLReadStructure::StructureRead()", "ReadError",
1062 FatalException, error_msg);
1063 }
1064 }
1065}
1066
1067// --------------------------------------------------------------------
1069{
1070 G4VPhysicalVolume* physvolPtr =
1072
1073 if(physvolPtr == nullptr)
1074 {
1075 G4String error_msg = "Referenced physvol '" + ref + "' was not found!";
1076 G4Exception("G4GDMLReadStructure::GetPhysvol()", "ReadError",
1077 FatalException, error_msg);
1078 }
1079
1080 return physvolPtr;
1081}
1082
1083// --------------------------------------------------------------------
1085{
1086 G4LogicalVolume* volumePtr =
1088
1089 if(volumePtr == nullptr)
1090 {
1091 G4String error_msg = "Referenced volume '" + ref + "' was not found!";
1092 G4Exception("G4GDMLReadStructure::GetVolume()", "ReadError", FatalException,
1093 error_msg);
1094 }
1095
1096 return volumePtr;
1097}
1098
1099// --------------------------------------------------------------------
1101{
1102 auto pos = assemblyMap.find(ref);
1103 if(pos != assemblyMap.cend())
1104 {
1105 return pos->second;
1106 }
1107 return nullptr;
1108}
1109
1110// --------------------------------------------------------------------
1112 G4LogicalVolume* logvol) const
1113{
1114 auto pos = auxMap.find(logvol);
1115 if(pos != auxMap.cend())
1116 {
1117 return pos->second;
1118 }
1119 else
1120 {
1121 return G4GDMLAuxListType();
1122 }
1123}
1124
1125// --------------------------------------------------------------------
1127 const G4String& setupName)
1128{
1129 G4String sname = GetSetup(setupName);
1130 if(sname == "")
1131 {
1132 return nullptr;
1133 }
1134
1135 G4LogicalVolume* volume = GetVolume(GenerateName(sname, dostrip));
1137
1138 G4VPhysicalVolume* pvWorld = nullptr;
1139
1140 if(setuptoPV[setupName])
1141 {
1142 pvWorld = setuptoPV[setupName];
1143 }
1144 else
1145 {
1146 pvWorld = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, 0), volume,
1147 volume->GetName() + "_PV", 0, 0, 0);
1148 setuptoPV[setupName] = pvWorld;
1149 }
1150 return pvWorld;
1151}
1152
1153// --------------------------------------------------------------------
1155{
1156 eval.Clear();
1157 setuptoPV.clear();
1158 auxMap.clear();
1159}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
std::map< G4LogicalVolume *, G4GDMLAuxListType > G4GDMLAuxMapType
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Scale3D G4Scale3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
void MakeImprint(G4LogicalVolume *pMotherLV, G4ThreeVector &translationInMother, G4RotationMatrix *pRotationInMother, G4int copyNumBase=0, G4bool surfCheck=false)
void AddPlacedAssembly(G4AssemblyVolume *pAssembly, G4Transform3D &transformation)
void AddPlacedVolume(G4LogicalVolume *pPlacedVolume, G4ThreeVector &translation, G4RotationMatrix *rotation)
G4double Evaluate(const G4String &)
G4int EvaluateInteger(const G4String &)
G4ThreeVector GetScale(const G4String &)
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
G4RotationMatrix GetRotationMatrix(const G4ThreeVector &)
G4ThreeVector GetPosition(const G4String &)
G4String RefRead(const xercesc::DOMElement *const)
G4ThreeVector GetRotation(const G4String &)
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
virtual void ParamvolRead(const xercesc::DOMElement *const, G4LogicalVolume *)
G4String GetSetup(const G4String &)
G4SurfaceProperty * GetSurfaceProperty(const G4String &) const
G4VSolid * GetSolid(const G4String &) const
const G4GDMLAuxMapType * GetAuxMap() const
void BorderSurfaceRead(const xercesc::DOMElement *const)
virtual void Volume_contentRead(const xercesc::DOMElement *const)
virtual void StructureRead(const xercesc::DOMElement *const)
void AssemblyRead(const xercesc::DOMElement *const)
G4LogicalVolume * GetVolume(const G4String &) const
void PhysvolRead(const xercesc::DOMElement *const, G4AssemblyVolume *assembly=0)
void ReplicaRead(const xercesc::DOMElement *const replicaElement, G4LogicalVolume *logvol, G4int number)
EAxis AxisRead(const xercesc::DOMElement *const axisElement)
G4VPhysicalVolume * GetWorldVolume(const G4String &)
G4GDMLAuxMapType auxMap
G4double QuantityRead(const xercesc::DOMElement *const readElement)
void DivisionvolRead(const xercesc::DOMElement *const)
std::map< std::string, G4VPhysicalVolume * > setuptoPV
void SkinSurfaceRead(const xercesc::DOMElement *const)
G4LogicalVolume * pMotherLogical
void ReplicavolRead(const xercesc::DOMElement *const, G4int number)
virtual void VolumeRead(const xercesc::DOMElement *const)
G4GDMLAssemblyMapType assemblyMap
G4VPhysicalVolume * GetPhysvol(const G4String &) const
G4AssemblyVolume * GetAssembly(const G4String &) const
G4LogicalVolume * FileRead(const xercesc::DOMElement *const)
G4GDMLAuxListType GetVolumeAuxiliaryInformation(G4LogicalVolume *) const
G4bool check
Definition: G4GDMLRead.hh:158
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:104
G4bool validate
Definition: G4GDMLRead.hh:157
G4GDMLAuxStructType AuxiliaryRead(const xercesc::DOMElement *const auxElem)
Definition: G4GDMLRead.cc:287
G4bool dostrip
Definition: G4GDMLRead.hh:159
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:87
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:70
void LoopRead(const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
Definition: G4GDMLRead.cc:186
virtual void Volume_contentRead(const xercesc::DOMElement *const)=0
void Read(const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
Definition: G4GDMLRead.cc:414
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
virtual void StructureRead(const xercesc::DOMElement *const)=0
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:156
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
static G4LogicalVolumeStore * GetInstance()
void SetVisAttributes(const G4VisAttributes *pVA)
const G4String & GetName() const
static G4PVDivisionFactory * GetInstance()
static G4PhysicalVolumeStore * GetInstance()
G4VPhysicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
G4PhysicalVolumesPair Divide(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofDivisions, G4double width, G4double offset)
static G4ReflectionFactory * Instance()
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)
G4PhysicalVolumesPair Replicate(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofReplicas, G4double width, G4double offset=0.)
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)
static const G4VisAttributes & GetInvisible()
EAxis
Definition: geomdefs.hh:54
@ kPhi
Definition: geomdefs.hh:60
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
@ kUndefined
Definition: geomdefs.hh:61
@ kRho
Definition: geomdefs.hh:58
Definition: xmlparse.cc:187
#define position
Definition: xmlparse.cc:622