Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes
G4GDMLReadDefine Class Referenceabstract

#include <G4GDMLReadDefine.hh>

Inheritance diagram for G4GDMLReadDefine:
G4GDMLRead G4GDMLReadMaterials G4GDMLReadSolids G4GDMLReadSetup G4GDMLReadParamvol G4GDMLReadStructure

Public Member Functions

virtual void DefineRead (const xercesc::DOMElement *const)
 
virtual void ExtensionRead (const xercesc::DOMElement *const)
 
const G4GDMLAuxListTypeGetAuxList () const
 
G4double GetConstant (const G4String &)
 
G4GDMLMatrix GetMatrix (const G4String &)
 
G4ThreeVector GetPosition (const G4String &)
 
G4double GetQuantity (const G4String &)
 
G4ThreeVector GetRotation (const G4String &)
 
G4ThreeVector GetScale (const G4String &)
 
virtual G4String GetSetup (const G4String &)=0
 
G4double GetVariable (const G4String &)
 
virtual G4LogicalVolumeGetVolume (const G4String &) const =0
 
G4bool IsValidID (const G4String &) const
 
virtual void MaterialsRead (const xercesc::DOMElement *const)=0
 
void OverlapCheck (G4bool)
 
virtual void Paramvol_contentRead (const xercesc::DOMElement *const)=0
 
void Read (const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
 
void SetReverseSearch (G4bool flag)
 
virtual void SetupRead (const xercesc::DOMElement *const)=0
 
virtual void SolidsRead (const xercesc::DOMElement *const)=0
 
void StripName (G4String &) const
 
void StripNames () const
 
virtual void StructureRead (const xercesc::DOMElement *const)=0
 
virtual void UserinfoRead (const xercesc::DOMElement *const)
 
virtual void Volume_contentRead (const xercesc::DOMElement *const)=0
 

Protected Member Functions

G4GDMLAuxStructType AuxiliaryRead (const xercesc::DOMElement *const auxElem)
 
void ConstantRead (const xercesc::DOMElement *const)
 
void ExpressionRead (const xercesc::DOMElement *const)
 
 G4GDMLReadDefine ()
 
G4String GenerateName (const G4String &name, G4bool strip=false)
 
void GeneratePhysvolName (const G4String &, G4VPhysicalVolume *)
 
G4RotationMatrix GetRotationMatrix (const G4ThreeVector &)
 
void LoopRead (const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
 
void MatrixRead (const xercesc::DOMElement *const)
 
void PositionRead (const xercesc::DOMElement *const)
 
void QuantityRead (const xercesc::DOMElement *const)
 
G4String RefRead (const xercesc::DOMElement *const)
 
void RotationRead (const xercesc::DOMElement *const)
 
void ScaleRead (const xercesc::DOMElement *const)
 
G4String Strip (const G4String &) const
 
G4String Transcode (const XMLCh *const)
 
void VariableRead (const xercesc::DOMElement *const)
 
void VectorRead (const xercesc::DOMElement *const, G4ThreeVector &)
 
virtual ~G4GDMLReadDefine ()
 

Protected Attributes

G4bool check = false
 
G4bool dostrip = true
 
G4GDMLEvaluator eval
 
std::map< G4String, G4GDMLMatrixmatrixMap
 
std::map< G4String, G4ThreeVectorpositionMap
 
std::map< G4String, G4doublequantityMap
 
G4bool reverseSearch = false
 
std::map< G4String, G4ThreeVectorrotationMap
 
std::map< G4String, G4ThreeVectorscaleMap
 
G4bool validate = true
 

Private Attributes

G4GDMLAuxListType auxGlobalList
 
G4int inLoop = 0
 
G4int loopCount = 0
 

Detailed Description

Definition at line 66 of file G4GDMLReadDefine.hh.

Constructor & Destructor Documentation

◆ G4GDMLReadDefine()

G4GDMLReadDefine::G4GDMLReadDefine ( )
protected

Definition at line 142 of file G4GDMLReadDefine.cc.

143 : G4GDMLRead()
144{
145}

◆ ~G4GDMLReadDefine()

G4GDMLReadDefine::~G4GDMLReadDefine ( )
protectedvirtual

Definition at line 148 of file G4GDMLReadDefine.cc.

149{
150}

Member Function Documentation

◆ AuxiliaryRead()

G4GDMLAuxStructType G4GDMLRead::AuxiliaryRead ( const xercesc::DOMElement *const  auxElem)
protectedinherited

Definition at line 295 of file G4GDMLRead.cc.

297{
298 G4GDMLAuxStructType auxstruct = { "", "", "", 0 };
299 G4GDMLAuxListType* auxList = nullptr;
300
301 const xercesc::DOMNamedNodeMap* const attributes =
302 auxiliaryElement->getAttributes();
303 XMLSize_t attributeCount = attributes->getLength();
304
305 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
306 ++attribute_index)
307 {
308 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
309
310 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
311 {
312 continue;
313 }
314
315 const xercesc::DOMAttr* const attribute =
316 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
317 if(!attribute)
318 {
319 G4Exception("G4GDMLRead::AuxiliaryRead()", "InvalidRead", FatalException,
320 "No attribute found!");
321 return auxstruct;
322 }
323 const G4String attName = Transcode(attribute->getName());
324 const G4String attValue = Transcode(attribute->getValue());
325
326 if(attName == "auxtype")
327 {
328 auxstruct.type = attValue;
329 }
330 else if(attName == "auxvalue")
331 {
332 auxstruct.value = attValue;
333 }
334 else if(attName == "auxunit")
335 {
336 auxstruct.unit = attValue;
337 }
338 }
339
340 for(xercesc::DOMNode* iter = auxiliaryElement->getFirstChild();
341 iter != nullptr; iter = iter->getNextSibling())
342 {
343 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
344 {
345 continue;
346 }
347
348 const xercesc::DOMElement* const child =
349 dynamic_cast<xercesc::DOMElement*>(iter);
350 if(!child)
351 {
352 G4Exception("G4GDMLRead::AuxiliaryRead()", "InvalidRead", FatalException,
353 "No child found!");
354 break;
355 }
356 const G4String tag = Transcode(child->getTagName());
357
358 if(tag == "auxiliary")
359 {
360 if(!auxList)
361 {
362 auxList = new G4GDMLAuxListType;
363 }
364 auxList->push_back(AuxiliaryRead(child));
365 }
366 }
367
368 if(auxList)
369 {
370 auxstruct.auxList = auxList;
371 }
372
373 return auxstruct;
374}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
G4GDMLAuxStructType AuxiliaryRead(const xercesc::DOMElement *const auxElem)
Definition: G4GDMLRead.cc:295
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
std::vector< G4GDMLAuxStructType > * auxList
Definition: xmlparse.cc:187

References G4GDMLRead::AuxiliaryRead(), G4GDMLAuxStructType::auxList, FatalException, G4Exception(), G4GDMLRead::Transcode(), G4GDMLAuxStructType::type, G4GDMLAuxStructType::unit, and G4GDMLAuxStructType::value.

Referenced by G4GDMLRead::AuxiliaryRead(), G4GDMLRead::UserinfoRead(), and G4GDMLReadStructure::VolumeRead().

◆ ConstantRead()

void G4GDMLReadDefine::ConstantRead ( const xercesc::DOMElement * const  constantElement)
protected

Definition at line 167 of file G4GDMLReadDefine.cc.

169{
170 G4String name = "";
171 G4double value = 0.0;
172
173 const xercesc::DOMNamedNodeMap* const attributes =
174 constantElement->getAttributes();
175 XMLSize_t attributeCount = attributes->getLength();
176
177 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
178 ++attribute_index)
179 {
180 xercesc::DOMNode* node = attributes->item(attribute_index);
181
182 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
183 {
184 continue;
185 }
186
187 const xercesc::DOMAttr* const attribute =
188 dynamic_cast<xercesc::DOMAttr*>(node);
189 if(attribute == nullptr)
190 {
191 G4Exception("G4GDMLRead::ConstantRead()", "InvalidRead", FatalException,
192 "No attribute found!");
193 return;
194 }
195 const G4String attName = Transcode(attribute->getName());
196 const G4String attValue = Transcode(attribute->getValue());
197
198 if(attName == "name")
199 {
200 name = attValue;
201 }
202 else if(attName == "value")
203 {
204 value = eval.Evaluate(attValue);
205 }
206 }
207
208 eval.DefineConstant(name, value);
209}
double G4double
Definition: G4Types.hh:83
G4double Evaluate(const G4String &)
void DefineConstant(const G4String &, G4double)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:156
const char * name(G4int ptype)

