Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLReadStructure Class Reference

#include <G4GDMLReadStructure.hh>

+ Inheritance diagram for G4GDMLReadStructure:

Public Member Functions

 G4GDMLReadStructure ()
 
virtual ~G4GDMLReadStructure ()
 
G4VPhysicalVolumeGetPhysvol (const G4String &) const
 
G4LogicalVolumeGetVolume (const G4String &) const
 
G4AssemblyVolumeGetAssembly (const G4String &) const
 
G4GDMLAuxListType GetVolumeAuxiliaryInformation (G4LogicalVolume *) const
 
G4VPhysicalVolumeGetWorldVolume (const G4String &)
 
const G4GDMLAuxMapTypeGetAuxMap () const
 
void Clear ()
 
virtual void VolumeRead (const xercesc::DOMElement *const)
 
virtual void Volume_contentRead (const xercesc::DOMElement *const)
 
virtual void StructureRead (const xercesc::DOMElement *const)
 
- Public Member Functions inherited from G4GDMLReadParamvol
virtual void ParamvolRead (const xercesc::DOMElement *const, G4LogicalVolume *)
 
virtual void Paramvol_contentRead (const xercesc::DOMElement *const)
 
- Public Member Functions inherited from G4GDMLReadSetup
G4String GetSetup (const G4String &)
 
virtual void SetupRead (const xercesc::DOMElement *const element)
 
- Public Member Functions inherited from G4GDMLReadSolids
G4VSolidGetSolid (const G4String &) const
 
G4SurfacePropertyGetSurfaceProperty (const G4String &) const
 
virtual void SolidsRead (const xercesc::DOMElement *const)
 
- Public Member Functions inherited from G4GDMLReadMaterials
G4ElementGetElement (const G4String &, G4bool verbose=true) const
 
G4IsotopeGetIsotope (const G4String &, G4bool verbose=true) const
 
G4MaterialGetMaterial (const G4String &, G4bool verbose=true) const
 
virtual void MaterialsRead (const xercesc::DOMElement *const)
 
- Public Member Functions inherited from G4GDMLReadDefine
G4bool IsValidID (const G4String &) const
 
G4double GetConstant (const G4String &)
 
G4double GetVariable (const G4String &)
 
G4double GetQuantity (const G4String &)
 
G4ThreeVector GetPosition (const G4String &)
 
G4ThreeVector GetRotation (const G4String &)
 
G4ThreeVector GetScale (const G4String &)
 
G4GDMLMatrix GetMatrix (const G4String &)
 
virtual void DefineRead (const xercesc::DOMElement *const)
 
void SetReverseSearch (G4bool flag)
 
- Public Member Functions inherited from G4GDMLRead
virtual void DefineRead (const xercesc::DOMElement *const)=0
 
virtual void MaterialsRead (const xercesc::DOMElement *const)=0
 
virtual void SetupRead (const xercesc::DOMElement *const)=0
 
virtual void SolidsRead (const xercesc::DOMElement *const)=0
 
virtual void Paramvol_contentRead (const xercesc::DOMElement *const)=0
 
virtual void Volume_contentRead (const xercesc::DOMElement *const)=0
 
virtual void StructureRead (const xercesc::DOMElement *const)=0
 
virtual void ExtensionRead (const xercesc::DOMElement *const)
 
virtual void UserinfoRead (const xercesc::DOMElement *const)
 
virtual G4LogicalVolumeGetVolume (const G4String &) const =0
 
virtual G4String GetSetup (const G4String &)=0
 
