Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4GDMLReadStructure Class Reference

#include <G4GDMLReadStructure.hh>

Inheritance diagram for G4GDMLReadStructure:
G4GDMLReadParamvol G4GDMLReadSetup G4GDMLReadSolids G4GDMLReadMaterials G4GDMLReadDefine G4GDMLRead G03ColorReader

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)
 
- Public Member Functions inherited from G4GDMLRead
virtual void ExtensionRead (const xercesc::DOMElement *const)
 
void Read (const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
 
void StripNames () const
 
void OverlapCheck (G4bool)
 

Protected Member Functions

G4GDMLAuxPairType AuxiliaryRead (const xercesc::DOMElement *const)
 
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 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)
 
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)
 
- 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 StripName (G4String &) const
 
void GeneratePhysvolName (const G4String &, G4VPhysicalVolume *)
 
void LoopRead (const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
 

Protected Attributes

G4GDMLAuxMapType auxMap
 
G4GDMLAssemblyMapType assemblyMap
 
G4LogicalVolumepMotherLogical
 
std::map< std::string,
G4VPhysicalVolume * > 
setuptoPV
 
- Protected Attributes inherited from G4GDMLReadParamvol
G4GDMLParameterisationparameterisation
 
- Protected Attributes inherited from G4GDMLReadSetup
std::map< G4String, G4StringsetupMap
 
- Protected Attributes inherited from G4GDMLReadDefine
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
 
G4bool check
 

Detailed Description

Definition at line 62 of file G4GDMLReadStructure.hh.

Constructor & Destructor Documentation

G4GDMLReadStructure::G4GDMLReadStructure ( )

Definition at line 48 of file G4GDMLReadStructure.cc.

50 {
51 }
G4LogicalVolume * pMotherLogical
G4GDMLReadStructure::~G4GDMLReadStructure ( )
virtual

Definition at line 53 of file G4GDMLReadStructure.cc.

54 {
55 }

Member Function Documentation

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

Definition at line 635 of file G4GDMLReadStructure.cc.

References assemblyMap, FatalException, G4cout, G4endl, G4Exception(), G4GDMLRead::GenerateName(), PhysvolRead(), and G4GDMLRead::Transcode().

Referenced by StructureRead().

636 {
637  XMLCh *name_attr = xercesc::XMLString::transcode("name");
638  const G4String name = Transcode(assemblyElement->getAttribute(name_attr));
639  xercesc::XMLString::release(&name_attr);
640 
641  G4AssemblyVolume* pAssembly = new G4AssemblyVolume();
642  assemblyMap.insert(std::make_pair(GenerateName(name), pAssembly));
643 
644  for (xercesc::DOMNode* iter = assemblyElement->getFirstChild();
645  iter != 0; iter = iter->getNextSibling())
646  {
647  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
648  const xercesc::DOMElement* const child
649  = dynamic_cast<xercesc::DOMElement*>(iter);
650  if (!child)
651  {
652  G4Exception("G4GDMLReadStructure::AssemblyRead()",
653  "InvalidRead", FatalException, "No child found!");
654  return;
655  }
656  const G4String tag = Transcode(child->getTagName());
657 
658  if (tag=="physvol")
659  {
660  PhysvolRead(child, pAssembly);
661  }
662  else
663  {
664  G4cout << "Unsupported GDML tag '" << tag
665  << "' for Geant4 assembly structure !" << G4endl;
666  }
667  }
668 }
Definition: xmlparse.cc:179
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
G4GLOB_DLL std::ostream G4cout
G4GDMLAssemblyMapType assemblyMap
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
#define G4endl
Definition: G4ios.hh:61
void PhysvolRead(const xercesc::DOMElement *const, G4AssemblyVolume *assembly=0)
G4GDMLAuxPairType G4GDMLReadStructure::AuxiliaryRead ( const xercesc::DOMElement * const  auxiliaryElement)
protected

Definition at line 58 of file G4GDMLReadStructure.cc.

References FatalException, G4Exception(), G4GDMLRead::Transcode(), G4GDMLAuxPairType::type, and G4GDMLAuxPairType::value.

Referenced by G03ColorReader::VolumeRead(), and VolumeRead().

59 {
60  G4GDMLAuxPairType auxpair = {"",""};
61 
62  const xercesc::DOMNamedNodeMap* const attributes
63  = auxiliaryElement->getAttributes();
64  XMLSize_t attributeCount = attributes->getLength();
65 
66  for (XMLSize_t attribute_index=0;
67  attribute_index<attributeCount; attribute_index++)
68  {
69  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
70 
71  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
72  { continue; }
73 
74  const xercesc::DOMAttr* const attribute
75  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
76  if (!attribute)
77  {
78  G4Exception("G4GDMLReadStructure::AuxiliaryRead()",
79  "InvalidRead", FatalException, "No attribute found!");
80  return auxpair;
81  }
82  const G4String attName = Transcode(attribute->getName());
83  const G4String attValue = Transcode(attribute->getValue());
84 
85  if (attName=="auxtype") { auxpair.type = attValue; } else
86  // if (attName=="auxvalue") { auxpair.value = eval.Evaluate(attValue); }
87  if (attName=="auxvalue") { auxpair.value = attValue; }
88  }
89 
90  return auxpair;
91 }
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
EAxis G4GDMLReadStructure::AxisRead ( const xercesc::DOMElement *const  axisElement)
protected

Definition at line 517 of file G4GDMLReadStructure.cc.

References G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), FatalException, G4Exception(), kPhi, kRho, kUndefined, kXAxis, kYAxis, kZAxis, and G4GDMLRead::Transcode().

Referenced by ReplicaRead().