References G4GDMLEvaluator::DefineConstant(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4InuclParticleNames::name(), and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ DefineRead()

void G4GDMLReadDefine::DefineRead ( const xercesc::DOMElement * const  defineElement)
virtual

Implements G4GDMLRead.

Definition at line 596 of file G4GDMLReadDefine.cc.

598{
599#ifdef G4VERBOSE
600 G4cout << "G4GDML: Reading definitions..." << G4endl;
601#endif
602 for(xercesc::DOMNode* iter = defineElement->getFirstChild(); iter != nullptr;
603 iter = iter->getNextSibling())
604 {
605 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
606 {
607 continue;
608 }
609
610 const xercesc::DOMElement* const child =
611 dynamic_cast<xercesc::DOMElement*>(iter);
612 if(child == nullptr)
613 {
614 G4Exception("G4GDMLRead::DefineRead()", "InvalidRead", FatalException,
615 "No child found!");
616 return;
617 }
618 const G4String tag = Transcode(child->getTagName());
619
620 if(tag == "constant")
621 {
622 ConstantRead(child);
623 }
624 else if(tag == "matrix")
625 {
626 MatrixRead(child);
627 }
628 else if(tag == "position")
629 {
630 PositionRead(child);
631 }
632 else if(tag == "rotation")
633 {
634 RotationRead(child);
635 }
636 else if(tag == "scale")
637 {
638 ScaleRead(child);
639 }
640 else if(tag == "variable")
641 {
642 VariableRead(child);
643 }
644 else if(tag == "quantity")
645 {
646 QuantityRead(child);
647 }
648 else if(tag == "expression")
649 {
650 ExpressionRead(child);
651 }
652 else
653 {
654 G4String error_msg = "Unknown tag in define: " + tag;
655 G4Exception("G4GDMLReadDefine::defineRead()", "ReadError", FatalException,
656 error_msg);
657 }
658 }
659}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void ExpressionRead(const xercesc::DOMElement *const)
void VariableRead(const xercesc::DOMElement *const)
void QuantityRead(const xercesc::DOMElement *const)
void PositionRead(const xercesc::DOMElement *const)
void MatrixRead(const xercesc::DOMElement *const)
void RotationRead(const xercesc::DOMElement *const)
void ScaleRead(const xercesc::DOMElement *const)
void ConstantRead(const xercesc::DOMElement *const)

References ConstantRead(), ExpressionRead(), FatalException, G4cout, G4endl, G4Exception(), MatrixRead(), PositionRead(), QuantityRead(), RotationRead(), ScaleRead(), G4GDMLRead::Transcode(), and VariableRead().

Referenced by G4GDMLReadMaterials::MaterialsRead(), and G4GDMLReadSolids::SolidsRead().

◆ ExpressionRead()

void G4GDMLReadDefine::ExpressionRead ( const xercesc::DOMElement * const  expElement)
protected

Definition at line 212 of file G4GDMLReadDefine.cc.

214{
215 G4String name = "";
216 G4double value = 0.0;
217
218 const xercesc::DOMNamedNodeMap* const attributes =
219 expElement->getAttributes();
220 XMLSize_t attributeCount = attributes->getLength();
221
222 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
223 ++attribute_index)
224 {
225 xercesc::DOMNode* node = attributes->item(attribute_index);
226
227 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
228 {
229 continue;
230 }
231
232 const xercesc::DOMAttr* const attribute =
233 dynamic_cast<xercesc::DOMAttr*>(node);
234 if(attribute == nullptr)
235 {
236 G4Exception("G4GDMLRead::ExpressionRead()", "InvalidRead", FatalException,
237 "No attribute found!");
238 return;
239 }
240 const G4String attName = Transcode(attribute->getName());
241 const G4String attValue = Transcode(attribute->getValue());
242
243 if(attName == "name")
244 {
245 name = attValue;
246 }
247 }
248
249 const G4String expValue = Transcode(expElement->getTextContent());
250 value = eval.Evaluate(expValue);
251 eval.DefineConstant(name, value);
252}