void Read (const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
 
void StripNames () const
 
void StripName (G4String &) const
 
void OverlapCheck (G4bool)
 
const G4GDMLAuxListTypeGetAuxList () const
 

Protected Member Functions

void AssemblyRead (const xercesc::DOMElement *const)
 
void DivisionvolRead (const xercesc::DOMElement *const)
 
G4LogicalVolumeFileRead (const xercesc::DOMElement *const)
 
void PhysvolRead (const xercesc::DOMElement *const, G4AssemblyVolume *assembly=0)
 
void ReplicavolRead (const xercesc::DOMElement *const, G4int number)
 
void ReplicaRead (const xercesc::DOMElement *const replicaElement, G4LogicalVolume *logvol, G4int number)
 
EAxis AxisRead (const xercesc::DOMElement *const axisElement)
 
G4double QuantityRead (const xercesc::DOMElement *const readElement)
 
void BorderSurfaceRead (const xercesc::DOMElement *const)
 
void SkinSurfaceRead (const xercesc::DOMElement *const)
 
- Protected Member Functions inherited from G4GDMLReadParamvol
 G4GDMLReadParamvol ()
 
virtual ~G4GDMLReadParamvol ()
 
void Box_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Trd_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Trap_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Tube_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Cone_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Sphere_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Orb_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Torus_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Ellipsoid_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Para_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Hype_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Polycone_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void Polyhedra_dimensionsRead (const xercesc::DOMElement *const, G4GDMLParameterisation::PARAMETER &)
 
void ParameterisedRead (const xercesc::DOMElement *const)
 
void ParametersRead (const xercesc::DOMElement *const)
 
- Protected Member Functions inherited from G4GDMLReadSetup
 G4GDMLReadSetup ()
 
virtual ~G4GDMLReadSetup ()
 
- Protected Member Functions inherited from G4GDMLReadSolids
 G4GDMLReadSolids ()
 
virtual ~G4GDMLReadSolids ()
 
void BooleanRead (const xercesc::DOMElement *const, const BooleanOp)
 
void BoxRead (const xercesc::DOMElement *const)
 
void ConeRead (const xercesc::DOMElement *const)
 
void ElconeRead (const xercesc::DOMElement *const)
 
void EllipsoidRead (const xercesc::DOMElement *const)
 
void EltubeRead (const xercesc::DOMElement *const)
 
void XtruRead (const xercesc::DOMElement *const)
 
void HypeRead (const xercesc::DOMElement *const)
 
void MultiUnionNodeRead (const xercesc::DOMElement *const, G4MultiUnion *const)
 
void MultiUnionRead (const xercesc::DOMElement *const)
 
void OrbRead (const xercesc::DOMElement *const)
 
void ParaRead (const xercesc::DOMElement *const)
 
void ParaboloidRead (const xercesc::DOMElement *const)
 
void PolyconeRead (const xercesc::DOMElement *const)
 
void GenericPolyconeRead (const xercesc::DOMElement *const)
 
void PolyhedraRead (const xercesc::DOMElement *const)
 
void GenericPolyhedraRead (const xercesc::DOMElement *const)
 
G4QuadrangularFacetQuadrangularRead (const xercesc::DOMElement *const)
 
void ReflectedSolidRead (const xercesc::DOMElement *const)
 
void ScaledSolidRead (const xercesc::DOMElement *const)
 
G4ExtrudedSolid::ZSection SectionRead (const xercesc::DOMElement *const, G4double)
 
void SphereRead (const xercesc::DOMElement *const)
 
void TessellatedRead (const xercesc::DOMElement *const)
 
void TetRead (const xercesc::DOMElement *const)
 
void TorusRead (const xercesc::DOMElement *const)
 
void GenTrapRead (const xercesc::DOMElement *const)
 
void TrapRead (const xercesc::DOMElement *const)
 
void TrdRead (const xercesc::DOMElement *const)
 
void TubeRead (const xercesc::DOMElement *const)
 
void CutTubeRead (const xercesc::DOMElement *const)
 
void TwistedboxRead (const xercesc::DOMElement *const)
 
void TwistedtrapRead (const xercesc::DOMElement *const)
 
void TwistedtrdRead (const xercesc::DOMElement *const)
 
void TwistedtubsRead (const xercesc::DOMElement *const)
 
G4TriangularFacetTriangularRead (const xercesc::DOMElement *const)
 
G4TwoVector TwoDimVertexRead (const xercesc::DOMElement *const, G4double)
 
zplaneType ZplaneRead (const xercesc::DOMElement *const)
 
rzPointType RZPointRead (const xercesc::DOMElement *const)
 
void OpticalSurfaceRead (const xercesc::DOMElement *const)
 
void PropertyRead (const xercesc::DOMElement *const, G4OpticalSurface *)
 
- Protected Member Functions inherited from G4GDMLReadMaterials
 G4GDMLReadMaterials ()
 
virtual ~G4GDMLReadMaterials ()
 
G4double AtomRead (const xercesc::DOMElement *const)
 
G4int CompositeRead (const xercesc::DOMElement *const, G4String &)
 
G4double DRead (const xercesc::DOMElement *const)
 
G4double PRead (const xercesc::DOMElement *const)
 
G4double TRead (const xercesc::DOMElement *const)
 
G4double MEERead (const xercesc::DOMElement *const)
 
void ElementRead (const xercesc::DOMElement *const)
 
G4double FractionRead (const xercesc::DOMElement *const, G4String &)
 
void IsotopeRead (const xercesc::DOMElement *const)
 
void MaterialRead (const xercesc::DOMElement *const)
 
void MixtureRead (const xercesc::DOMElement *const, G4Element *)
 
void MixtureRead (const xercesc::DOMElement *const, G4Material *)
 
void PropertyRead (const xercesc::DOMElement *const, G4Material *)
 
- Protected Member Functions inherited from G4GDMLReadDefine
 G4GDMLReadDefine ()
 
virtual ~G4GDMLReadDefine ()
 
G4RotationMatrix GetRotationMatrix (const G4ThreeVector &)
 
void VectorRead (const xercesc::DOMElement *const, G4ThreeVector &)
 
G4String RefRead (const xercesc::DOMElement *const)
 
void ConstantRead (const xercesc::DOMElement *const)
 
void MatrixRead (const xercesc::DOMElement *const)
 
void PositionRead (const xercesc::DOMElement *const)
 
void RotationRead (const xercesc::DOMElement *const)
 
void ScaleRead (const xercesc::DOMElement *const)
 
void VariableRead (const xercesc::DOMElement *const)
 
void QuantityRead (const xercesc::DOMElement *const)
 
void ExpressionRead (const xercesc::DOMElement *const)
 
- Protected Member Functions inherited from G4GDMLRead
 G4GDMLRead ()
 
virtual ~G4GDMLRead ()
 
G4String Transcode (const XMLCh *const)
 
G4String GenerateName (const G4String &name, G4bool strip=false)
 
G4String Strip (const G4String &) const
 
void GeneratePhysvolName (const G4String &, G4VPhysicalVolume *)
 
void LoopRead (const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
 
G4GDMLAuxStructType AuxiliaryRead (const xercesc::DOMElement *const auxElem)
 

Protected Attributes

G4GDMLAuxMapType auxMap
 
G4GDMLAssemblyMapType assemblyMap
 
G4LogicalVolumepMotherLogical = nullptr
 
std::map< std::string, G4VPhysicalVolume * > setuptoPV
 
G4bool strip = false
 
- Protected Attributes inherited from G4GDMLReadParamvol
G4GDMLParameterisationparameterisation = nullptr
 
- Protected Attributes inherited from G4GDMLReadSetup
std::map< G4String, G4StringsetupMap
 
- Protected Attributes inherited from G4GDMLReadDefine
G4bool reverseSearch = false
 
std::map< G4String, G4doublequantityMap
 
std::map< G4String, G4ThreeVectorpositionMap
 
std::map< G4String, G4ThreeVectorrotationMap
 
std::map< G4String, G4ThreeVectorscaleMap
 
std::map< G4String, G4GDMLMatrixmatrixMap
 
- Protected Attributes inherited from G4GDMLRead
G4GDMLEvaluator eval
 
G4bool validate = true
 
G4bool check = false
 
G4bool dostrip = true
 

Detailed Description

Definition at line 51 of file G4GDMLReadStructure.hh.

Constructor & Destructor Documentation

◆ G4GDMLReadStructure()

G4GDMLReadStructure::G4GDMLReadStructure ( )

Definition at line 47 of file G4GDMLReadStructure.cc.

◆ ~G4GDMLReadStructure()

G4GDMLReadStructure::~G4GDMLReadStructure ( )
virtual

Definition at line 53 of file G4GDMLReadStructure.cc.

54{
55}

Member Function Documentation

◆ AssemblyRead()

void G4GDMLReadStructure::AssemblyRead ( const xercesc::DOMElement * const  assemblyElement)
protected

Definition at line 818 of file G4GDMLReadStructure.cc.

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}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void PhysvolRead(const xercesc::DOMElement *const, G4AssemblyVolume *assembly=0)
G4GDMLAssemblyMapType assemblyMap
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:70
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const char * name(G4int ptype)
Definition: xmlparse.c:284

Referenced by StructureRead().

◆ AxisRead()

EAxis G4GDMLReadStructure::AxisRead ( const xercesc::DOMElement *const  axisElement)
protected

Definition at line 643 of file G4GDMLReadStructure.cc.

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}
G4double Evaluate(const G4String &)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:156
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