518 {
519  EAxis axis = kUndefined;
520 
521  const xercesc::DOMNamedNodeMap* const attributes
522  = axisElement->getAttributes();
523  XMLSize_t attributeCount = attributes->getLength();
524 
525  for (XMLSize_t attribute_index=0;
526  attribute_index<attributeCount; attribute_index++)
527  {
528  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
529 
530  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
531  { continue; }
532 
533  const xercesc::DOMAttr* const attribute
534  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
535  if (!attribute)
536  {
537  G4Exception("G4GDMLReadStructure::AxisRead()",
538  "InvalidRead", FatalException, "No attribute found!");
539  return axis;
540  }
541  const G4String attName = Transcode(attribute->getName());
542  const G4String attValue = Transcode(attribute->getValue());
543  if (attName=="x")
544  { if( eval.Evaluate(attValue)==1.) {axis=kXAxis;} }
545  else if (attName=="y")
546  { if( eval.Evaluate(attValue)==1.) {axis=kYAxis;} }
547  else if (attName=="z")
548  { if( eval.Evaluate(attValue)==1.) {axis=kZAxis;} }
549  else if (attName=="rho")
550  { if( eval.Evaluate(attValue)==1.) {axis=kRho;} }
551  else if (attName=="phi")
552  { if( eval.Evaluate(attValue)==1.) {axis=kPhi;} }
553  }
554 
555  return axis;
556 }
Definition: geomdefs.hh:54
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:144
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
EAxis
Definition: geomdefs.hh:54
Definition: geomdefs.hh:54
G4double Evaluate(const G4String &)
void G4GDMLReadStructure::BorderSurfaceRead ( const xercesc::DOMElement * const  bordersurfaceElement)
protected

Definition at line 94 of file G4GDMLReadStructure.cc.

References FatalException, G4Exception(), G4GDMLRead::GenerateName(), GetPhysvol(), G4GDMLReadSolids::GetSurfaceProperty(), G4GDMLReadDefine::RefRead(), G4GDMLRead::Strip(), and G4GDMLRead::Transcode().

Referenced by StructureRead().

95 {
96  G4String name;
97  G4VPhysicalVolume* pv1 = 0;
98  G4VPhysicalVolume* pv2 = 0;
99  G4SurfaceProperty* prop = 0;
100  G4int index = 0;
101 
102  const xercesc::DOMNamedNodeMap* const attributes
103  = bordersurfaceElement->getAttributes();
104  XMLSize_t attributeCount = attributes->getLength();
105 
106  for (XMLSize_t attribute_index=0;
107  attribute_index<attributeCount; attribute_index++)
108  {
109  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
110 
111  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
112  { continue; }
113 
114  const xercesc::DOMAttr* const attribute
115  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
116  if (!attribute)
117  {
118  G4Exception("G4GDMLReadStructure::BorderSurfaceRead()",
119  "InvalidRead", FatalException, "No attribute found!");
120  return;
121  }
122  const G4String attName = Transcode(attribute->getName());
123  const G4String attValue = Transcode(attribute->getValue());
124 
125  if (attName=="name")
126  { name = GenerateName(attValue); } else
127  if (attName=="surfaceproperty")
128  { prop = GetSurfaceProperty(GenerateName(attValue)); }
129  }
130 
131  for (xercesc::DOMNode* iter = bordersurfaceElement->getFirstChild();
132  iter != 0; iter = iter->getNextSibling())
133  {
134  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
135 
136  const xercesc::DOMElement* const child
137  = dynamic_cast<xercesc::DOMElement*>(iter);
138  if (!child)
139  {
140  G4Exception("G4GDMLReadStructure::BorderSurfaceRead()",
141  "InvalidRead", FatalException, "No child found!");
142  return;
143  }
144  const G4String tag = Transcode(child->getTagName());
145 
146  if (tag != "physvolref") { continue; }
147 
148  if (index==0)
149  { pv1 = GetPhysvol(GenerateName(RefRead(child))); index++; } else
150  if (index==1)
151  { pv2 = GetPhysvol(GenerateName(RefRead(child))); index++; } else
152  break;
153  }
154 
155  new G4LogicalBorderSurface(Strip(name),pv1,pv2,prop);
156 }
Definition: xmlparse.cc:179
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
G4SurfaceProperty * GetSurfaceProperty(const G4String &) const
int G4int
Definition: G4Types.hh:78
G4VPhysicalVolume * GetPhysvol(const G4String &) const
G4String RefRead(const xercesc::DOMElement *const)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:100
void G4GDMLReadStructure::Clear ( )

Definition at line 923 of file G4GDMLReadStructure.cc.

References G4GDMLEvaluator::Clear(), G4GDMLRead::eval, and setuptoPV.

924 {
925  eval.Clear();
926  setuptoPV.clear();
927 }
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:144
std::map< std::string, G4VPhysicalVolume * > setuptoPV
void G4GDMLReadStructure::DivisionvolRead ( const xercesc::DOMElement * const  divisionvolElement)
protected

Definition at line 159 of file G4GDMLReadStructure.cc.

References G4ReflectionFactory::Divide(), G4GDMLRead::eval, G4GDMLEvaluator::Evaluate(), G4GDMLEvaluator::EvaluateInteger(), FatalException, G4String::first(), G4Exception(), G4GDMLRead::GenerateName(), G4GDMLRead::GeneratePhysvolName(), G4PVDivisionFactory::GetInstance(), G4LogicalVolume::GetName(), GetVolume(), G4ReflectionFactory::Instance(), kPhi, kRho, kUndefined, kXAxis, kYAxis, kZAxis, pMotherLogical, G4GDMLReadDefine::RefRead(), G4GDMLRead::Transcode(), and width.