References G4GDMLEvaluator::DefineConstant(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4InuclParticleNames::name(), and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ ExtensionRead()

void G4GDMLRead::ExtensionRead ( const xercesc::DOMElement * const  )
virtualinherited

Definition at line 414 of file G4GDMLRead.cc.

415{
416 G4String error_msg = "No handle to user-code for parsing extensions!";
417 G4Exception("G4GDMLRead::ExtensionRead()", "NotImplemented", JustWarning,
418 error_msg);
419}
@ JustWarning

References G4Exception(), and JustWarning.

Referenced by G4GDMLRead::Read().

◆ GenerateName()

G4String G4GDMLRead::GenerateName ( const G4String name,
G4bool  strip = false 
)
protectedinherited

Definition at line 70 of file G4GDMLRead.cc.

71{
72 G4String nameOut(nameIn);
73
74 if(inLoop > 0)
75 {
76 nameOut = eval.SolveBrackets(nameOut);
77 }
78 if(strip)
79 {
80 StripName(nameOut);
81 }
82
83 return nameOut;
84}
G4String SolveBrackets(const G4String &)
void StripName(G4String &) const
Definition: G4GDMLRead.cc:112
G4int inLoop
Definition: G4GDMLRead.hh:163
void strip(G4String &str, char c=' ')
Remove leading and trailing characters from string.

References G4GDMLRead::eval, G4GDMLRead::inLoop, G4GDMLEvaluator::SolveBrackets(), G4StrUtil::strip(), and G4GDMLRead::StripName().

Referenced by G4GDMLReadStructure::AssemblyRead(), G4GDMLReadSolids::BooleanRead(), G4GDMLReadStructure::BorderSurfaceRead(), G4GDMLReadSolids::BoxRead(), G4GDMLReadSolids::ConeRead(), G4GDMLReadSolids::CutTubeRead(), G4GDMLReadStructure::DivisionvolRead(), G4GDMLReadSolids::ElconeRead(), G4GDMLReadMaterials::ElementRead(), G4GDMLReadSolids::EllipsoidRead(), G4GDMLReadSolids::EltubeRead(), G4GDMLReadStructure::FileRead(), G4GDMLReadSolids::GenericPolyconeRead(), G4GDMLReadSolids::GenericPolyhedraRead(), G4GDMLReadSolids::GenTrapRead(), G4GDMLReadStructure::GetWorldVolume(), G4GDMLReadSolids::HypeRead(), G4GDMLReadMaterials::IsotopeRead(), G4GDMLReadMaterials::MaterialRead(), MatrixRead(), G4GDMLReadMaterials::MixtureRead(), G4GDMLReadSolids::MultiUnionNodeRead(), G4GDMLReadSolids::MultiUnionRead(), G4GDMLReadSolids::OpticalSurfaceRead(), G4GDMLReadSolids::OrbRead(), G4GDMLReadSolids::ParaboloidRead(), G4GDMLReadParamvol::ParametersRead(), G4GDMLReadParamvol::ParamvolRead(), G4GDMLReadSolids::ParaRead(), G4GDMLReadStructure::PhysvolRead(), G4GDMLReadSolids::PolyconeRead(), G4GDMLReadSolids::PolyhedraRead(), PositionRead(), G4GDMLReadMaterials::PropertyRead(), G4GDMLReadSolids::PropertyRead(), G4GDMLReadSolids::QuadrangularRead(), G4GDMLReadSolids::ReflectedSolidRead(), G4GDMLReadStructure::ReplicaRead(), G4GDMLReadStructure::ReplicavolRead(), RotationRead(), G4GDMLReadSolids::ScaledSolidRead(), ScaleRead(), G4GDMLReadSetup::SetupRead(), G4GDMLReadStructure::SkinSurfaceRead(), G4GDMLReadSolids::SphereRead(), G4GDMLReadSolids::TessellatedRead(), G4GDMLReadSolids::TetRead(), G4GDMLReadSolids::TorusRead(), G4GDMLReadSolids::TrapRead(), G4GDMLReadSolids::TrdRead(), G4GDMLReadSolids::TriangularRead(), G4GDMLReadSolids::TubeRead(), G4GDMLReadSolids::TwistedboxRead(), G4GDMLReadSolids::TwistedtrapRead(), G4GDMLReadSolids::TwistedtrdRead(), G4GDMLReadSolids::TwistedtubsRead(), G4GDMLReadStructure::VolumeRead(), and G4GDMLReadSolids::XtruRead().

◆ GeneratePhysvolName()

void G4GDMLRead::GeneratePhysvolName ( const G4String nameIn,
G4VPhysicalVolume physvol 
)
protectedinherited

Definition at line 87 of file G4GDMLRead.cc.

89{
90 G4String nameOut(nameIn);
91
92 if(nameIn.empty())
93 {
94 std::stringstream stream;
95 stream << physvol->GetLogicalVolume()->GetName() << "_PV";
96 nameOut = stream.str();
97 }
98 nameOut = eval.SolveBrackets(nameOut);
99
100 physvol->SetName(nameOut);
101}
const G4String & GetName() const
G4LogicalVolume * GetLogicalVolume() const
void SetName(const G4String &pName)

References G4GDMLRead::eval, G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetName(), G4VPhysicalVolume::SetName(), and G4GDMLEvaluator::SolveBrackets().

Referenced by G4GDMLReadStructure::DivisionvolRead(), G4GDMLReadStructure::PhysvolRead(), and G4GDMLReadStructure::ReplicaRead().

◆ GetAuxList()

const G4GDMLAuxListType * G4GDMLRead::GetAuxList ( ) const
inherited

Definition at line 563 of file G4GDMLRead.cc.

564{
565 return &auxGlobalList;
566}
G4GDMLAuxListType auxGlobalList
Definition: G4GDMLRead.hh:164

References G4GDMLRead::auxGlobalList.

◆ GetConstant()

G4double G4GDMLReadDefine::GetConstant ( const G4String ref)

Definition at line 758 of file G4GDMLReadDefine.cc.

759{
760 return eval.GetConstant(ref);
761}
G4double GetConstant(const G4String &)

References G4GDMLRead::eval, and G4GDMLEvaluator::GetConstant().

◆ GetMatrix()

G4GDMLMatrix G4GDMLReadDefine::GetMatrix ( const G4String ref)

Definition at line 818 of file G4GDMLReadDefine.cc.

819{
820 if(matrixMap.find(ref) == matrixMap.end())
821 {
822 G4String error_msg = "Matrix '" + ref + "' was not found!";
823 G4Exception("G4GDMLReadDefine::getMatrix()", "ReadError", FatalException,
824 error_msg);
825 }
826 return matrixMap[ref];
827}
std::map< G4String, G4GDMLMatrix > matrixMap

References FatalException, G4Exception(), and matrixMap.

Referenced by G4GDMLReadMaterials::PropertyRead(), and G4GDMLReadSolids::PropertyRead().

◆ GetPosition()

G4ThreeVector G4GDMLReadDefine::GetPosition ( const G4String ref)

Definition at line 782 of file G4GDMLReadDefine.cc.

783{
784 if(positionMap.find(ref) == positionMap.cend())
785 {
786 G4String error_msg = "Position '" + ref + "' was not found!";
787 G4Exception("G4GDMLReadDefine::getPosition()", "ReadError", FatalException,
788 error_msg);
789 }
790 return positionMap[ref];
791}
std::map< G4String, G4ThreeVector > positionMap

References FatalException, G4Exception(), and positionMap.

Referenced by G4GDMLReadSolids::BooleanRead(), G4GDMLReadSolids::MultiUnionNodeRead(), G4GDMLReadParamvol::ParametersRead(), G4GDMLReadStructure::PhysvolRead(), G4GDMLReadSolids::QuadrangularRead(), G4GDMLReadStructure::ReplicaRead(), G4GDMLReadSolids::TetRead(), and G4GDMLReadSolids::TriangularRead().

◆ GetQuantity()

G4double G4GDMLReadDefine::GetQuantity ( const G4String ref)

Definition at line 770 of file G4GDMLReadDefine.cc.

771{
772 if(quantityMap.find(ref) == quantityMap.cend())
773 {
774 G4String error_msg = "Quantity '" + ref + "' was not found!";
775 G4Exception("G4GDMLReadDefine::getQuantity()", "ReadError", FatalException,
776 error_msg);
777 }
778 return quantityMap[ref];
779}
std::map< G4String, G4double > quantityMap

References FatalException, G4Exception(), and quantityMap.

Referenced by G4GDMLReadMaterials::MaterialRead().

◆ GetRotation()

G4ThreeVector G4GDMLReadDefine::GetRotation ( const G4String ref)

Definition at line 794 of file G4GDMLReadDefine.cc.

795{
796 if(rotationMap.find(ref) == rotationMap.cend())
797 {
798 G4String error_msg = "Rotation '" + ref + "' was not found!";
799 G4Exception("G4GDMLReadDefine::getRotation()", "ReadError", FatalException,
800 error_msg);
801 }
802 return rotationMap[ref];
803}
std::map< G4String, G4ThreeVector > rotationMap

References FatalException, G4Exception(), and rotationMap.

Referenced by G4GDMLReadSolids::BooleanRead(), G4GDMLReadSolids::MultiUnionNodeRead(), G4GDMLReadParamvol::ParametersRead(), G4GDMLReadStructure::PhysvolRead(), and G4GDMLReadStructure::ReplicaRead().

◆ GetRotationMatrix()

G4RotationMatrix G4GDMLReadDefine::GetRotationMatrix ( const G4ThreeVector angles)
protected

Definition at line 153 of file G4GDMLReadDefine.cc.

155{
157
158 rot.rotateX(angles.x());
159 rot.rotateY(angles.y());
160 rot.rotateZ(angles.z());
161 rot.rectify(); // Rectify matrix from possible roundoff errors
162
163 return rot;
164}
double z() const
double x() const
double y() const
HepRotation & rotateX(double delta)
Definition: Rotation.cc:61
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74

References CLHEP::HepRotation::rectify(), CLHEP::HepRotation::rotateX(), CLHEP::HepRotation::rotateY(), CLHEP::HepRotation::rotateZ(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by G4GDMLReadSolids::BooleanRead(), G4GDMLReadSolids::MultiUnionNodeRead(), G4GDMLReadStructure::PhysvolRead(), and G4GDMLReadSolids::ReflectedSolidRead().

◆ GetScale()

G4ThreeVector G4GDMLReadDefine::GetScale ( const G4String ref)

Definition at line 806 of file G4GDMLReadDefine.cc.

807{
808 if(scaleMap.find(ref) == scaleMap.end())
809 {
810 G4String error_msg = "Scale '" + ref + "' was not found!";
811 G4Exception("G4GDMLReadDefine::getScale()", "ReadError", FatalException,
812 error_msg);
813 }
814 return scaleMap[ref];
815}
std::map< G4String, G4ThreeVector > scaleMap

References FatalException, G4Exception(), and scaleMap.

Referenced by G4GDMLReadStructure::PhysvolRead(), and G4GDMLReadSolids::ScaledSolidRead().

◆ GetSetup()

virtual G4String G4GDMLRead::GetSetup ( const G4String )
pure virtualinherited

Implemented in G4GDMLReadSetup.

◆ GetVariable()

G4double G4GDMLReadDefine::GetVariable ( const G4String ref)

Definition at line 764 of file G4GDMLReadDefine.cc.

765{
766 return eval.GetVariable(ref);
767}
G4double GetVariable(const G4String &)

References G4GDMLRead::eval, and G4GDMLEvaluator::GetVariable().

◆ GetVolume()

virtual G4LogicalVolume * G4GDMLRead::GetVolume ( const G4String ) const
pure virtualinherited

◆ IsValidID()

G4bool G4GDMLReadDefine::IsValidID ( const G4String ref) const

Definition at line 752 of file G4GDMLReadDefine.cc.

753{
754 return eval.IsVariable(ref);
755}
G4bool IsVariable(const G4String &) const

References G4GDMLRead::eval, and G4GDMLEvaluator::IsVariable().

◆ LoopRead()

void G4GDMLRead::LoopRead ( const xercesc::DOMElement * const  element,
void(G4GDMLRead::*)(const xercesc::DOMElement *const)  func 
)
protectedinherited

Definition at line 194 of file G4GDMLRead.cc.

196{
197 G4String var;
198 G4String from;
199 G4String to;
200 G4String step;
201
202 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
203 XMLSize_t attributeCount = attributes->getLength();
204
205 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
206 ++attribute_index)
207 {
208 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
209
210 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
211 {
212 continue;
213 }
214
215 const xercesc::DOMAttr* const attribute =
216 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
217 if(!attribute)
218 {
219 G4Exception("G4GDMLRead::LoopRead()", "InvalidRead", FatalException,
220 "No attribute found!");
221 return;
222 }
223 const G4String attribute_name = Transcode(attribute->getName());
224 const G4String attribute_value = Transcode(attribute->getValue());
225
226 if(attribute_name == "for")
227 {
228 var = attribute_value;
229 }
230 else if(attribute_name == "from")
231 {
232 from = attribute_value;
233 }
234 else if(attribute_name == "to")
235 {
236 to = attribute_value;
237 }
238 else if(attribute_name == "step")
239 {
240 step = attribute_value;
241 }
242 }
243
244 if(var.empty())
245 {
246 G4Exception("G4GDMLRead::loopRead()", "InvalidRead", FatalException,
247 "No variable is determined for loop!");
248 }
249
250 if(!eval.IsVariable(var))
251 {
252 G4Exception("G4GDMLRead::loopRead()", "InvalidRead", FatalException,
253 "Variable is not defined in loop!");
254 }
255
256 G4int _var = eval.EvaluateInteger(var);
257 G4int _from = eval.EvaluateInteger(from);
258 G4int _to = eval.EvaluateInteger(to);
259 G4int _step = eval.EvaluateInteger(step);
260
261 if(!from.empty())
262 {
263 _var = _from;
264 }
265
266 if((_from < _to) && (_step <= 0))
267 {
268 G4Exception("G4GDMLRead::loopRead()", "InvalidRead", FatalException,
269 "Infinite loop!");
270 }
271 if((_from > _to) && (_step >= 0))
272 {
273 G4Exception("G4GDMLRead::loopRead()", "InvalidRead", FatalException,
274 "Infinite loop!");
275 }
276
277 ++inLoop;
278
279 while(_var <= _to)
280 {
281 eval.SetVariable(var, _var);
282 (this->*func)(element);
283 _var += _step;
284 ++loopCount;
285 }
286
287 --inLoop;
288 if(!inLoop)
289 {
290 loopCount = 0;
291 }
292}
int G4int
Definition: G4Types.hh:85
void SetVariable(const G4String &, G4double)
G4int EvaluateInteger(const G4String &)
G4int loopCount
Definition: G4GDMLRead.hh:163

References G4GDMLRead::eval, G4GDMLEvaluator::EvaluateInteger(), FatalException, G4Exception(), G4GDMLRead::inLoop, G4GDMLEvaluator::IsVariable(), G4GDMLRead::loopCount, G4GDMLEvaluator::SetVariable(), and G4GDMLRead::Transcode().

Referenced by G4GDMLReadParamvol::ParameterisedRead(), G4GDMLReadParamvol::Paramvol_contentRead(), G4GDMLReadSolids::SolidsRead(), G4GDMLReadStructure::StructureRead(), and G4GDMLReadStructure::Volume_contentRead().

◆ MaterialsRead()

virtual void G4GDMLRead::MaterialsRead ( const xercesc::DOMElement * const  )
pure virtualinherited

Implemented in G4GDMLReadMaterials.

Referenced by G4GDMLRead::Read().

◆ MatrixRead()

void G4GDMLReadDefine::MatrixRead ( const xercesc::DOMElement * const  matrixElement)
protected

Definition at line 255 of file G4GDMLReadDefine.cc.

257{
258 G4String name = "";
259 G4int coldim = 0;
260 G4String values = "";
261
262 const xercesc::DOMNamedNodeMap* const attributes =
263 matrixElement->getAttributes();
264 XMLSize_t attributeCount = attributes->getLength();
265
266 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
267 ++attribute_index)
268 {
269 xercesc::DOMNode* node = attributes->item(attribute_index);
270
271 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
272 {
273 continue;
274 }
275
276 const xercesc::DOMAttr* const attribute =
277 dynamic_cast<xercesc::DOMAttr*>(node);
278 if(attribute == nullptr)
279 {
280 G4Exception("G4GDMLRead::MatrixRead()", "InvalidRead", FatalException,
281 "No attribute found!");
282 return;
283 }
284 const G4String attName = Transcode(attribute->getName());
285 const G4String attValue = Transcode(attribute->getValue());
286
287 if(attName == "name")
288 {
289 name = GenerateName(attValue);
290 }
291 else if(attName == "coldim")
292 {
293 coldim = eval.EvaluateInteger(attValue);
294 }
295 else if(attName == "values")
296 {
297 values = attValue;
298 }
299 }
300
301 std::stringstream MatrixValueStream(values);
302 std::vector<G4double> valueList;
303
304 while(!MatrixValueStream.eof())
305 {
306 G4String MatrixValue;
307 MatrixValueStream >> MatrixValue;
308 valueList.push_back(eval.Evaluate(MatrixValue));
309 }
310
311 eval.DefineMatrix(name, coldim, valueList);
312
313 G4GDMLMatrix matrix(valueList.size() / coldim, coldim);
314
315 for(std::size_t i = 0; i < valueList.size(); ++i)
316 {
317 matrix.Set(i / coldim, i % coldim, valueList[i]);
318 }
319
320 matrixMap[name] = matrix;
321}
void DefineMatrix(const G4String &, G4int, std::vector< G4double >)
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:70