Referenced by ReplicaRead().

◆ BorderSurfaceRead()

void G4GDMLReadStructure::BorderSurfaceRead ( const xercesc::DOMElement * const  bordersurfaceElement)
protected

Definition at line 58 of file G4GDMLReadStructure.cc.

60{
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}
int G4int
Definition: G4Types.hh:85
G4String RefRead(const xercesc::DOMElement *const)
G4SurfaceProperty * GetSurfaceProperty(const G4String &) const
G4VPhysicalVolume * GetPhysvol(const G4String &) const
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:104

Referenced by StructureRead().

◆ Clear()

void G4GDMLReadStructure::Clear ( )

Definition at line 1162 of file G4GDMLReadStructure.cc.

1163{
1164 eval.Clear();
1165 setuptoPV.clear();
1166 auxMap.clear();
1167}
G4GDMLAuxMapType auxMap
std::map< std::string, G4VPhysicalVolume * > setuptoPV

◆ DivisionvolRead()

void G4GDMLReadStructure::DivisionvolRead ( const xercesc::DOMElement * const  divisionvolElement)
protected

Definition at line 143 of file G4GDMLReadStructure.cc.

145{
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}
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
double G4double
Definition: G4Types.hh:83
G4int EvaluateInteger(const G4String &)
G4LogicalVolume * GetVolume(const G4String &) const
G4LogicalVolume * pMotherLogical
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:87
const G4String & GetName() const
static G4PVDivisionFactory * GetInstance()
G4PhysicalVolumesPair Divide(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofDivisions, G4double width, G4double offset)
static G4ReflectionFactory * Instance()
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)