Referenced by Volume_contentRead().

160 {
161  G4String name;
162  G4double unit = 1.0;
163  G4double width = 0.0;
164  G4double offset = 0.0;
165  G4int number = 0;
166  EAxis axis = kUndefined;
167  G4LogicalVolume* logvol = 0;
168 
169  const xercesc::DOMNamedNodeMap* const attributes
170  = divisionvolElement->getAttributes();
171  XMLSize_t attributeCount = attributes->getLength();
172 
173  for (XMLSize_t attribute_index=0;
174  attribute_index<attributeCount; attribute_index++)
175  {
176  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
177 
178  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
179  { continue; }
180 
181  const xercesc::DOMAttr* const attribute
182  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
183  if (!attribute)
184  {
185  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
186  "InvalidRead", FatalException, "No attribute found!");
187  return;
188  }
189  const G4String attName = Transcode(attribute->getName());
190  const G4String attValue = Transcode(attribute->getValue());
191 
192  if (attName=="name") { name = attValue; } else
193  if (attName=="unit") { unit = eval.Evaluate(attValue); } else
194  if (attName=="width") { width = eval.Evaluate(attValue); } else
195  if (attName=="offset") { offset = eval.Evaluate(attValue); } else
196  if (attName=="number") { number = eval.EvaluateInteger(attValue); } else
197  if (attName=="axis")
198  {
199  if (attValue=="kXAxis") { axis = kXAxis; } else
200  if (attValue=="kYAxis") { axis = kYAxis; } else
201  if (attValue=="kZAxis") { axis = kZAxis; } else
202  if (attValue=="kRho") { axis = kRho; } else
203  if (attValue=="kPhi") { axis = kPhi; }
204  }
205  }
206 
207  width *= unit;
208  offset *= unit;
209 
210  for (xercesc::DOMNode* iter = divisionvolElement->getFirstChild();
211  iter != 0;iter = iter->getNextSibling())
212  {
213  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
214 
215  const xercesc::DOMElement* const child
216  = dynamic_cast<xercesc::DOMElement*>(iter);
217  if (!child)
218  {
219  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
220  "InvalidRead", FatalException, "No child found!");
221  return;
222  }
223  const G4String tag = Transcode(child->getTagName());
224 
225  if (tag=="volumeref") { logvol = GetVolume(GenerateName(RefRead(child))); }
226  }
227 
228  if (!logvol) { return; }
229 
232 
233  G4String pv_name = logvol->GetName() + "_div";
234  if ((number != 0) && (width == 0.0))
235  {
237  ->Divide(pv_name,logvol,pMotherLogical,axis,number,offset);
238  }
239  else if ((number == 0) && (width != 0.0))
240  {
242  ->Divide(pv_name,logvol,pMotherLogical,axis,width,offset);
243  }
244  else
245  {
247  ->Divide(pv_name,logvol,pMotherLogical,axis,number,width,offset);
248  }
249 
250  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
251  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
252 }
G4PhysicalVolumesPair Divide(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofDivisions, G4double width, G4double offset)
Definition: geomdefs.hh:54
G4int EvaluateInteger(const G4String &)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:144
G4String GetName() const
G4int first(char) const
Definition: xmlparse.cc:179
#define width
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:84
int G4int
Definition: G4Types.hh:78
static G4ReflectionFactory * Instance()
G4String RefRead(const xercesc::DOMElement *const)
static G4PVDivisionFactory * GetInstance()
G4LogicalVolume * GetVolume(const G4String &) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
EAxis
Definition: geomdefs.hh:54
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
double G4double
Definition: G4Types.hh:76
Definition: geomdefs.hh:54
G4double Evaluate(const G4String &)
G4LogicalVolume * G4GDMLReadStructure::FileRead ( const xercesc::DOMElement * const  fileElement)
protected

Definition at line 255 of file G4GDMLReadStructure.cc.

References auxMap, FatalException, G4Exception(), G4GDMLRead::GenerateName(), GetAuxMap(), G4GDMLReadSetup::GetSetup(), GetVolume(), G4GDMLRead::Read(), G4GDMLRead::Transcode(), and G4GDMLRead::validate.

Referenced by PhysvolRead().