References G4GDMLEvaluator::DefineMatrix(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), G4GDMLEvaluator::EvaluateInteger(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), matrixMap, G4InuclParticleNames::name(), G4GDMLMatrix::Set(), and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ OverlapCheck()

void G4GDMLRead::OverlapCheck ( G4bool  flag)
inherited

Definition at line 64 of file G4GDMLRead.cc.

65{
66 check = flag;
67}
G4bool check
Definition: G4GDMLRead.hh:158

References G4GDMLRead::check.

◆ Paramvol_contentRead()

virtual void G4GDMLRead::Paramvol_contentRead ( const xercesc::DOMElement * const  )
pure virtualinherited

◆ PositionRead()

void G4GDMLReadDefine::PositionRead ( const xercesc::DOMElement * const  positionElement)
protected

Definition at line 324 of file G4GDMLReadDefine.cc.

326{
327 G4String name = "";
328 G4double unit = 1.0;
329 G4ThreeVector position(0., 0., 0.);
330
331 const xercesc::DOMNamedNodeMap* const attributes =
332 positionElement->getAttributes();
333 XMLSize_t attributeCount = attributes->getLength();
334
335 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
336 ++attribute_index)
337 {
338 xercesc::DOMNode* node = attributes->item(attribute_index);
339
340 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
341 {
342 continue;
343 }
344
345 const xercesc::DOMAttr* const attribute =
346 dynamic_cast<xercesc::DOMAttr*>(node);
347 if(attribute == nullptr)
348 {
349 G4Exception("G4GDMLRead::PositionRead()", "InvalidRead", FatalException,
350 "No attribute found!");
351 return;
352 }
353 const G4String attName = Transcode(attribute->getName());
354 const G4String attValue = Transcode(attribute->getValue());
355
356 if(attName == "name")
357 {
358 name = GenerateName(attValue);
359 }
360 else if(attName == "unit")
361 {
362 unit = G4UnitDefinition::GetValueOf(attValue);
363 if(G4UnitDefinition::GetCategory(attValue) != "Length")
364 {
365 G4Exception("G4GDMLReadDefine::PositionRead()", "InvalidRead",
366 FatalException, "Invalid unit for length!");
367 }
368 }
369 else if(attName == "x")
370 {
371 position.setX(eval.Evaluate(attValue));
372 }
373 else if(attName == "y")
374 {
375 position.setY(eval.Evaluate(attValue));
376 }
377 else if(attName == "z")
378 {
379 position.setZ(eval.Evaluate(attValue));
380 }
381 }
382
383 positionMap[name] = position * unit;
384}
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)
#define position
Definition: xmlparse.cc:622

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4UnitDefinition::GetCategory(), G4UnitDefinition::GetValueOf(), G4InuclParticleNames::name(), position, positionMap, and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ QuantityRead()

void G4GDMLReadDefine::QuantityRead ( const xercesc::DOMElement * const  element)
protected

Definition at line 547 of file G4GDMLReadDefine.cc.

548{
549 G4String name = "";
550 G4double unit = 1.0;
551 G4double value = 0.0;
552
553 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
554 XMLSize_t attributeCount = attributes->getLength();
555
556 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
557 ++attribute_index)
558 {
559 xercesc::DOMNode* node = attributes->item(attribute_index);
560
561 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
562 {
563 continue;
564 }
565
566 const xercesc::DOMAttr* const attribute =
567 dynamic_cast<xercesc::DOMAttr*>(node);
568 if(attribute == nullptr)
569 {
570 G4Exception("G4GDMLRead::QuantityRead()", "InvalidRead", FatalException,
571 "No attribute found!");
572 return;
573 }
574 const G4String attName = Transcode(attribute->getName());
575 const G4String attValue = Transcode(attribute->getValue());
576
577 if(attName == "name")
578 {
579 name = attValue;
580 }
581 else if(attName == "value")
582 {
583 value = eval.Evaluate(attValue);
584 }
585 else if(attName == "unit")
586 {
587 unit = G4UnitDefinition::GetValueOf(attValue);
588 }
589 }
590
591 quantityMap[name] = value * unit;
592 eval.DefineConstant(name, value * unit);
593}

References G4GDMLEvaluator::DefineConstant(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4UnitDefinition::GetValueOf(), G4InuclParticleNames::name(), quantityMap, and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ Read()

void G4GDMLRead::Read ( const G4String fileName,
G4bool  validation,
G4bool  isModule,
G4bool  strip = true 
)
inherited

Definition at line 422 of file G4GDMLRead.cc.

424{
425 dostrip = strip;
426#ifdef G4VERBOSE
427 if(isModule)
428 {
429 G4cout << "G4GDML: Reading module '" << fileName << "'..." << G4endl;
430 }
431 else
432 {
433 G4cout << "G4GDML: Reading '" << fileName << "'..." << G4endl;
434 }
435#endif
436 inLoop = 0;
437 validate = validation;
438
439 xercesc::ErrorHandler* handler = new G4GDMLErrorHandler(!validate);
440 xercesc::XercesDOMParser* parser = new xercesc::XercesDOMParser;
441
442 if(validate)
443 {
444 parser->setValidationScheme(xercesc::XercesDOMParser::Val_Always);
445 }
446 parser->setValidationSchemaFullChecking(validate);
447 parser->setCreateEntityReferenceNodes(false);
448 // Entities will be automatically resolved by Xerces
449
450 parser->setDoNamespaces(true);
451 parser->setDoSchema(validate);
452 parser->setErrorHandler(handler);
453
454 try
455 {
456 parser->parse(fileName.c_str());
457 } catch(const xercesc::XMLException& e)
458 {
459 G4cout << "G4GDML: " << Transcode(e.getMessage()) << G4endl;
460 } catch(const xercesc::DOMException& e)
461 {
462 G4cout << "G4GDML: " << Transcode(e.getMessage()) << G4endl;
463 }
464
465 xercesc::DOMDocument* doc = parser->getDocument();
466
467 if(doc == nullptr)
468 {
469 G4String error_msg = "Unable to open document: " + fileName;
470 G4Exception("G4GDMLRead::Read()", "InvalidRead", FatalException, error_msg);
471 return;
472 }
473 xercesc::DOMElement* element = doc->getDocumentElement();
474
475 if(element == nullptr )
476 {
477 std::ostringstream message;
478 message << "ERROR - Empty document or unable to validate schema!" << G4endl
479 << " Check Internet connection is ON in case of schema"
480 << G4endl
481 << " validation enabled and location defined as URL in"
482 << G4endl << " the GDML file - " << fileName
483 << " - being imported!" << G4endl
484 << " Otherwise, verify GDML schema server is reachable!";
485 G4Exception("G4GDMLRead::Read()", "InvalidRead", FatalException, message);
486 return;
487 }
488
489 for(xercesc::DOMNode* iter = element->getFirstChild(); iter != nullptr;
490 iter = iter->getNextSibling())
491 {
492 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
493 {
494 continue;
495 }
496
497 const xercesc::DOMElement* const child =
498 dynamic_cast<xercesc::DOMElement*>(iter);
499 if(child == nullptr)
500 {
501 G4Exception("G4GDMLRead::Read()", "InvalidRead", FatalException,
502 "No child found!");
503 return;
504 }
505 const G4String tag = Transcode(child->getTagName());
506
507 if(tag == "define")
508 {
509 DefineRead(child);
510 }
511 else if(tag == "materials")
512 {
513 MaterialsRead(child);
514 }
515 else if(tag == "solids")
516 {
517 SolidsRead(child);
518 }
519 else if(tag == "setup")
520 {
521 SetupRead(child);
522 }
523 else if(tag == "structure")
524 {
525 StructureRead(child);
526 }
527 else if(tag == "userinfo")
528 {
529 UserinfoRead(child);
530 }
531 else if(tag == "extension")
532 {
533 ExtensionRead(child);
534 }
535 else
536 {
537 G4String error_msg = "Unknown tag in gdml: " + tag;
538 G4Exception("G4GDMLRead::Read()", "InvalidRead", FatalException,
539 error_msg);
540 }
541 }
542
543 delete parser;
544 delete handler;
545
546 if(isModule)
547 {
548#ifdef G4VERBOSE
549 G4cout << "G4GDML: Reading module '" << fileName << "' done!" << G4endl;
550#endif
551 }
552 else
553 {
554 G4cout << "G4GDML: Reading '" << fileName << "' done!" << G4endl;
555 if(strip)
556 {
557 StripNames();
558 }
559 }
560}
G4bool validate
Definition: G4GDMLRead.hh:157
G4bool dostrip
Definition: G4GDMLRead.hh:159
virtual void SolidsRead(const xercesc::DOMElement *const)=0
virtual void MaterialsRead(const xercesc::DOMElement *const)=0
virtual void UserinfoRead(const xercesc::DOMElement *const)
Definition: G4GDMLRead.cc:377
virtual void SetupRead(const xercesc::DOMElement *const)=0
virtual void ExtensionRead(const xercesc::DOMElement *const)
Definition: G4GDMLRead.cc:414
virtual void DefineRead(const xercesc::DOMElement *const)=0
virtual void StructureRead(const xercesc::DOMElement *const)=0
void StripNames() const
Definition: G4GDMLRead.cc:122

References G4GDMLRead::DefineRead(), G4GDMLRead::dostrip, G4GDMLRead::ExtensionRead(), FatalException, G4cout, G4endl, G4Exception(), G4GDMLRead::inLoop, G4GDMLRead::MaterialsRead(), geant4_check_module_cycles::parser, G4GDMLRead::SetupRead(), G4GDMLRead::SolidsRead(), G4StrUtil::strip(), G4GDMLRead::StripNames(), G4GDMLRead::StructureRead(), G4GDMLRead::Transcode(), G4GDMLRead::UserinfoRead(), and G4GDMLRead::validate.

Referenced by G4GDMLReadStructure::FileRead().

◆ RefRead()

G4String G4GDMLReadDefine::RefRead ( const xercesc::DOMElement * const  element)
protected

Definition at line 714 of file G4GDMLReadDefine.cc.

715{
716 G4String ref;
717
718 const xercesc::DOMNamedNodeMap* const attributes = element->getAttributes();
719 XMLSize_t attributeCount = attributes->getLength();
720
721 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
722 ++attribute_index)
723 {
724 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
725
726 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
727 {
728 continue;
729 }
730
731 const xercesc::DOMAttr* const attribute =
732 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
733 if(attribute == nullptr)
734 {
735 G4Exception("G4GDMLRead::Read()", "InvalidRead", FatalException,
736 "No attribute found!");
737 return ref;
738 }
739 const G4String attName = Transcode(attribute->getName());
740 const G4String attValue = Transcode(attribute->getValue());
741
742 if(attName == "ref")
743 {
744 ref = attValue;
745 }
746 }
747
748 return ref;
749}

References FatalException, G4Exception(), and G4GDMLRead::Transcode().

Referenced by G4GDMLReadSolids::BooleanRead(), G4GDMLReadStructure::BorderSurfaceRead(), G4GDMLReadStructure::DivisionvolRead(), G4GDMLReadMaterials::MaterialRead(), G4GDMLReadSolids::MultiUnionNodeRead(), G4GDMLReadParamvol::ParametersRead(), G4GDMLReadParamvol::ParamvolRead(), G4GDMLReadStructure::PhysvolRead(), G4GDMLReadStructure::ReplicaRead(), G4GDMLReadStructure::ReplicavolRead(), G4GDMLReadSolids::ScaledSolidRead(), G4GDMLReadSetup::SetupRead(), G4GDMLReadStructure::SkinSurfaceRead(), and G4GDMLReadStructure::VolumeRead().

◆ RotationRead()

void G4GDMLReadDefine::RotationRead ( const xercesc::DOMElement * const  rotationElement)
protected

Definition at line 387 of file G4GDMLReadDefine.cc.

389{
390 G4String name = "";
391 G4double unit = 1.0;
392 G4ThreeVector rotation(0., 0., 0.);
393
394 const xercesc::DOMNamedNodeMap* const attributes =
395 rotationElement->getAttributes();
396 XMLSize_t attributeCount = attributes->getLength();
397
398 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
399 ++attribute_index)
400 {
401 xercesc::DOMNode* node = attributes->item(attribute_index);
402
403 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
404 {
405 continue;
406 }
407
408 const xercesc::DOMAttr* const attribute =
409 dynamic_cast<xercesc::DOMAttr*>(node);
410 if(attribute == nullptr)
411 {
412 G4Exception("G4GDMLRead::RotationRead()", "InvalidRead", FatalException,
413 "No attribute found!");
414 return;
415 }
416 const G4String attName = Transcode(attribute->getName());
417 const G4String attValue = Transcode(attribute->getValue());
418
419 if(attName == "name")
420 {
421 name = GenerateName(attValue);
422 }
423 else if(attName == "unit")
424 {
425 unit = G4UnitDefinition::GetValueOf(attValue);
426 if(G4UnitDefinition::GetCategory(attValue) != "Angle")
427 {
428 G4Exception("G4GDMLReadDefine::RotationRead()", "InvalidRead",
429 FatalException, "Invalid unit for angle!");
430 }
431 }
432 else if(attName == "x")
433 {
434 rotation.setX(eval.Evaluate(attValue));
435 }
436 else if(attName == "y")
437 {
438 rotation.setY(eval.Evaluate(attValue));
439 }
440 else if(attName == "z")
441 {
442 rotation.setZ(eval.Evaluate(attValue));
443 }
444 }
445
446 rotationMap[name] = rotation * unit;
447}

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4UnitDefinition::GetCategory(), G4UnitDefinition::GetValueOf(), G4InuclParticleNames::name(), rotationMap, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ ScaleRead()

void G4GDMLReadDefine::ScaleRead ( const xercesc::DOMElement * const  scaleElement)
protected

Definition at line 450 of file G4GDMLReadDefine.cc.

451{
452 G4String name = "";
453 G4ThreeVector scale(1.0, 1.0, 1.0);
454
455 const xercesc::DOMNamedNodeMap* const attributes =
456 scaleElement->getAttributes();
457 XMLSize_t attributeCount = attributes->getLength();
458
459 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
460 ++attribute_index)
461 {
462 xercesc::DOMNode* node = attributes->item(attribute_index);
463
464 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
465 {
466 continue;
467 }
468
469 const xercesc::DOMAttr* const attribute =
470 dynamic_cast<xercesc::DOMAttr*>(node);
471 if(attribute == nullptr)
472 {
473 G4Exception("G4GDMLRead::ScaleRead()", "InvalidRead", FatalException,
474 "No attribute found!");
475 return;
476 }
477 const G4String attName = Transcode(attribute->getName());
478 const G4String attValue = Transcode(attribute->getValue());
479
480 if(attName == "name")
481 {
482 name = GenerateName(attValue);
483 }
484 else if(attName == "x")
485 {
486 scale.setX(eval.Evaluate(attValue));
487 }
488 else if(attName == "y")
489 {
490 scale.setY(eval.Evaluate(attValue));
491 }
492 else if(attName == "z")
493 {
494 scale.setZ(eval.Evaluate(attValue));
495 }
496 }
497
498 scaleMap[name] = scale;
499}

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4InuclParticleNames::name(), scaleMap, CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ SetReverseSearch()

void G4GDMLReadDefine::SetReverseSearch ( G4bool  flag)
inline

Definition at line 81 of file G4GDMLReadDefine.hh.

82 { reverseSearch = flag; }

References reverseSearch.

◆ SetupRead()

virtual void G4GDMLRead::SetupRead ( const xercesc::DOMElement * const  )
pure virtualinherited

Implemented in G4GDMLReadSetup.

Referenced by G4GDMLRead::Read().

◆ SolidsRead()

virtual void G4GDMLRead::SolidsRead ( const xercesc::DOMElement * const  )
pure virtualinherited

◆ Strip()

G4String G4GDMLRead::Strip ( const G4String name) const
protectedinherited

◆ StripName()

void G4GDMLRead::StripName ( G4String name) const
inherited

Definition at line 112 of file G4GDMLRead.cc.

113{
114 auto idx = name.find("0x");
115 if(idx != G4String::npos)
116 {
117 name.erase(idx);
118 }
119}

References G4InuclParticleNames::name().

Referenced by G4GDMLParser::ExportRegions(), G4GDMLRead::GenerateName(), G4GDMLParser::ImportRegions(), G4GDMLRead::Strip(), and G4GDMLRead::StripNames().

◆ StripNames()

void G4GDMLRead::StripNames ( ) const
inherited

Definition at line 122 of file G4GDMLRead.cc.

123{
124 // Strips off names of volumes, solids elements and materials from possible
125 // reference pointers or IDs attached to their original identifiers.
126
130 const G4ElementTable* elements = G4Element::GetElementTable();
132
133 G4cout << "Stripping off GDML names of materials, solids and volumes ..."
134 << G4endl;
135
136 G4String sname;
137 std::size_t i;
138
139 // Solids...
140 //
141 for(i = 0; i < solids->size(); ++i)
142 {
143 G4VSolid* psol = (*solids)[i];
144 sname = psol->GetName();
145 StripName(sname);
146 psol->SetName(sname);
147 }
148 solids->UpdateMap();
149
150 // Logical volumes...
151 //
152 for(i = 0; i < lvols->size(); ++i)
153 {
154 G4LogicalVolume* lvol = (*lvols)[i];
155 sname = lvol->GetName();
156 StripName(sname);
157 lvol->SetName(sname);
158 }
159 lvols->UpdateMap();
160
161 // Physical volumes...
162 //
163 for(i = 0; i < pvols->size(); ++i)
164 {
165 G4VPhysicalVolume* pvol = (*pvols)[i];
166 sname = pvol->GetName();
167 StripName(sname);
168 pvol->SetName(sname);
169 }
170 pvols->UpdateMap();
171
172 // Materials...
173 //
174 for(i = 0; i < materials->size(); ++i)
175 {
176 G4Material* pmat = (*materials)[i];
177 sname = pmat->GetName();
178 StripName(sname);
179 pmat->SetName(sname);
180 }
181
182 // Elements...
183 //
184 for(i = 0; i < elements->size(); ++i)
185 {
186 G4Element* pelm = (*elements)[i];
187 sname = pelm->GetName();
188 StripName(sname);
189 pelm->SetName(sname);
190 }
191}
std::vector< G4Element * > G4ElementTable
std::vector< G4Material * > G4MaterialTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
void SetName(const G4String &name)
Definition: G4Element.hh:215
const G4String & GetName() const
Definition: G4Element.hh:127
static G4LogicalVolumeStore * GetInstance()
void SetName(const G4String &pName)
void SetName(const G4String &name)
Definition: G4Material.hh:285
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:672
const G4String & GetName() const
Definition: G4Material.hh:173
static G4PhysicalVolumeStore * GetInstance()
void UpdateMap()
static G4SolidStore * GetInstance()
const G4String & GetName() const
G4String GetName() const
void SetName(const G4String &name)
Definition: G4VSolid.cc:127

References G4cout, G4endl, G4Element::GetElementTable(), G4LogicalVolumeStore::GetInstance(), G4PhysicalVolumeStore::GetInstance(), G4SolidStore::GetInstance(), G4Material::GetMaterialTable(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), G4VSolid::GetName(), G4Element::GetName(), G4Material::GetName(), G4VSolid::SetName(), G4Element::SetName(), G4Material::SetName(), G4LogicalVolume::SetName(), G4VPhysicalVolume::SetName(), G4GDMLRead::StripName(), G4LogicalVolumeStore::UpdateMap(), G4PhysicalVolumeStore::UpdateMap(), and G4SolidStore::UpdateMap().

Referenced by G4GDMLRead::Read().

◆ StructureRead()

virtual void G4GDMLRead::StructureRead ( const xercesc::DOMElement * const  )
pure virtualinherited

◆ Transcode()

G4String G4GDMLRead::Transcode ( const XMLCh * const  toTranscode)
protectedinherited

Definition at line 55 of file G4GDMLRead.cc.

56{
57 char* char_str = xercesc::XMLString::transcode(toTranscode);
58 G4String my_str(char_str);
59 xercesc::XMLString::release(&char_str);
60 return my_str;
61}

Referenced by G4GDMLReadStructure::AssemblyRead(), G4GDMLReadMaterials::AtomRead(), G4GDMLRead::AuxiliaryRead(), G4GDMLReadStructure::AxisRead(), G4GDMLReadSolids::BooleanRead(), G4GDMLReadStructure::BorderSurfaceRead(), G4GDMLReadParamvol::Box_dimensionsRead(), G4GDMLReadSolids::BoxRead(), G4GDMLReadMaterials::CompositeRead(), G4GDMLReadParamvol::Cone_dimensionsRead(), G4GDMLReadSolids::ConeRead(), ConstantRead(), G4GDMLReadSolids::CutTubeRead(), DefineRead(), G4GDMLReadStructure::DivisionvolRead(), G4GDMLReadMaterials::DRead(), G4GDMLReadSolids::ElconeRead(), G4GDMLReadMaterials::ElementRead(), G4GDMLReadParamvol::Ellipsoid_dimensionsRead(), G4GDMLReadSolids::EllipsoidRead(), G4GDMLReadSolids::EltubeRead(), ExpressionRead(), G4GDMLReadStructure::FileRead(), G4GDMLReadMaterials::FractionRead(), G4GDMLReadSolids::GenericPolyconeRead(), G4GDMLReadSolids::GenericPolyhedraRead(), G4GDMLReadSolids::GenTrapRead(), G4GDMLReadParamvol::Hype_dimensionsRead(), G4GDMLReadSolids::HypeRead(), G4GDMLReadMaterials::IsotopeRead(), G4GDMLRead::LoopRead(), G4GDMLReadMaterials::MaterialRead(), G4GDMLReadMaterials::MaterialsRead(), MatrixRead(), G4GDMLReadMaterials::MEERead(), G4GDMLReadMaterials::MixtureRead(), G4GDMLReadSolids::MultiUnionNodeRead(), G4GDMLReadSolids::MultiUnionRead(), G4GDMLReadSolids::OpticalSurfaceRead(), G4GDMLReadParamvol::Orb_dimensionsRead(), G4GDMLReadSolids::OrbRead(), G4GDMLReadParamvol::Para_dimensionsRead(), G4GDMLReadSolids::ParaboloidRead(), G4GDMLReadParamvol::ParameterisedRead(), G4GDMLReadParamvol::ParametersRead(), G4GDMLReadParamvol::Paramvol_contentRead(), G4GDMLReadParamvol::ParamvolRead(), G4GDMLReadSolids::ParaRead(), G4GDMLReadStructure::PhysvolRead(), G4GDMLReadParamvol::Polycone_dimensionsRead(), G4GDMLReadSolids::PolyconeRead(), G4GDMLReadParamvol::Polyhedra_dimensionsRead(), G4GDMLReadSolids::PolyhedraRead(), PositionRead(), G4GDMLReadMaterials::PRead(), G4GDMLReadMaterials::PropertyRead(), G4GDMLReadSolids::PropertyRead(), G4GDMLReadSolids::QuadrangularRead(), G4GDMLReadStructure::QuantityRead(), QuantityRead(), G4GDMLRead::Read(), G4GDMLReadSolids::ReflectedSolidRead(), RefRead(), G4GDMLReadStructure::ReplicaRead(), G4GDMLReadStructure::ReplicavolRead(), RotationRead(), G4GDMLReadSolids::RZPointRead(), G4GDMLReadSolids::ScaledSolidRead(), ScaleRead(), G4GDMLReadSolids::SectionRead(), G4GDMLReadSetup::SetupRead(), G4GDMLReadStructure::SkinSurfaceRead(), G4GDMLReadSolids::SolidsRead(), G4GDMLReadParamvol::Sphere_dimensionsRead(), G4GDMLReadSolids::SphereRead(), G4GDMLReadStructure::StructureRead(), G4GDMLReadSolids::TessellatedRead(), G4GDMLReadSolids::TetRead(), G4GDMLReadParamvol::Torus_dimensionsRead(), G4GDMLReadSolids::TorusRead(), G4GDMLReadParamvol::Trap_dimensionsRead(), G4GDMLReadSolids::TrapRead(), G4GDMLReadParamvol::Trd_dimensionsRead(), G4GDMLReadSolids::TrdRead(), G4GDMLReadMaterials::TRead(), G4GDMLReadSolids::TriangularRead(), G4GDMLReadParamvol::Tube_dimensionsRead(), G4GDMLReadSolids::TubeRead(), G4GDMLReadSolids::TwistedboxRead(), G4GDMLReadSolids::TwistedtrapRead(), G4GDMLReadSolids::TwistedtrdRead(), G4GDMLReadSolids::TwistedtubsRead(), G4GDMLReadSolids::TwoDimVertexRead(), G4GDMLRead::UserinfoRead(), VariableRead(), VectorRead(), G4GDMLReadStructure::Volume_contentRead(), G4GDMLReadStructure::VolumeRead(), G4GDMLReadSolids::XtruRead(), and G4GDMLReadSolids::ZplaneRead().

◆ UserinfoRead()

void G4GDMLRead::UserinfoRead ( const xercesc::DOMElement * const  userinfoElement)
virtualinherited

Definition at line 377 of file G4GDMLRead.cc.

378{
379#ifdef G4VERBOSE
380 G4cout << "G4GDML: Reading userinfo..." << G4endl;
381#endif
382 for(xercesc::DOMNode* iter = userinfoElement->getFirstChild();
383 iter != nullptr; iter = iter->getNextSibling())
384 {
385 if(iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE)
386 {
387 continue;
388 }
389
390 const xercesc::DOMElement* const child =
391 dynamic_cast<xercesc::DOMElement*>(iter);
392 if(!child)
393 {
394 G4Exception("G4GDMLRead::UserinfoRead()", "InvalidRead", FatalException,
395 "No child found!");
396 return;
397 }
398 const G4String tag = Transcode(child->getTagName());
399
400 if(tag == "auxiliary")
401 {
402 auxGlobalList.push_back(AuxiliaryRead(child));
403 }
404 else
405 {
406 G4String error_msg = "Unknown tag in structure: " + tag;
407 G4Exception("G4GDMLRead::UserinfoRead()", "ReadError", FatalException,
408 error_msg);
409 }
410 }
411}

References G4GDMLRead::auxGlobalList, G4GDMLRead::AuxiliaryRead(), FatalException, G4cout, G4endl, G4Exception(), and G4GDMLRead::Transcode().

Referenced by G4GDMLRead::Read().

◆ VariableRead()

void G4GDMLReadDefine::VariableRead ( const xercesc::DOMElement * const  variableElement)
protected

Definition at line 502 of file G4GDMLReadDefine.cc.

504{
505 G4String name = "";
506 G4double value = 0.0;
507
508 const xercesc::DOMNamedNodeMap* const attributes =
509 variableElement->getAttributes();
510 XMLSize_t attributeCount = attributes->getLength();
511
512 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
513 ++attribute_index)
514 {
515 xercesc::DOMNode* node = attributes->item(attribute_index);
516
517 if(node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
518 {
519 continue;
520 }
521
522 const xercesc::DOMAttr* const attribute =
523 dynamic_cast<xercesc::DOMAttr*>(node);
524 if(attribute == nullptr)
525 {
526 G4Exception("G4GDMLRead::VariableRead()", "InvalidRead", FatalException,
527 "No attribute found!");
528 return;
529 }
530 const G4String attName = Transcode(attribute->getName());
531 const G4String attValue = Transcode(attribute->getValue());
532
533 if(attName == "name")
534 {
535 name = attValue;
536 }
537 else if(attName == "value")
538 {
539 value = eval.Evaluate(attValue);
540 }
541 }
542
543 eval.DefineVariable(name, value);
544}
void DefineVariable(const G4String &, G4double)

References G4GDMLEvaluator::DefineVariable(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4InuclParticleNames::name(), and G4GDMLRead::Transcode().

Referenced by DefineRead().

◆ VectorRead()

void G4GDMLReadDefine::VectorRead ( const xercesc::DOMElement * const  vectorElement,
G4ThreeVector vec 
)
protected

Definition at line 662 of file G4GDMLReadDefine.cc.

664{
665 G4double unit = 1.0;
666
667 const xercesc::DOMNamedNodeMap* const attributes =
668 vectorElement->getAttributes();
669 XMLSize_t attributeCount = attributes->getLength();
670
671 for(XMLSize_t attribute_index = 0; attribute_index < attributeCount;
672 ++attribute_index)
673 {
674 xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
675
676 if(attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
677 {
678 continue;
679 }
680
681 const xercesc::DOMAttr* const attribute =
682 dynamic_cast<xercesc::DOMAttr*>(attribute_node);
683 if(attribute == nullptr)
684 {
685 G4Exception("G4GDMLRead::VectorRead()", "InvalidRead", FatalException,
686 "No attribute found!");
687 return;
688 }
689 const G4String attName = Transcode(attribute->getName());
690 const G4String attValue = Transcode(attribute->getValue());
691
692 if(attName == "unit")
693 {
694 unit = G4UnitDefinition::GetValueOf(attValue);
695 }
696 else if(attName == "x")
697 {
698 vec.setX(eval.Evaluate(attValue));
699 }
700 else if(attName == "y")
701 {
702 vec.setY(eval.Evaluate(attValue));
703 }
704 else if(attName == "z")
705 {
706 vec.setZ(eval.Evaluate(attValue));
707 }
708 }
709
710 vec *= unit;
711}
void setY(double)
void setZ(double)
void setX(double)

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), G4UnitDefinition::GetValueOf(), CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), and G4GDMLRead::Transcode().

Referenced by G4GDMLReadSolids::BooleanRead(), G4GDMLReadSolids::MultiUnionNodeRead(), G4GDMLReadParamvol::ParametersRead(), G4GDMLReadStructure::PhysvolRead(), G4GDMLReadStructure::ReplicaRead(), and G4GDMLReadSolids::ScaledSolidRead().

◆ Volume_contentRead()

virtual void G4GDMLRead::Volume_contentRead ( const xercesc::DOMElement * const  )
pure virtualinherited

Field Documentation

◆ auxGlobalList

G4GDMLAuxListType G4GDMLRead::auxGlobalList
privateinherited

Definition at line 164 of file G4GDMLRead.hh.

Referenced by G4GDMLRead::GetAuxList(), and G4GDMLRead::UserinfoRead().

◆ check

G4bool G4GDMLRead::check = false
protectedinherited

◆ dostrip

G4bool G4GDMLRead::dostrip = true
protectedinherited

Definition at line 159 of file G4GDMLRead.hh.

Referenced by G4GDMLReadStructure::GetWorldVolume(), and G4GDMLRead::Read().

◆ eval

G4GDMLEvaluator G4GDMLRead::eval
protectedinherited

Definition at line 156 of file G4GDMLRead.hh.

Referenced by G4GDMLReadMaterials::AtomRead(), G4GDMLReadStructure::AxisRead(), G4GDMLReadParamvol::Box_dimensionsRead(), G4GDMLReadSolids::BoxRead(), G4GDMLReadStructure::Clear(), G4GDMLReadMaterials::CompositeRead(), G4GDMLReadParamvol::Cone_dimensionsRead(), G4GDMLReadSolids::ConeRead(), ConstantRead(), G4GDMLReadSolids::CutTubeRead(), G4GDMLReadStructure::DivisionvolRead(), G4GDMLReadMaterials::DRead(), G4GDMLReadSolids::ElconeRead(), G4GDMLReadMaterials::ElementRead(), G4GDMLReadParamvol::Ellipsoid_dimensionsRead(), G4GDMLReadSolids::EllipsoidRead(), G4GDMLReadSolids::EltubeRead(), ExpressionRead(), G4GDMLReadMaterials::FractionRead(), G4GDMLRead::GenerateName(), G4GDMLRead::GeneratePhysvolName(), G4GDMLReadSolids::GenericPolyconeRead(), G4GDMLReadSolids::GenericPolyhedraRead(), G4GDMLReadSolids::GenTrapRead(), GetConstant(), GetVariable(), G4GDMLReadParamvol::Hype_dimensionsRead(), G4GDMLReadSolids::HypeRead(), G4GDMLReadMaterials::IsotopeRead(), IsValidID(), G4GDMLRead::LoopRead(), G4GDMLReadMaterials::MaterialRead(), MatrixRead(), G4GDMLReadMaterials::MEERead(), G4GDMLReadSolids::OpticalSurfaceRead(), G4GDMLReadParamvol::Orb_dimensionsRead(), G4GDMLReadSolids::OrbRead(), G4GDMLReadParamvol::Para_dimensionsRead(), G4GDMLReadSolids::ParaboloidRead(), G4GDMLReadParamvol::ParameterisedRead(), G4GDMLReadSolids::ParaRead(), G4GDMLReadStructure::PhysvolRead(), G4GDMLReadParamvol::Polycone_dimensionsRead(), G4GDMLReadSolids::PolyconeRead(), G4GDMLReadParamvol::Polyhedra_dimensionsRead(), G4GDMLReadSolids::PolyhedraRead(), PositionRead(), G4GDMLReadMaterials::PRead(), G4GDMLReadStructure::QuantityRead(), QuantityRead(), G4GDMLReadSolids::ReflectedSolidRead(), RotationRead(), G4GDMLReadSolids::RZPointRead(), ScaleRead(), G4GDMLReadSolids::SectionRead(), G4GDMLReadParamvol::Sphere_dimensionsRead(), G4GDMLReadSolids::SphereRead(), G4GDMLReadParamvol::Torus_dimensionsRead(), G4GDMLReadSolids::TorusRead(), G4GDMLReadParamvol::Trap_dimensionsRead(), G4GDMLReadSolids::TrapRead(), G4GDMLReadParamvol::Trd_dimensionsRead(), G4GDMLReadSolids::TrdRead(), G4GDMLReadMaterials::TRead(), G4GDMLReadParamvol::Tube_dimensionsRead(), G4GDMLReadSolids::TubeRead(), G4GDMLReadSolids::TwistedboxRead(), G4GDMLReadSolids::TwistedtrapRead(), G4GDMLReadSolids::TwistedtrdRead(), G4GDMLReadSolids::TwistedtubsRead(), G4GDMLReadSolids::TwoDimVertexRead(), VariableRead(), VectorRead(), G4GDMLReadStructure::Volume_contentRead(), and G4GDMLReadSolids::ZplaneRead().

◆ inLoop

G4int G4GDMLRead::inLoop = 0
privateinherited

Definition at line 163 of file G4GDMLRead.hh.

Referenced by G4GDMLRead::GenerateName(), G4GDMLRead::LoopRead(), and G4GDMLRead::Read().

◆ loopCount

G4int G4GDMLRead::loopCount = 0
privateinherited

Definition at line 163 of file G4GDMLRead.hh.

Referenced by G4GDMLRead::LoopRead().

◆ matrixMap

std::map<G4String, G4GDMLMatrix> G4GDMLReadDefine::matrixMap
protected

Definition at line 108 of file G4GDMLReadDefine.hh.

Referenced by GetMatrix(), and MatrixRead().

◆ positionMap

std::map<G4String, G4ThreeVector> G4GDMLReadDefine::positionMap
protected

Definition at line 105 of file G4GDMLReadDefine.hh.

Referenced by GetPosition(), and PositionRead().

◆ quantityMap

std::map<G4String, G4double> G4GDMLReadDefine::quantityMap
protected

Definition at line 104 of file G4GDMLReadDefine.hh.

Referenced by GetQuantity(), and QuantityRead().

◆ reverseSearch

G4bool G4GDMLReadDefine::reverseSearch = false
protected

◆ rotationMap

std::map<G4String, G4ThreeVector> G4GDMLReadDefine::rotationMap
protected

Definition at line 106 of file G4GDMLReadDefine.hh.

Referenced by GetRotation(), and RotationRead().

◆ scaleMap

std::map<G4String, G4ThreeVector> G4GDMLReadDefine::scaleMap
protected

Definition at line 107 of file G4GDMLReadDefine.hh.

Referenced by GetScale(), and ScaleRead().

◆ validate

G4bool G4GDMLRead::validate = true
protectedinherited

Definition at line 157 of file G4GDMLRead.hh.

Referenced by G4GDMLReadStructure::FileRead(), and G4GDMLRead::Read().


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