Referenced by Volume_contentRead().

◆ FileRead()

G4LogicalVolume * G4GDMLReadStructure::FileRead ( const xercesc::DOMElement * const  fileElement)
protected

Definition at line 297 of file G4GDMLReadStructure.cc.

299{
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}
std::map< G4LogicalVolume *, G4GDMLAuxListType > G4GDMLAuxMapType
bool G4bool
Definition: G4Types.hh:86
G4String GetSetup(const G4String &)
const G4GDMLAuxMapType * GetAuxMap() const
G4bool validate
Definition: G4GDMLRead.hh:157
void Read(const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
Definition: G4GDMLRead.cc:422

Referenced by PhysvolRead().

◆ GetAssembly()

G4AssemblyVolume * G4GDMLReadStructure::GetAssembly ( const G4String ref) const

Definition at line 1108 of file G4GDMLReadStructure.cc.

1109{
1110 auto pos = assemblyMap.find(ref);
1111 if(pos != assemblyMap.cend())
1112 {
1113 return pos->second;
1114 }
1115 return nullptr;
1116}

Referenced by PhysvolRead().

◆ GetAuxMap()

const G4GDMLAuxMapType * G4GDMLReadStructure::GetAuxMap ( ) const
inline

Definition at line 63 of file G4GDMLReadStructure.hh.

63{ return &auxMap; }

Referenced by FileRead().

◆ GetPhysvol()

G4VPhysicalVolume * G4GDMLReadStructure::GetPhysvol ( const G4String ref) const

Definition at line 1076 of file G4GDMLReadStructure.cc.

1077{
1078 G4VPhysicalVolume* physvolPtr
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}
static G4PhysicalVolumeStore * GetInstance()
G4VPhysicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const

Referenced by BorderSurfaceRead().

◆ GetVolume()

G4LogicalVolume * G4GDMLReadStructure::GetVolume ( const G4String ref) const
virtual

Implements G4GDMLRead.

Definition at line 1092 of file G4GDMLReadStructure.cc.

1093{
1094 G4LogicalVolume* volumePtr
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}
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const
static G4LogicalVolumeStore * GetInstance()

Referenced by DivisionvolRead(), FileRead(), GetWorldVolume(), PhysvolRead(), ReplicavolRead(), and SkinSurfaceRead().

◆ GetVolumeAuxiliaryInformation()

G4GDMLAuxListType G4GDMLReadStructure::GetVolumeAuxiliaryInformation ( G4LogicalVolume logvol) const

Definition at line 1119 of file G4GDMLReadStructure.cc.

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}
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType

◆ GetWorldVolume()

G4VPhysicalVolume * G4GDMLReadStructure::GetWorldVolume ( const G4String setupName)

Definition at line 1134 of file G4GDMLReadStructure.cc.

1136{
1137 G4String sname = GetSetup(setupName);
1138 if(sname == "")
1139 {
1140 return nullptr;
1141 }
1142
1143 G4LogicalVolume* volume = GetVolume(GenerateName(sname, dostrip));
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}
CLHEP::Hep3Vector G4ThreeVector
G4bool dostrip
Definition: G4GDMLRead.hh:159
void SetVisAttributes(const G4VisAttributes *pVA)
static const G4VisAttributes & GetInvisible()

◆ PhysvolRead()

void G4GDMLReadStructure::PhysvolRead ( const xercesc::DOMElement * const  physvolElement,
G4AssemblyVolume assembly = 0 
)
protected

Definition at line 366 of file G4GDMLReadStructure.cc.

368{
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}
HepGeom::Scale3D G4Scale3D
void MakeImprint(G4LogicalVolume *pMotherLV, G4ThreeVector &translationInMother, G4RotationMatrix *pRotationInMother, G4int copyNumBase=0, G4bool surfCheck=false)
G4ThreeVector GetScale(const G4String &)
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
G4RotationMatrix GetRotationMatrix(const G4ThreeVector &)
G4ThreeVector GetPosition(const G4String &)
G4ThreeVector GetRotation(const G4String &)
G4AssemblyVolume * GetAssembly(const G4String &) const
G4LogicalVolume * FileRead(const xercesc::DOMElement *const)
G4bool check
Definition: G4GDMLRead.hh:158
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)

Referenced by AssemblyRead(), and Volume_contentRead().

◆ QuantityRead()

G4double G4GDMLReadStructure::QuantityRead ( const xercesc::DOMElement *const  readElement)
protected

Definition at line 713 of file G4GDMLReadStructure.cc.

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}

Referenced by ReplicaRead().

◆ ReplicaRead()

void G4GDMLReadStructure::ReplicaRead ( const xercesc::DOMElement *const  replicaElement,
G4LogicalVolume logvol,
G4int  number 
)
protected

