Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/persistency/gdml/src/G4GDMLReadStructure.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 31 #include "G4GDMLReadStructure.hh"
 32 
 33 #include "G4UnitsTable.hh"
 34 #include "G4LogicalVolume.hh"
 35 #include "G4VPhysicalVolume.hh"
 36 #include "G4PVPlacement.hh"
 37 #include "G4LogicalVolumeStore.hh"
 38 #include "G4PhysicalVolumeStore.hh"
 39 #include "G4AssemblyVolume.hh"
 40 #include "G4ReflectionFactory.hh"
 41 #include "G4PVDivisionFactory.hh"
 42 #include "G4LogicalBorderSurface.hh"
 43 #include "G4LogicalSkinSurface.hh"
 44 #include "G4VisAttributes.hh"
 45 
 46 // --------------------------------------------------------------------
 47 G4GDMLReadStructure::G4GDMLReadStructure()
 48   : G4GDMLReadParamvol()
 49 {
 50 }
 51 
 52 // --------------------------------------------------------------------
 53 G4GDMLReadStructure::~G4GDMLReadStructure()
 54 {
 55 }
 56 
 57 // --------------------------------------------------------------------
 58 void G4GDMLReadStructure::BorderSurfaceRead(
 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 // --------------------------------------------------------------------
143 void G4GDMLReadStructure::DivisionvolRead(
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 
266   G4PVDivisionFactory::GetInstance();
267   G4PhysicalVolumesPair pair;
268 
269   G4String pv_name = logvol->GetName() + "_div";
270   if((number != 0) && (width == 0.0))
271   {
272     pair = G4ReflectionFactory::Instance()->Divide(
273       pv_name, logvol, pMotherLogical, axis, number, offset);
274   }
275   else if((number == 0) && (width != 0.0))
276   {
277     pair = G4ReflectionFactory::Instance()->Divide(
278       pv_name, logvol, pMotherLogical, axis, width, offset);
279   }
280   else
281   {
282     pair = G4ReflectionFactory::Instance()->Divide(
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 // --------------------------------------------------------------------
297 G4LogicalVolume* G4GDMLReadStructure::FileRead(
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 // --------------------------------------------------------------------
366 void G4GDMLReadStructure::PhysvolRead(
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     {
457       position = GetPosition(GenerateName(RefRead(child)));
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";
504       G4PhysicalVolumesPair pair = G4ReflectionFactory::Instance()->Place(
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 // --------------------------------------------------------------------
520 void G4GDMLReadStructure::ReplicavolRead(
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 // --------------------------------------------------------------------
563 void G4GDMLReadStructure::ReplicaRead(
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     {
602       position = GetPosition(GenerateName(RefRead(child)));
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";
629   G4PhysicalVolumesPair pair = G4ReflectionFactory::Instance()->Replicate(
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 // --------------------------------------------------------------------
643 EAxis G4GDMLReadStructure::AxisRead(
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 // --------------------------------------------------------------------
713 G4double G4GDMLReadStructure::QuantityRead(
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",
749                     FatalException,
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 // --------------------------------------------------------------------
763 void G4GDMLReadStructure::VolumeRead(
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 
806   pMotherLogical =
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 // --------------------------------------------------------------------
818 void G4GDMLReadStructure::AssemblyRead(
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   auto aName = GenerateName(name);
827   if(reverseSearch)
828   {
829     assemblyMap.insert_or_assign(aName, pAssembly);
830   }
831   else
832   { 
833     assemblyMap.insert(std::make_pair(aName, pAssembly));
834   }
835 
836   for(xercesc::DOMNode* iter = assemblyElement->getFirstChild();
837                         iter != nullptr; iter = iter->getNextSibling())
838   {
839     if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
840     {
841       continue;
842     }
843     const xercesc::DOMElement* const child =
844       dynamic_cast<xercesc::DOMElement*>(iter);
845     if(child == nullptr)
846     {
847       G4Exception("G4GDMLReadStructure::AssemblyRead()", "InvalidRead",
848                   FatalException, "No child found!");
849       return;
850     }
851     const G4String tag = Transcode(child->getTagName());
852 
853     if(tag == "physvol")
854     {
855       PhysvolRead(child, pAssembly);
856     }
857     else
858     {
859       G4cout << "Unsupported GDML tag '" << tag
860              << "' for Geant4 assembly structure !" << G4endl;
861     }
862   }
863 }
864 
865 // --------------------------------------------------------------------
866 void G4GDMLReadStructure::SkinSurfaceRead(
867   const xercesc::DOMElement* const skinsurfaceElement)
868 {
869   G4String name;
870   G4LogicalVolume* logvol = nullptr;
871   G4SurfaceProperty* prop = nullptr;
872 
873   const xercesc::DOMNamedNodeMap* const attributes =
874     skinsurfaceElement->getAttributes();
875   XMLSize_t attributeCount = attributes->getLength();
876 
877   for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
878       ++attribute_index)
879   {
880     xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
881 
882     if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
883     {
884       continue;
885     }
886 
887     const xercesc::DOMAttr* const attribute =
888       dynamic_cast<xercesc::DOMAttr*>(attribute_node);
889     if(attribute == nullptr)
890     {
891       G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "InvalidRead",
892                   FatalException, "No attribute found!");
893       return;
894     }
895     const G4String attName  = Transcode(attribute->getName());
896     const G4String attValue = Transcode(attribute->getValue());
897 
898     if(attName == "name")
899     {
900       name = GenerateName(attValue);
901     }
902     else if(attName == "surfaceproperty")
903     {
904       prop = GetSurfaceProperty(GenerateName(attValue));
905     }
906   }
907 
908   for(xercesc::DOMNode* iter = skinsurfaceElement->getFirstChild();
909                         iter != nullptr; iter = iter->getNextSibling())
910   {
911     if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
912     {
913       continue;
914     }
915 
916     const xercesc::DOMElement* const child =
917       dynamic_cast<xercesc::DOMElement*>(iter);
918     if(child == nullptr)
919     {
920       G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "InvalidRead",
921                   FatalException, "No child found!");
922       return;
923     }
924     const G4String tag = Transcode(child->getTagName());
925 
926     if(tag == "volumeref")
927     {
928       logvol = GetVolume(GenerateName(RefRead(child)));
929     }
930     else
931     {
932       G4String error_msg = "Unknown tag in skinsurface: " + tag;
933       G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "ReadError",
934                   FatalException, error_msg);
935     }
936   }
937 
938   new G4LogicalSkinSurface(Strip(name), logvol, prop);
939 }
940 
941 // --------------------------------------------------------------------
942 void G4GDMLReadStructure::Volume_contentRead(
943   const xercesc::DOMElement* const volumeElement)
944 {
945   for(xercesc::DOMNode* iter = volumeElement->getFirstChild(); iter != nullptr;
946       iter                   = iter->getNextSibling())
947   {
948     if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
949     {
950       continue;
951     }
952 
953     const xercesc::DOMElement* const child =
954       dynamic_cast<xercesc::DOMElement*>(iter);
955     if(child == nullptr)
956     {
957       G4Exception("G4GDMLReadStructure::Volume_contentRead()", "InvalidRead",
958                   FatalException, "No child found!");
959       return;
960     }
961     const G4String tag = Transcode(child->getTagName());
962 
963     if((tag == "auxiliary") || (tag == "materialref") || (tag == "solidref"))
964     {
965       // These are already processed in VolumeRead()
966     }
967     else if(tag == "paramvol")
968     {
969       ParamvolRead(child, pMotherLogical);
970     }
971     else if(tag == "physvol")
972     {
973       PhysvolRead(child);
974     }
975     else if(tag == "replicavol")
976     {
977       G4int number = 1;
978       const xercesc::DOMNamedNodeMap* const attributes = child->getAttributes();
979       XMLSize_t attributeCount = attributes->getLength();
980       for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
981           ++attribute_index)
982       {
983         xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
984         if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
985         {
986           continue;
987         }
988         const xercesc::DOMAttr* const attribute =
989           dynamic_cast<xercesc::DOMAttr*>(attribute_node);
990         if(attribute == nullptr)
991         {
992           G4Exception("G4GDMLReadStructure::Volume_contentRead()",
993                       "InvalidRead", FatalException, "No attribute found!");
994           return;
995         }
996         const G4String attName  = Transcode(attribute->getName());
997         const G4String attValue = Transcode(attribute->getValue());
998         if(attName == "number")
999         {
1000           number = eval.EvaluateInteger(attValue);
1001         }
1002       }
1003       ReplicavolRead(child, number);
1004     }
1005     else if(tag == "divisionvol")
1006     {
1007       DivisionvolRead(child);
1008     }
1009     else if(tag == "loop")
1010     {
1011       LoopRead(child, &G4GDMLRead::Volume_contentRead);
1012     }
1013     else
1014     {
1015       G4cout << "Treating unknown GDML tag in volume '" << tag
1016              << "' as GDML extension..." << G4endl;
1017     }
1018   }
1019 }
1020 
1021 // --------------------------------------------------------------------
1022 void G4GDMLReadStructure::StructureRead(
1023   const xercesc::DOMElement* const structureElement)
1024 {
1025 #ifdef G4VERBOSE
1026   G4cout << "G4GDML: Reading structure..." << G4endl;
1027 #endif
1028   for(xercesc::DOMNode* iter = structureElement->getFirstChild();
1029                         iter != nullptr; iter = iter->getNextSibling())
1030   {
1031     if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
1032     {
1033       continue;
1034     }
1035 
1036     const xercesc::DOMElement* const child =
1037       dynamic_cast<xercesc::DOMElement*>(iter);
1038     if(child == nullptr)
1039     {
1040       G4Exception("G4GDMLReadStructure::StructureRead()", "InvalidRead",
1041                   FatalException, "No child found!");
1042       return;
1043     }
1044     const G4String tag = Transcode(child->getTagName());
1045 
1046     if(tag == "bordersurface")
1047     {
1048       BorderSurfaceRead(child);
1049     }
1050     else if(tag == "skinsurface")
1051     {
1052       SkinSurfaceRead(child);
1053     }
1054     else if(tag == "volume")
1055     {
1056       VolumeRead(child);
1057     }
1058     else if(tag == "assembly")
1059     {
1060       AssemblyRead(child);
1061     }
1062     else if(tag == "loop")
1063     {
1064       LoopRead(child, &G4GDMLRead::StructureRead);
1065     }
1066     else
1067     {
1068       G4String error_msg = "Unknown tag in structure: " + tag;
1069       G4Exception("G4GDMLReadStructure::StructureRead()", "ReadError",
1070                   FatalException, error_msg);
1071     }
1072   }
1073 }
1074 
1075 // --------------------------------------------------------------------
1076 G4VPhysicalVolume* G4GDMLReadStructure::GetPhysvol(const G4String& ref) const
1077 {
1078   G4VPhysicalVolume* physvolPtr
1079     = G4PhysicalVolumeStore::GetInstance()->GetVolume(ref,false,reverseSearch);
1080 
1081   if(physvolPtr == nullptr)
1082   {
1083     G4String error_msg = "Referenced physvol '" + ref + "' was not found!";
1084     G4Exception("G4GDMLReadStructure::GetPhysvol()", "ReadError",
1085                 FatalException, error_msg);
1086   }
1087 
1088   return physvolPtr;
1089 }
1090 
1091 // --------------------------------------------------------------------
1092 G4LogicalVolume* G4GDMLReadStructure::GetVolume(const G4String& ref) const
1093 {
1094   G4LogicalVolume* volumePtr
1095     = G4LogicalVolumeStore::GetInstance()->GetVolume(ref,false,reverseSearch);
1096 
1097   if(volumePtr == nullptr)
1098   {
1099     G4String error_msg = "Referenced volume '" + ref + "' was not found!";
1100     G4Exception("G4GDMLReadStructure::GetVolume()", "ReadError", FatalException,
1101                 error_msg);
1102   }
1103 
1104   return volumePtr;
1105 }
1106 
1107 // --------------------------------------------------------------------
1108 G4AssemblyVolume* G4GDMLReadStructure::GetAssembly(const G4String& ref) const
1109 {
1110   auto pos = assemblyMap.find(ref);
1111   if(pos != assemblyMap.cend())
1112   {
1113     return pos->second;
1114   }
1115   return nullptr;
1116 }
1117 
1118 // --------------------------------------------------------------------
1119 G4GDMLAuxListType G4GDMLReadStructure::GetVolumeAuxiliaryInformation(
1120   G4LogicalVolume* logvol) const
1121 {
1122   auto pos = auxMap.find(logvol);
1123   if(pos != auxMap.cend())
1124   {
1125     return pos->second;
1126   }
1127   else
1128   {
1129     return G4GDMLAuxListType();
1130   }
1131 }
1132 
1133 // --------------------------------------------------------------------
1134 G4VPhysicalVolume* G4GDMLReadStructure::GetWorldVolume(
1135   const G4String& setupName)
1136 {
1137   G4String sname = GetSetup(setupName);
1138   if(sname == "")
1139   {
1140     return nullptr;
1141   }
1142 
1143   G4LogicalVolume* volume = GetVolume(GenerateName(sname, dostrip));
1144   volume->SetVisAttributes(G4VisAttributes::GetInvisible());
1145 
1146   G4VPhysicalVolume* pvWorld = nullptr;
1147 
1148   if(setuptoPV[setupName])
1149   {
1150     pvWorld = setuptoPV[setupName];
1151   }
1152   else
1153   {
1154     pvWorld = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, 0), volume,
1155                                 volume->GetName() + "_PV", 0, 0, 0);
1156     setuptoPV[setupName] = pvWorld;
1157   }
1158   return pvWorld;
1159 }
1160 
1161 // --------------------------------------------------------------------
1162 void G4GDMLReadStructure::Clear()
1163 {
1164   eval.Clear();
1165   setuptoPV.clear();
1166   auxMap.clear();
1167 }
1168