256 {
257  G4String name;
258  G4String volname;
259 
260  const xercesc::DOMNamedNodeMap* const attributes
261  = fileElement->getAttributes();
262  XMLSize_t attributeCount = attributes->getLength();
263 
264  for (XMLSize_t attribute_index=0;
265  attribute_index<attributeCount; attribute_index++)
266  {
267  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
268 
269  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
270  { continue; }
271 
272  const xercesc::DOMAttr* const attribute
273  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
274  if (!attribute)
275  {
276  G4Exception("G4GDMLReadStructure::FileRead()",
277  "InvalidRead", FatalException, "No attribute found!");
278  return 0;
279  }
280  const G4String attName = Transcode(attribute->getName());
281  const G4String attValue = Transcode(attribute->getValue());
282 
283  if (attName=="name") { name = attValue; } else
284  if (attName=="volname") { volname = attValue; }
285  }
286 
287  const G4bool isModule = true;
288  G4GDMLReadStructure structure;
289  structure.Read(name,validate,isModule);
290 
291  // Register existing auxiliar information defined in child module
292  //
293  const G4GDMLAuxMapType* aux = structure.GetAuxMap();
294  if (!aux->empty())
295  {
296  G4GDMLAuxMapType::const_iterator pos;
297  for (pos = aux->begin(); pos != aux->end(); ++pos)
298  {
299  auxMap.insert(std::make_pair(pos->first,pos->second));
300  }
301  }
302 
303  // Return volume structure from child module
304  //
305  if (volname.empty())
306  {
307  return structure.GetVolume(structure.GetSetup("Default"));
308  }
309  else
310  {
311  return structure.GetVolume(structure.GenerateName(volname));
312  }
313 }
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
G4LogicalVolume * GetVolume(const G4String &) const
bool G4bool
Definition: G4Types.hh:79
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4GDMLAuxMapType * GetAuxMap() const
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
G4bool validate
Definition: G4GDMLRead.hh:145
std::map< G4LogicalVolume *, G4GDMLAuxListType > G4GDMLAuxMapType
G4String GetSetup(const G4String &)
void Read(const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
Definition: G4GDMLRead.cc:271
G4AssemblyVolume * G4GDMLReadStructure::GetAssembly ( const G4String ref) const

Definition at line 881 of file G4GDMLReadStructure.cc.

References assemblyMap.

Referenced by PhysvolRead().

882 {
883  G4GDMLAssemblyMapType::const_iterator pos = assemblyMap.find(ref);
884  if (pos != assemblyMap.end()) { return pos->second; }
885  return 0;
886 }
G4GDMLAssemblyMapType assemblyMap
const G4GDMLAuxMapType * G4GDMLReadStructure::GetAuxMap ( ) const

Definition at line 897 of file G4GDMLReadStructure.cc.

References auxMap.

Referenced by FileRead().

898 {
899  return &auxMap;
900 }
G4VPhysicalVolume * G4GDMLReadStructure::GetPhysvol ( const G4String ref) const

Definition at line 849 of file G4GDMLReadStructure.cc.

References FatalException, G4Exception(), G4PhysicalVolumeStore::GetInstance(), and G4PhysicalVolumeStore::GetVolume().

Referenced by BorderSurfaceRead().

850 {
851  G4VPhysicalVolume* physvolPtr =
853 
854  if (!physvolPtr)
855  {
856  G4String error_msg = "Referenced physvol '" + ref + "' was not found!";
857  G4Exception("G4GDMLReadStructure::GetPhysvol()", "ReadError",
858  FatalException, error_msg);
859  }
860 
861  return physvolPtr;
862 }
static G4PhysicalVolumeStore * GetInstance()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4VPhysicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
G4LogicalVolume * G4GDMLReadStructure::GetVolume ( const G4String ref) const
virtual

Implements G4GDMLRead.

Definition at line 865 of file G4GDMLReadStructure.cc.

References FatalException, G4Exception(), G4LogicalVolumeStore::GetInstance(), and G4LogicalVolumeStore::GetVolume().

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

866 {
867  G4LogicalVolume *volumePtr
869 
870  if (!volumePtr)
871  {
872  G4String error_msg = "Referenced volume '" + ref + "' was not found!";
873  G4Exception("G4GDMLReadStructure::GetVolume()", "ReadError",
874  FatalException, error_msg);
875  }
876 
877  return volumePtr;
878 }
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
static G4LogicalVolumeStore * GetInstance()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4GDMLAuxListType G4GDMLReadStructure::GetVolumeAuxiliaryInformation ( G4LogicalVolume logvol) const

Definition at line 889 of file G4GDMLReadStructure.cc.

References auxMap.

890 {
891  G4GDMLAuxMapType::const_iterator pos = auxMap.find(logvol);
892  if (pos != auxMap.end()) { return pos->second; }
893  else { return G4GDMLAuxListType(); }
894 }
std::vector< G4GDMLAuxPairType > G4GDMLAuxListType
G4VPhysicalVolume * G4GDMLReadStructure::GetWorldVolume ( const G4String setupName)

Definition at line 903 of file G4GDMLReadStructure.cc.

References G4LogicalVolume::GetName(), G4GDMLReadSetup::GetSetup(), GetVolume(), G4VisAttributes::Invisible, setuptoPV, G4LogicalVolume::SetVisAttributes(), and G4GDMLRead::Strip().

904 {
905  G4LogicalVolume* volume = GetVolume(Strip(GetSetup(setupName)));
907 
908  G4VPhysicalVolume* pvWorld = 0;
909 
910  if(setuptoPV[setupName])
911  {
912  pvWorld = setuptoPV[setupName];
913  }
914  else
915  {
916  pvWorld = new G4PVPlacement(0,G4ThreeVector(0,0,0),volume,
917  volume->GetName()+"_PV",0,0,0);
918  setuptoPV[setupName] = pvWorld;
919  }
920  return pvWorld;
921 }
G4String GetName() const
CLHEP::Hep3Vector G4ThreeVector
G4LogicalVolume * GetVolume(const G4String &) const
std::map< std::string, G4VPhysicalVolume * > setuptoPV
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:100
static const G4VisAttributes Invisible
void SetVisAttributes(const G4VisAttributes *pVA)
G4String GetSetup(const G4String &)
void G4GDMLReadStructure::PhysvolRead ( const xercesc::DOMElement * const  physvolElement,
G4AssemblyVolume assembly = 0 
)
protected

Definition at line 316 of file G4GDMLReadStructure.cc.

References G4AssemblyVolume::AddPlacedVolume(), G4GDMLRead::check, FatalException, FileRead(), G4String::first(), G4Exception(), G4GDMLRead::GenerateName(), G4GDMLRead::GeneratePhysvolName(), GetAssembly(), G4LogicalVolume::GetName(), G4GDMLReadDefine::GetPosition(), G4GDMLReadDefine::GetRotation(), G4GDMLReadDefine::GetRotationMatrix(), G4GDMLReadDefine::GetScale(), GetVolume(), G4ReflectionFactory::Instance(), G4AssemblyVolume::MakeImprint(), G4ReflectionFactory::Place(), pMotherLogical, position, G4GDMLReadDefine::RefRead(), G4GDMLRead::Transcode(), G4GDMLReadDefine::VectorRead(), CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

Referenced by AssemblyRead(), and Volume_contentRead().

318 {
319  G4String name;
320  G4LogicalVolume* logvol = 0;
321  G4AssemblyVolume* assembly = 0;
322  G4ThreeVector position(0.0,0.0,0.0);
323  G4ThreeVector rotation(0.0,0.0,0.0);
324  G4ThreeVector scale(1.0,1.0,1.0);
325 
326  const xercesc::DOMNamedNodeMap* const attributes
327  = physvolElement->getAttributes();
328  XMLSize_t attributeCount = attributes->getLength();
329 
330  for (XMLSize_t attribute_index=0;
331  attribute_index<attributeCount; attribute_index++)
332  {
333  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
334 
335  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
336  { continue; }
337 
338  const xercesc::DOMAttr* const attribute
339  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
340  if (!attribute)
341  {
342  G4Exception("G4GDMLReadStructure::PhysvolRead()",
343  "InvalidRead", FatalException, "No attribute found!");
344  return;
345  }
346  const G4String attName = Transcode(attribute->getName());
347  const G4String attValue = Transcode(attribute->getValue());
348 
349  if (attName=="name") { name = attValue; }
350  }
351 
352  for (xercesc::DOMNode* iter = physvolElement->getFirstChild();
353  iter != 0; iter = iter->getNextSibling())
354  {
355  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
356 
357  const xercesc::DOMElement* const child
358  = dynamic_cast<xercesc::DOMElement*>(iter);
359  if (!child)
360  {
361  G4Exception("G4GDMLReadStructure::PhysvolRead()",
362  "InvalidRead", FatalException, "No child found!");
363  return;
364  }
365  const G4String tag = Transcode(child->getTagName());
366 
367  if (tag=="volumeref")
368  {
369  const G4String& child_name = GenerateName(RefRead(child));
370  assembly = GetAssembly(child_name);
371  if (!assembly) { logvol = GetVolume(child_name); }
372  }
373  else if (tag=="file")
374  { logvol = FileRead(child); }
375  else if (tag=="position")
376  { VectorRead(child,position); }
377  else if (tag=="rotation")
378  { VectorRead(child,rotation); }
379  else if (tag=="scale")
380  { VectorRead(child,scale); }
381  else if (tag=="positionref")
382  { position = GetPosition(GenerateName(RefRead(child))); }
383  else if (tag=="rotationref")
384  { rotation = GetRotation(GenerateName(RefRead(child))); }
385  else if (tag=="scaleref")
386  { scale = GetScale(GenerateName(RefRead(child))); }
387  else
388  {
389  G4String error_msg = "Unknown tag in physvol: " + tag;
390  G4Exception("G4GDMLReadStructure::PhysvolRead()", "ReadError",
391  FatalException, error_msg);
392  return;
393  }
394  }
395 
396  G4Transform3D transform(GetRotationMatrix(rotation).inverse(),position);
397  transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
398 
399  if (pAssembly) // Fill assembly structure
400  {
401  if (!logvol) { return; }
402  pAssembly->AddPlacedVolume(logvol, transform);
403  }
404  else // Generate physical volume tree or do assembly imprint
405  {
406  if (assembly)
407  {
408  assembly->MakeImprint(pMotherLogical, transform, 0, check);
409  }
410  else
411  {
412  if (!logvol) { return; }
413  G4String pv_name = logvol->GetName() + "_PV";
415  ->Place(transform,pv_name,logvol,pMotherLogical,false,0,check);
416 
417  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
418  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
419  }
420  }
421 }
G4AssemblyVolume * GetAssembly(const G4String &) const
G4LogicalVolume * FileRead(const xercesc::DOMElement *const)
void MakeImprint(G4LogicalVolume *pMotherLV, G4ThreeVector &translationInMother, G4RotationMatrix *pRotationInMother, G4int copyNumBase=0, G4bool surfCheck=false)
G4String GetName() const
G4int first(char) const
Definition: xmlparse.cc:179
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
G4ThreeVector GetScale(const G4String &)
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:84
static G4ReflectionFactory * Instance()
G4String RefRead(const xercesc::DOMElement *const)
G4PhysicalVolumesPair Place(const G4Transform3D &transform3D, const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, G4bool isMany, G4int copyNo, G4bool surfCheck=false)
G4LogicalVolume * GetVolume(const G4String &) const
HepGeom::Scale3D G4Scale3D
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
int position
Definition: filter.cc:7
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
G4ThreeVector GetRotation(const G4String &)
G4RotationMatrix GetRotationMatrix(const G4ThreeVector &)
G4bool check
Definition: G4GDMLRead.hh:146
G4ThreeVector GetPosition(const G4String &)
G4double G4GDMLReadStructure::QuantityRead ( const xercesc::DOMElement *const  readElement)
protected

Definition at line 559 of file G4GDMLReadStructure.cc.

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

Referenced by ReplicaRead().

560 {
561  G4double value = 0.0;
562  G4double unit = 0.0;
563  const xercesc::DOMNamedNodeMap* const attributes
564  = readElement->getAttributes();
565  XMLSize_t attributeCount = attributes->getLength();
566 
567  for (XMLSize_t attribute_index=0;
568  attribute_index<attributeCount; attribute_index++)
569  {
570  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
571 
572  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
573  { continue; }
574  const xercesc::DOMAttr* const attribute
575  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
576  if (!attribute)
577  {
578  G4Exception("G4GDMLReadStructure::QuantityRead()",
579  "InvalidRead", FatalException, "No attribute found!");
580  return value;
581  }
582  const G4String attName = Transcode(attribute->getName());
583  const G4String attValue = Transcode(attribute->getValue());
584 
585  if (attName=="unit") { unit = eval.Evaluate(attValue); } else
586  if (attName=="value"){ value= eval.Evaluate(attValue); }
587  }
588 
589  return value*unit;
590 }
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:144
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const XML_Char int const XML_Char * value
double G4double
Definition: G4Types.hh:76
G4double Evaluate(const G4String &)
void G4GDMLReadStructure::ReplicaRead ( const xercesc::DOMElement *const  replicaElement,
G4LogicalVolume logvol,
G4int  number 
)
protected

Definition at line 460 of file G4GDMLReadStructure.cc.

References AxisRead(), FatalException, G4String::first(), G4Exception(), G4GDMLRead::GenerateName(), G4GDMLRead::GeneratePhysvolName(), G4LogicalVolume::GetName(), G4GDMLReadDefine::GetPosition(), G4GDMLReadDefine::GetRotation(), G4ReflectionFactory::Instance(), kUndefined, pMotherLogical, position, QuantityRead(), G4GDMLReadDefine::RefRead(), G4ReflectionFactory::Replicate(), G4GDMLRead::Transcode(), G4GDMLReadDefine::VectorRead(), and width.

Referenced by ReplicavolRead().

462 {
463  G4double width = 0.0;
464  G4double offset = 0.0;
465  G4ThreeVector position(0.0,0.0,0.0);
466  G4ThreeVector rotation(0.0,0.0,0.0);
467  EAxis axis = kUndefined;
468  G4String name;
469 
470  for (xercesc::DOMNode* iter = replicaElement->getFirstChild();
471  iter != 0; iter = iter->getNextSibling())
472  {
473  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
474 
475  const xercesc::DOMElement* const child
476  = dynamic_cast<xercesc::DOMElement*>(iter);
477  if (!child)
478  {
479  G4Exception("G4GDMLReadStructure::ReplicaRead()",
480  "InvalidRead", FatalException, "No child found!");
481  return;
482  }
483  const G4String tag = Transcode(child->getTagName());
484 
485  if (tag=="position")
486  { VectorRead(child,position); } else
487  if (tag=="rotation")
488  { VectorRead(child,rotation); } else
489  if (tag=="positionref")
490  { position = GetPosition(GenerateName(RefRead(child))); } else
491  if (tag=="rotationref")
492  { rotation = GetRotation(GenerateName(RefRead(child))); } else
493  if (tag=="direction")
494  { axis=AxisRead(child); } else
495  if (tag=="width")
496  { width=QuantityRead(child); } else
497  if (tag=="offset")
498  { offset=QuantityRead(child); }
499  else
500  {
501  G4String error_msg = "Unknown tag in ReplicaRead: " + tag;
502  G4Exception("G4GDMLReadStructure::ReplicaRead()", "ReadError",
503  FatalException, error_msg);
504  }
505  }
506 
507  G4String pv_name = logvol->GetName() + "_PV";
509  ->Replicate(pv_name,logvol,pMotherLogical,axis,number,width,offset);
510 
511  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
512  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
513 
514 }
G4String GetName() const
G4int first(char) const
Definition: xmlparse.cc:179
G4double QuantityRead(const xercesc::DOMElement *const readElement)
#define width
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
G4PhysicalVolumesPair Replicate(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofReplicas, G4double width, G4double offset=0)
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:84
static G4ReflectionFactory * Instance()
G4String RefRead(const xercesc::DOMElement *const)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
EAxis
Definition: geomdefs.hh:54
EAxis AxisRead(const xercesc::DOMElement *const axisElement)
int position
Definition: filter.cc:7
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
G4ThreeVector GetRotation(const G4String &)
double G4double
Definition: G4Types.hh:76
G4ThreeVector GetPosition(const G4String &)
void G4GDMLReadStructure::ReplicavolRead ( const xercesc::DOMElement * const  replicavolElement,
G4int  number 
)
protected

Definition at line 424 of file G4GDMLReadStructure.cc.

References FatalException, G4Exception(), G4GDMLRead::GenerateName(), GetVolume(), G4GDMLReadDefine::RefRead(), ReplicaRead(), and G4GDMLRead::Transcode().

Referenced by Volume_contentRead().

425 {
426  G4LogicalVolume* logvol = 0;
427  for (xercesc::DOMNode* iter = replicavolElement->getFirstChild();
428  iter != 0; iter = iter->getNextSibling())
429  {
430  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
431 
432  const xercesc::DOMElement* const child
433  = dynamic_cast<xercesc::DOMElement*>(iter);
434  if (!child)
435  {
436  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
437  "InvalidRead", FatalException, "No child found!");
438  return;
439  }
440  const G4String tag = Transcode(child->getTagName());
441 
442  if (tag=="volumeref")
443  {
444  logvol = GetVolume(GenerateName(RefRead(child)));
445  }
446  else if (tag=="replicate_along_axis")
447  {
448  if (logvol) { ReplicaRead(child,logvol,number); }
449  }
450  else
451  {
452  G4String error_msg = "Unknown tag in ReplicavolRead: " + tag;
453  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
454  "ReadError", FatalException, error_msg);
455  }
456  }
457 }
Definition: xmlparse.cc:179
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
G4String RefRead(const xercesc::DOMElement *const)
G4LogicalVolume * GetVolume(const G4String &) const
void ReplicaRead(const xercesc::DOMElement *const replicaElement, G4LogicalVolume *logvol, G4int number)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
void G4GDMLReadStructure::SkinSurfaceRead ( const xercesc::DOMElement * const  skinsurfaceElement)
protected

Definition at line 671 of file G4GDMLReadStructure.cc.

References FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4GDMLReadSolids::GetSurfaceProperty(), GetVolume(), G4GDMLReadDefine::RefRead(), G4GDMLRead::Strip(), and G4GDMLRead::Transcode().

Referenced by StructureRead().

672 {
673  G4String name;
674  G4LogicalVolume* logvol = 0;
675  G4SurfaceProperty* prop = 0;
676 
677  const xercesc::DOMNamedNodeMap* const attributes
678  = skinsurfaceElement->getAttributes();
679  XMLSize_t attributeCount = attributes->getLength();
680 
681  for (XMLSize_t attribute_index=0;
682  attribute_index<attributeCount; attribute_index++)
683  {
684  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
685 
686  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
687  { continue; }
688 
689  const xercesc::DOMAttr* const attribute
690  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
691  if (!attribute)
692  {
693  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()",
694  "InvalidRead", FatalException, "No attribute found!");
695  return;
696  }
697  const G4String attName = Transcode(attribute->getName());
698  const G4String attValue = Transcode(attribute->getValue());
699 
700  if (attName=="name")
701  { name = GenerateName(attValue); } else
702  if (attName=="surfaceproperty")
703  { prop = GetSurfaceProperty(GenerateName(attValue)); }
704  }
705 
706  for (xercesc::DOMNode* iter = skinsurfaceElement->getFirstChild();
707  iter != 0; iter = iter->getNextSibling())
708  {
709  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
710 
711  const xercesc::DOMElement* const child
712  = dynamic_cast<xercesc::DOMElement*>(iter);
713  if (!child)
714  {
715  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()",
716  "InvalidRead", FatalException, "No child found!");
717  return;
718  }
719  const G4String tag = Transcode(child->getTagName());
720 
721  if (tag=="volumeref")
722  {
723  logvol = GetVolume(GenerateName(RefRead(child)));
724  }
725  else
726  {
727  G4String error_msg = "Unknown tag in skinsurface: " + tag;
728  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "ReadError",
729  FatalException, error_msg);
730  }
731  }
732 
733  new G4LogicalSkinSurface(Strip(name),logvol,prop);
734 }
Definition: xmlparse.cc:179
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
G4SurfaceProperty * GetSurfaceProperty(const G4String &) const
G4String RefRead(const xercesc::DOMElement *const)
G4LogicalVolume * GetVolume(const G4String &) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
G4String Strip(const G4String &) const
Definition: G4GDMLRead.cc:100
void G4GDMLReadStructure::StructureRead ( const xercesc::DOMElement * const  structureElement)
virtual

Implements G4GDMLRead.

Definition at line 815 of file G4GDMLReadStructure.cc.

References AssemblyRead(), BorderSurfaceRead(), FatalException, G4cout, G4endl, G4Exception(), G4GDMLRead::LoopRead(), SkinSurfaceRead(), G4GDMLRead::StructureRead(), G4GDMLRead::Transcode(), and VolumeRead().

816 {
817  G4cout << "G4GDML: Reading structure..." << G4endl;
818 
819  for (xercesc::DOMNode* iter = structureElement->getFirstChild();
820  iter != 0; iter = iter->getNextSibling())
821  {
822  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
823 
824  const xercesc::DOMElement* const child
825  = dynamic_cast<xercesc::DOMElement*>(iter);
826  if (!child)
827  {
828  G4Exception("G4GDMLReadStructure::StructureRead()",
829  "InvalidRead", FatalException, "No child found!");
830  return;
831  }
832  const G4String tag = Transcode(child->getTagName());
833 
834  if (tag=="bordersurface") { BorderSurfaceRead(child); } else
835  if (tag=="skinsurface") { SkinSurfaceRead(child); } else
836  if (tag=="volume") { VolumeRead(child); } else
837  if (tag=="assembly") { AssemblyRead(child); } else
838  if (tag=="loop") { LoopRead(child,&G4GDMLRead::StructureRead); }
839  else
840  {
841  G4String error_msg = "Unknown tag in structure: " + tag;
842  G4Exception("G4GDMLReadStructure::StructureRead()",
843  "ReadError", FatalException, error_msg);
844  }
845  }
846 }
Definition: xmlparse.cc:179
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
virtual void StructureRead(const xercesc::DOMElement *const)=0
void LoopRead(const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
Definition: G4GDMLRead.cc:179
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void BorderSurfaceRead(const xercesc::DOMElement *const)
void SkinSurfaceRead(const xercesc::DOMElement *const)
#define G4endl
Definition: G4ios.hh:61
void AssemblyRead(const xercesc::DOMElement *const)
virtual void VolumeRead(const xercesc::DOMElement *const)
void G4GDMLReadStructure::Volume_contentRead ( const xercesc::DOMElement * const  volumeElement)
virtual

Implements G4GDMLRead.

Definition at line 737 of file G4GDMLReadStructure.cc.

References DivisionvolRead(), G4GDMLRead::eval, G4GDMLEvaluator::EvaluateInteger(), FatalException, G4cout, G4endl, G4Exception(), G4GDMLRead::LoopRead(), G4GDMLReadParamvol::ParamvolRead(), PhysvolRead(), pMotherLogical, ReplicavolRead(), G4GDMLRead::Transcode(), and G4GDMLRead::Volume_contentRead().

Referenced by G03ColorReader::VolumeRead(), and VolumeRead().

738 {
739  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
740  iter != 0; iter = iter->getNextSibling())
741  {
742  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
743 
744  const xercesc::DOMElement* const child
745  = dynamic_cast<xercesc::DOMElement*>(iter);
746  if (!child)
747  {
748  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
749  "InvalidRead", FatalException, "No child found!");
750  return;
751  }
752  const G4String tag = Transcode(child->getTagName());
753 
754  if ((tag=="auxiliary") || (tag=="materialref") || (tag=="solidref"))
755  {
756  // These are already processed in VolumeRead()
757  }
758  else if (tag=="paramvol")
759  {
761  }
762  else if (tag=="physvol")
763  {
764  PhysvolRead(child);
765  }
766  else if (tag=="replicavol")
767  {
768  G4int number = 1;
769  const xercesc::DOMNamedNodeMap* const attributes
770  = child->getAttributes();
771  XMLSize_t attributeCount = attributes->getLength();
772  for (XMLSize_t attribute_index=0;
773  attribute_index<attributeCount; attribute_index++)
774  {
775  xercesc::DOMNode* attribute_node
776  = attributes->item(attribute_index);
777  if (attribute_node->getNodeType()!=xercesc::DOMNode::ATTRIBUTE_NODE)
778  {
779  continue;
780  }
781  const xercesc::DOMAttr* const attribute
782  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
783  if (!attribute)
784  {
785  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
786  "InvalidRead", FatalException, "No attribute found!");
787  return;
788  }
789  const G4String attName = Transcode(attribute->getName());
790  const G4String attValue = Transcode(attribute->getValue());
791  if (attName=="number")
792  {
793  number = eval.EvaluateInteger(attValue);
794  }
795  }
796  ReplicavolRead(child,number);
797  }
798  else if (tag=="divisionvol")
799  {
800  DivisionvolRead(child);
801  }
802  else if (tag=="loop")
803  {
805  }
806  else
807  {
808  G4cout << "Treating unknown GDML tag in volume '" << tag
809  << "' as GDML extension..." << G4endl;
810  }
811  }
812 }
G4int EvaluateInteger(const G4String &)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:144
virtual void ParamvolRead(const xercesc::DOMElement *const, G4LogicalVolume *)
virtual void Volume_contentRead(const xercesc::DOMElement *const)=0
Definition: xmlparse.cc:179
void DivisionvolRead(const xercesc::DOMElement *const)
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
void ReplicavolRead(const xercesc::DOMElement *const, G4int number)
int G4int
Definition: G4Types.hh:78
void LoopRead(const xercesc::DOMElement *const, void(G4GDMLRead::*)(const xercesc::DOMElement *const))
Definition: G4GDMLRead.cc:179
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
#define G4endl
Definition: G4ios.hh:61
void PhysvolRead(const xercesc::DOMElement *const, G4AssemblyVolume *assembly=0)
void G4GDMLReadStructure::VolumeRead ( const xercesc::DOMElement * const  volumeElement)
virtual

Reimplemented in G03ColorReader.

Definition at line 593 of file G4GDMLReadStructure.cc.

References AuxiliaryRead(), auxMap, FatalException, G4Exception(), G4GDMLRead::GenerateName(), G4GDMLReadMaterials::GetMaterial(), G4GDMLReadSolids::GetSolid(), pMotherLogical, G4GDMLReadDefine::RefRead(), G4GDMLRead::Transcode(), and Volume_contentRead().

Referenced by StructureRead().

594 {
595  G4VSolid* solidPtr = 0;
596  G4Material* materialPtr = 0;
597  G4GDMLAuxListType auxList;
598 
599  XMLCh *name_attr = xercesc::XMLString::transcode("name");
600  const G4String name = Transcode(volumeElement->getAttribute(name_attr));
601  xercesc::XMLString::release(&name_attr);
602 
603  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
604  iter != 0; iter = iter->getNextSibling())
605  {
606  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
607 
608  const xercesc::DOMElement* const child
609  = dynamic_cast<xercesc::DOMElement*>(iter);
610  if (!child)
611  {
612  G4Exception("G4GDMLReadStructure::VolumeRead()",
613  "InvalidRead", FatalException, "No child found!");
614  return;
615  }
616  const G4String tag = Transcode(child->getTagName());
617 
618  if (tag=="auxiliary")
619  { auxList.push_back(AuxiliaryRead(child)); } else
620  if (tag=="materialref")
621  { materialPtr = GetMaterial(GenerateName(RefRead(child),true)); } else
622  if (tag=="solidref")
623  { solidPtr = GetSolid(GenerateName(RefRead(child))); }
624  }
625 
626  pMotherLogical = new G4LogicalVolume(solidPtr,materialPtr,
627  GenerateName(name),0,0,0);
628 
629  if (!auxList.empty()) { auxMap[pMotherLogical] = auxList; }
630 
631  Volume_contentRead(volumeElement);
632 }
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
Definition: xmlparse.cc:179
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
const XML_Char * name
G4GDMLAuxPairType AuxiliaryRead(const xercesc::DOMElement *const)
G4VSolid * GetSolid(const G4String &) const
G4String RefRead(const xercesc::DOMElement *const)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
std::vector< G4GDMLAuxPairType > G4GDMLAuxListType
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
virtual void Volume_contentRead(const xercesc::DOMElement *const)

Field Documentation

G4GDMLAssemblyMapType G4GDMLReadStructure::assemblyMap
protected

Definition at line 101 of file G4GDMLReadStructure.hh.

Referenced by AssemblyRead(), and GetAssembly().

G4GDMLAuxMapType G4GDMLReadStructure::auxMap
protected
G4LogicalVolume* G4GDMLReadStructure::pMotherLogical
protected
std::map<std::string, G4VPhysicalVolume*> G4GDMLReadStructure::setuptoPV
protected

Definition at line 103 of file G4GDMLReadStructure.hh.

Referenced by Clear(), and GetWorldVolume().


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