Definition at line 563 of file G4GDMLReadStructure.cc.

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;
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}
EAxis AxisRead(const xercesc::DOMElement *const axisElement)
G4double QuantityRead(const xercesc::DOMElement *const readElement)
G4PhysicalVolumesPair Replicate(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofReplicas, G4double width, G4double offset=0.)

Referenced by ReplicavolRead().

◆ ReplicavolRead()

void G4GDMLReadStructure::ReplicavolRead ( const xercesc::DOMElement * const  replicavolElement,
G4int  number 
)
protected

Definition at line 520 of file G4GDMLReadStructure.cc.

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}
void ReplicaRead(const xercesc::DOMElement *const replicaElement, G4LogicalVolume *logvol, G4int number)

Referenced by Volume_contentRead().

◆ SkinSurfaceRead()

void G4GDMLReadStructure::SkinSurfaceRead ( const xercesc::DOMElement * const  skinsurfaceElement)
protected

Definition at line 866 of file G4GDMLReadStructure.cc.

868{
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}

Referenced by StructureRead().

◆ StructureRead()

void G4GDMLReadStructure::StructureRead ( const xercesc::DOMElement * const  structureElement)
virtual

Implements G4GDMLRead.

Definition at line 1022 of file G4GDMLReadStructure.cc.

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 {
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}
void BorderSurfaceRead(const xercesc::DOMElement *const)
void AssemblyRead(const xercesc::DOMElement *const)
void SkinSurfaceRead(const xercesc::DOMElement *const)
virtual void VolumeRead(const xercesc::DOMElement *const)
void LoopRead(const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
Definition: G4GDMLRead.cc:194
virtual void StructureRead(const xercesc::DOMElement *const)=0

◆ Volume_contentRead()

void G4GDMLReadStructure::Volume_contentRead ( const xercesc::DOMElement * const  volumeElement)
virtual

Implements G4GDMLRead.

Definition at line 942 of file G4GDMLReadStructure.cc.

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 {
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 {
1012 }
1013 else
1014 {
1015 G4cout << "Treating unknown GDML tag in volume '" << tag
1016 << "' as GDML extension..." << G4endl;
1017 }
1018 }
1019}
virtual void ParamvolRead(const xercesc::DOMElement *const, G4LogicalVolume *)
void DivisionvolRead(const xercesc::DOMElement *const)
void ReplicavolRead(const xercesc::DOMElement *const, G4int number)
virtual void Volume_contentRead(const xercesc::DOMElement *const)=0

Referenced by VolumeRead().

◆ VolumeRead()

void G4GDMLReadStructure::VolumeRead ( const xercesc::DOMElement * const  volumeElement)
virtual

Definition at line 763 of file G4GDMLReadStructure.cc.

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}
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
G4VSolid * GetSolid(const G4String &) const
virtual void Volume_contentRead(const xercesc::DOMElement *const)
G4GDMLAuxStructType AuxiliaryRead(const xercesc::DOMElement *const auxElem)
Definition: G4GDMLRead.cc:295

Referenced by StructureRead().

Member Data Documentation

◆ assemblyMap

G4GDMLAssemblyMapType G4GDMLReadStructure::assemblyMap
protected

Definition at line 88 of file G4GDMLReadStructure.hh.

Referenced by AssemblyRead(), and GetAssembly().

◆ auxMap

G4GDMLAuxMapType G4GDMLReadStructure::auxMap
protected

◆ pMotherLogical

G4LogicalVolume* G4GDMLReadStructure::pMotherLogical = nullptr
protected

◆ setuptoPV

std::map<std::string, G4VPhysicalVolume*> G4GDMLReadStructure::setuptoPV
protected

Definition at line 90 of file G4GDMLReadStructure.hh.

Referenced by Clear(), and GetWorldVolume().

◆ strip

G4bool G4GDMLReadStructure::strip = false
protected

Definition at line 91 of file G4GDMLReadStructure.hh.


The documentation for this class was generated from the following files: