Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4GDMLReadStructure Class Reference

#include <G4GDMLReadStructure.hh>

Inheritance diagram for G4GDMLReadStructure:
Collaboration diagram for G4GDMLReadStructure:

Public Member Functions

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

Protected Member Functions

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

Protected Attributes

G4GDMLAuxMapType auxMap
 
G4GDMLAssemblyMapType assemblyMap
 
G4LogicalVolumepMotherLogical
 
std::map< std::string,
G4VPhysicalVolume * > 
setuptoPV
 
G4bool strip
 
- 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
 
G4bool dostrip
 

Detailed Description

Definition at line 55 of file G4GDMLReadStructure.hh.

Constructor & Destructor Documentation

G4GDMLReadStructure::G4GDMLReadStructure ( )

Definition at line 49 of file G4GDMLReadStructure.cc.

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

Definition at line 54 of file G4GDMLReadStructure.cc.

55 {
56 }

Member Function Documentation

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

Definition at line 615 of file G4GDMLReadStructure.cc.

616 {
617  XMLCh *name_attr = xercesc::XMLString::transcode("name");
618  const G4String name = Transcode(assemblyElement->getAttribute(name_attr));
619  xercesc::XMLString::release(&name_attr);
620 
621  G4AssemblyVolume* pAssembly = new G4AssemblyVolume();
622  assemblyMap.insert(std::make_pair(GenerateName(name), pAssembly));
623 
624  for (xercesc::DOMNode* iter = assemblyElement->getFirstChild();
625  iter != 0; iter = iter->getNextSibling())
626  {
627  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
628  const xercesc::DOMElement* const child
629  = dynamic_cast<xercesc::DOMElement*>(iter);
630  if (!child)
631  {
632  G4Exception("G4GDMLReadStructure::AssemblyRead()",
633  "InvalidRead", FatalException, "No child found!");
634  return;
635  }
636  const G4String tag = Transcode(child->getTagName());
637 
638  if (tag=="physvol")
639  {
640  PhysvolRead(child, pAssembly);
641  }
642  else
643  {
644  G4cout << "Unsupported GDML tag '" << tag
645  << "' for Geant4 assembly structure !" << G4endl;
646  }
647  }
648 }
const XML_Char * name
Definition: expat.h:151
Definition: xmlparse.cc:187
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
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)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 493 of file G4GDMLReadStructure.cc.

494 {
495  EAxis axis = kUndefined;
496 
497  const xercesc::DOMNamedNodeMap* const attributes
498  = axisElement->getAttributes();
499  XMLSize_t attributeCount = attributes->getLength();
500 
501  for (XMLSize_t attribute_index=0;
502  attribute_index<attributeCount; attribute_index++)
503  {
504  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
505 
506  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
507  { continue; }
508 
509  const xercesc::DOMAttr* const attribute
510  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
511  if (!attribute)
512  {
513  G4Exception("G4GDMLReadStructure::AxisRead()",
514  "InvalidRead", FatalException, "No attribute found!");
515  return axis;
516  }
517  const G4String attName = Transcode(attribute->getName());
518  const G4String attValue = Transcode(attribute->getValue());
519  if (attName=="x")
520  { if( eval.Evaluate(attValue)==1.) {axis=kXAxis;} }
521  else if (attName=="y")
522  { if( eval.Evaluate(attValue)==1.) {axis=kYAxis;} }
523  else if (attName=="z")
524  { if( eval.Evaluate(attValue)==1.) {axis=kZAxis;} }
525  else if (attName=="rho")
526  { if( eval.Evaluate(attValue)==1.) {axis=kRho;} }
527  else if (attName=="phi")
528  { if( eval.Evaluate(attValue)==1.) {axis=kPhi;} }
529  }
530 
531  return axis;
532 }
Definition: geomdefs.hh:54
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
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 &)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 60 of file G4GDMLReadStructure.cc.

61 {
62  G4String name;
63  G4VPhysicalVolume* pv1 = 0;
64  G4VPhysicalVolume* pv2 = 0;
65  G4SurfaceProperty* prop = 0;
66  G4int index = 0;
67 
68  const xercesc::DOMNamedNodeMap* const attributes
69  = bordersurfaceElement->getAttributes();
70  XMLSize_t attributeCount = attributes->getLength();
71 
72  for (XMLSize_t attribute_index=0;
73  attribute_index<attributeCount; attribute_index++)
74  {
75  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
76 
77  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
78  { continue; }
79 
80  const xercesc::DOMAttr* const attribute
81  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
82  if (!attribute)
83  {
84  G4Exception("G4GDMLReadStructure::BorderSurfaceRead()",
85  "InvalidRead", FatalException, "No attribute found!");
86  return;
87  }
88  const G4String attName = Transcode(attribute->getName());
89  const G4String attValue = Transcode(attribute->getValue());
90 
91  if (attName=="name")
92  { name = GenerateName(attValue); } else
93  if (attName=="surfaceproperty")
94  { prop = GetSurfaceProperty(GenerateName(attValue)); }
95  }
96 
97  for (xercesc::DOMNode* iter = bordersurfaceElement->getFirstChild();
98  iter != 0; iter = iter->getNextSibling())
99  {
100  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
101 
102  const xercesc::DOMElement* const child
103  = dynamic_cast<xercesc::DOMElement*>(iter);
104  if (!child)
105  {
106  G4Exception("G4GDMLReadStructure::BorderSurfaceRead()",
107  "InvalidRead", FatalException, "No child found!");
108  return;
109  }
110  const G4String tag = Transcode(child->getTagName());
111 
112  if (tag != "physvolref") { continue; }
113 
114  if (index==0)
115  { pv1 = GetPhysvol(GenerateName(RefRead(child))); index++; } else
116  if (index==1)
117  { pv2 = GetPhysvol(GenerateName(RefRead(child))); index++; } else
118  break;
119  }
120 
121  new G4LogicalBorderSurface(Strip(name),pv1,pv2,prop);
122 }
const XML_Char * name
Definition: expat.h:151
Definition: xmlparse.cc:187
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
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:97

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLReadStructure::Clear ( )

Definition at line 900 of file G4GDMLReadStructure.cc.

901 {
902  eval.Clear();
903  setuptoPV.clear();
904 }
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
std::map< std::string, G4VPhysicalVolume * > setuptoPV

Here is the call graph for this function:

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

Definition at line 125 of file G4GDMLReadStructure.cc.

126 {
127  G4String name;
128  G4double unit = 1.0;
129  G4double width = 0.0;
130  G4double offset = 0.0;
131  G4int number = 0;
132  EAxis axis = kUndefined;
133  G4LogicalVolume* logvol = 0;
134 
135  const xercesc::DOMNamedNodeMap* const attributes
136  = divisionvolElement->getAttributes();
137  XMLSize_t attributeCount = attributes->getLength();
138  G4String unitname;
139 
140  for (XMLSize_t attribute_index=0;
141  attribute_index<attributeCount; attribute_index++)
142  {
143  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
144 
145  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
146  { continue; }
147 
148  const xercesc::DOMAttr* const attribute
149  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
150  if (!attribute)
151  {
152  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
153  "InvalidRead", FatalException, "No attribute found!");
154  return;
155  }
156  const G4String attName = Transcode(attribute->getName());
157  const G4String attValue = Transcode(attribute->getValue());
158 
159  if (attName=="name") { name = attValue; } else
160  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue);
161  unitname = G4UnitDefinition::GetCategory(attValue);
162  } else
163  if (attName=="width") { width = eval.Evaluate(attValue); } else
164  if (attName=="offset") { offset = eval.Evaluate(attValue); } else
165  if (attName=="number") { number = eval.EvaluateInteger(attValue); } else
166  if (attName=="axis")
167  {
168  if (attValue=="kXAxis") { axis = kXAxis; } else
169  if (attValue=="kYAxis") { axis = kYAxis; } else
170  if (attValue=="kZAxis") { axis = kZAxis; } else
171  if (attValue=="kRho") { axis = kRho; } else
172  if (attValue=="kPhi") { axis = kPhi; }
173  }
174  }
175 
176  if ( ((axis == kXAxis || axis == kYAxis || axis == kZAxis) && unitname!="Length") ||
177  ((axis == kRho || axis == kPhi) && unitname!="Angle")) {
178  G4Exception("G4GDMLReadStructure::DivisionvolRead()", "InvalidRead",
179  FatalException, "Invalid unit!"); }
180 
181  width *= unit;
182  offset *= unit;
183 
184  for (xercesc::DOMNode* iter = divisionvolElement->getFirstChild();
185  iter != 0;iter = iter->getNextSibling())
186  {
187  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
188 
189  const xercesc::DOMElement* const child
190  = dynamic_cast<xercesc::DOMElement*>(iter);
191  if (!child)
192  {
193  G4Exception("G4GDMLReadStructure::DivisionvolRead()",
194  "InvalidRead", FatalException, "No child found!");
195  return;
196  }
197  const G4String tag = Transcode(child->getTagName());
198 
199  if (tag=="volumeref") { logvol = GetVolume(GenerateName(RefRead(child))); }
200  }
201 
202  if (!logvol) { return; }
203 
206 
207  G4String pv_name = logvol->GetName() + "_div";
208  if ((number != 0) && (width == 0.0))
209  {
211  ->Divide(pv_name,logvol,pMotherLogical,axis,number,offset);
212  }
213  else if ((number == 0) && (width != 0.0))
214  {
216  ->Divide(pv_name,logvol,pMotherLogical,axis,width,offset);
217  }
218  else
219  {
221  ->Divide(pv_name,logvol,pMotherLogical,axis,number,width,offset);
222  }
223 
224  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
225  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
226 }
G4PhysicalVolumesPair Divide(const G4String &name, G4LogicalVolume *LV, G4LogicalVolume *motherLV, EAxis axis, G4int nofDivisions, G4double width, G4double offset)
Definition: geomdefs.hh:54
const XML_Char * name
Definition: expat.h:151
G4int EvaluateInteger(const G4String &)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
G4int first(char) const
Definition: xmlparse.cc:187
#define width
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:81
int G4int
Definition: G4Types.hh:78
static G4ReflectionFactory * Instance()
G4String RefRead(const xercesc::DOMElement *const)
static G4PVDivisionFactory * GetInstance()
static G4double GetValueOf(const G4String &)
G4LogicalVolume * GetVolume(const G4String &) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4LogicalVolume * pMotherLogical
static G4String GetCategory(const G4String &)
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
EAxis
Definition: geomdefs.hh:54
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
Definition: geomdefs.hh:54
G4double Evaluate(const G4String &)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 229 of file G4GDMLReadStructure.cc.

230 {
231  G4String name;
232  G4String volname;
233 
234  const xercesc::DOMNamedNodeMap* const attributes
235  = fileElement->getAttributes();
236  XMLSize_t attributeCount = attributes->getLength();
237 
238  for (XMLSize_t attribute_index=0;
239  attribute_index<attributeCount; attribute_index++)
240  {
241  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
242 
243  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
244  { continue; }
245 
246  const xercesc::DOMAttr* const attribute
247  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
248  if (!attribute)
249  {
250  G4Exception("G4GDMLReadStructure::FileRead()",
251  "InvalidRead", FatalException, "No attribute found!");
252  return 0;
253  }
254  const G4String attName = Transcode(attribute->getName());
255  const G4String attValue = Transcode(attribute->getValue());
256 
257  if (attName=="name") { name = attValue; } else
258  if (attName=="volname") { volname = attValue; }
259  }
260 
261  const G4bool isModule = true;
262  G4GDMLReadStructure structure;
263  structure.Read(name,validate,isModule);
264 
265  // Register existing auxiliar information defined in child module
266  //
267  const G4GDMLAuxMapType* aux = structure.GetAuxMap();
268  if (!aux->empty())
269  {
270  G4GDMLAuxMapType::const_iterator pos;
271  for (pos = aux->begin(); pos != aux->end(); ++pos)
272  {
273  auxMap.insert(std::make_pair(pos->first,pos->second));
274  }
275  }
276 
277  // Return volume structure from child module
278  //
279  if (volname.empty())
280  {
281  return structure.GetVolume(structure.GetSetup("Default"));
282  }
283  else
284  {
285  return structure.GetVolume(structure.GenerateName(volname));
286  }
287 }
const XML_Char * name
Definition: expat.h:151
const G4GDMLAuxMapType * GetAuxMap() const
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
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
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
G4GDMLAuxMapType auxMap
G4bool validate
Definition: G4GDMLRead.hh:156
std::map< G4LogicalVolume *, G4GDMLAuxListType > G4GDMLAuxMapType
static const G4double pos
G4String GetSetup(const G4String &)
void Read(const G4String &, G4bool validation, G4bool isModule, G4bool strip=true)
Definition: G4GDMLRead.cc:361

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 861 of file G4GDMLReadStructure.cc.

862 {
863  G4GDMLAssemblyMapType::const_iterator pos = assemblyMap.find(ref);
864  if (pos != assemblyMap.end()) { return pos->second; }
865  return 0;
866 }
G4GDMLAssemblyMapType assemblyMap
static const G4double pos

Here is the caller graph for this function:

const G4GDMLAuxMapType* G4GDMLReadStructure::GetAuxMap ( ) const
inline

Definition at line 68 of file G4GDMLReadStructure.hh.

68 {return &auxMap;}
G4GDMLAuxMapType auxMap

Here is the caller graph for this function:

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

Definition at line 829 of file G4GDMLReadStructure.cc.

830 {
831  G4VPhysicalVolume* physvolPtr =
833 
834  if (!physvolPtr)
835  {
836  G4String error_msg = "Referenced physvol '" + ref + "' was not found!";
837  G4Exception("G4GDMLReadStructure::GetPhysvol()", "ReadError",
838  FatalException, error_msg);
839  }
840 
841  return physvolPtr;
842 }
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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4GDMLRead.

Definition at line 845 of file G4GDMLReadStructure.cc.

846 {
847  G4LogicalVolume *volumePtr
849 
850  if (!volumePtr)
851  {
852  G4String error_msg = "Referenced volume '" + ref + "' was not found!";
853  G4Exception("G4GDMLReadStructure::GetVolume()", "ReadError",
854  FatalException, error_msg);
855  }
856 
857  return volumePtr;
858 }
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

Here is the call graph for this function:

Here is the caller graph for this function:

G4GDMLAuxListType G4GDMLReadStructure::GetVolumeAuxiliaryInformation ( G4LogicalVolume logvol) const

Definition at line 869 of file G4GDMLReadStructure.cc.

870 {
871  G4GDMLAuxMapType::const_iterator pos = auxMap.find(logvol);
872  if (pos != auxMap.end()) { return pos->second; }
873  else { return G4GDMLAuxListType(); }
874 }
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
G4GDMLAuxMapType auxMap
static const G4double pos
G4VPhysicalVolume * G4GDMLReadStructure::GetWorldVolume ( const G4String setupName)

Definition at line 877 of file G4GDMLReadStructure.cc.

878 {
879  G4String sname = GetSetup(setupName);
880  if (sname == "") { return 0; }
881 
882  G4LogicalVolume* volume = GetVolume(GenerateName(sname, dostrip));
884 
885  G4VPhysicalVolume* pvWorld = 0;
886 
887  if(setuptoPV[setupName])
888  {
889  pvWorld = setuptoPV[setupName];
890  }
891  else
892  {
893  pvWorld = new G4PVPlacement(0,G4ThreeVector(0,0,0),volume,
894  volume->GetName()+"_PV",0,0,0);
895  setuptoPV[setupName] = pvWorld;
896  }
897  return pvWorld;
898 }
CLHEP::Hep3Vector G4ThreeVector
G4bool dostrip
Definition: G4GDMLRead.hh:158
G4LogicalVolume * GetVolume(const G4String &) const
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
std::map< std::string, G4VPhysicalVolume * > setuptoPV
const G4String & GetName() const
static const G4VisAttributes & GetInvisible()
void SetVisAttributes(const G4VisAttributes *pVA)
G4String GetSetup(const G4String &)

Here is the call graph for this function:

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

Definition at line 290 of file G4GDMLReadStructure.cc.

292 {
293  G4String name;
294  G4LogicalVolume* logvol = 0;
295  G4AssemblyVolume* assembly = 0;
296  G4ThreeVector position(0.0,0.0,0.0);
297  G4ThreeVector rotation(0.0,0.0,0.0);
298  G4ThreeVector scale(1.0,1.0,1.0);
299  G4int copynumber = 0;
300 
301  const xercesc::DOMNamedNodeMap* const attributes
302  = physvolElement->getAttributes();
303  XMLSize_t attributeCount = attributes->getLength();
304 
305  for (XMLSize_t attribute_index=0;
306  attribute_index<attributeCount; attribute_index++)
307  {
308  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
309 
310  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
311  { continue; }
312 
313  const xercesc::DOMAttr* const attribute
314  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
315  if (!attribute)
316  {
317  G4Exception("G4GDMLReadStructure::PhysvolRead()",
318  "InvalidRead", FatalException, "No attribute found!");
319  return;
320  }
321  const G4String attName = Transcode(attribute->getName());
322  const G4String attValue = Transcode(attribute->getValue());
323 
324  if (attName=="name") { name = attValue; }
325  if (attName=="copynumber") { copynumber = eval.EvaluateInteger(attValue); }
326  }
327 
328  for (xercesc::DOMNode* iter = physvolElement->getFirstChild();
329  iter != 0; iter = iter->getNextSibling())
330  {
331  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
332 
333  const xercesc::DOMElement* const child
334  = dynamic_cast<xercesc::DOMElement*>(iter);
335  if (!child)
336  {
337  G4Exception("G4GDMLReadStructure::PhysvolRead()",
338  "InvalidRead", FatalException, "No child found!");
339  return;
340  }
341  const G4String tag = Transcode(child->getTagName());
342 
343  if (tag=="volumeref")
344  {
345  const G4String& child_name = GenerateName(RefRead(child));
346  assembly = GetAssembly(child_name);
347  if (!assembly) { logvol = GetVolume(child_name); }
348  }
349  else if (tag=="file")
350  { logvol = FileRead(child); }
351  else if (tag=="position")
352  { VectorRead(child,position); }
353  else if (tag=="rotation")
354  { VectorRead(child,rotation); }
355  else if (tag=="scale")
356  { VectorRead(child,scale); }
357  else if (tag=="positionref")
358  { position = GetPosition(GenerateName(RefRead(child))); }
359  else if (tag=="rotationref")
360  { rotation = GetRotation(GenerateName(RefRead(child))); }
361  else if (tag=="scaleref")
362  { scale = GetScale(GenerateName(RefRead(child))); }
363  else
364  {
365  G4String error_msg = "Unknown tag in physvol: " + tag;
366  G4Exception("G4GDMLReadStructure::PhysvolRead()", "ReadError",
367  FatalException, error_msg);
368  return;
369  }
370  }
371 
372  G4Transform3D transform(GetRotationMatrix(rotation).inverse(),position);
373  transform = transform*G4Scale3D(scale.x(),scale.y(),scale.z());
374 
375  if (pAssembly) // Fill assembly structure
376  {
377  if (!logvol) { return; }
378  pAssembly->AddPlacedVolume(logvol, transform);
379  }
380  else // Generate physical volume tree or do assembly imprint
381  {
382  if (assembly)
383  {
384  assembly->MakeImprint(pMotherLogical, transform, 0, check);
385  }
386  else
387  {
388  if (!logvol) { return; }
389  G4String pv_name = logvol->GetName() + "_PV";
391  ->Place(transform,pv_name,logvol,pMotherLogical,false,copynumber,check);
392 
393  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
394  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
395  }
396  }
397 }
const XML_Char * name
Definition: expat.h:151
G4AssemblyVolume * GetAssembly(const G4String &) const
G4LogicalVolume * FileRead(const xercesc::DOMElement *const)
G4int EvaluateInteger(const G4String &)
void MakeImprint(G4LogicalVolume *pMotherLV, G4ThreeVector &translationInMother, G4RotationMatrix *pRotationInMother, G4int copyNumBase=0, G4bool surfCheck=false)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
G4int first(char) const
Definition: xmlparse.cc:187
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
G4ThreeVector GetScale(const G4String &)
void GeneratePhysvolName(const G4String &, G4VPhysicalVolume *)
Definition: G4GDMLRead.cc:81
int G4int
Definition: G4Types.hh:78
static G4ReflectionFactory * Instance()
G4String RefRead(const xercesc::DOMElement *const)
#define position
Definition: xmlparse.cc:622
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
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
G4ThreeVector GetRotation(const G4String &)
G4RotationMatrix GetRotationMatrix(const G4ThreeVector &)
const G4String & GetName() const
G4bool check
Definition: G4GDMLRead.hh:157
G4ThreeVector GetPosition(const G4String &)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 535 of file G4GDMLReadStructure.cc.

536 {
537  G4double value = 0.0;
538  G4double unit = 0.0;
539  const xercesc::DOMNamedNodeMap* const attributes
540  = readElement->getAttributes();
541  XMLSize_t attributeCount = attributes->getLength();
542 
543  for (XMLSize_t attribute_index=0;
544  attribute_index<attributeCount; attribute_index++)
545  {
546  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
547 
548  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
549  { continue; }
550  const xercesc::DOMAttr* const attribute
551  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
552  if (!attribute)
553  {
554  G4Exception("G4GDMLReadStructure::QuantityRead()",
555  "InvalidRead", FatalException, "No attribute found!");
556  return value;
557  }
558  const G4String attName = Transcode(attribute->getName());
559  const G4String attValue = Transcode(attribute->getValue());
560 
561  if (attName=="unit") { unit = G4UnitDefinition::GetValueOf(attValue);
562  if (G4UnitDefinition::GetCategory(attValue)!="Length" && G4UnitDefinition::GetCategory(attValue)!="Angle") {
563  G4Exception("G4GDMLReadStructure::QuantityRead()", "InvalidRead",
564  FatalException, "Invalid unit for lenght or angle (width, offset)!"); }
565  } else
566  if (attName=="value"){ value= eval.Evaluate(attValue); }
567  }
568 
569  return value*unit;
570 }
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
static G4double GetValueOf(const G4String &)
const XML_Char int const XML_Char * value
Definition: expat.h:331
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4String GetCategory(const G4String &)
double G4double
Definition: G4Types.hh:76
G4double Evaluate(const G4String &)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 436 of file G4GDMLReadStructure.cc.

438 {
439  G4double width = 0.0;
440  G4double offset = 0.0;
441  G4ThreeVector position(0.0,0.0,0.0);
442  G4ThreeVector rotation(0.0,0.0,0.0);
443  EAxis axis = kUndefined;
444  G4String name;
445 
446  for (xercesc::DOMNode* iter = replicaElement->getFirstChild();
447  iter != 0; iter = iter->getNextSibling())
448  {
449  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
450 
451  const xercesc::DOMElement* const child
452  = dynamic_cast<xercesc::DOMElement*>(iter);
453  if (!child)
454  {
455  G4Exception("G4GDMLReadStructure::ReplicaRead()",
456  "InvalidRead", FatalException, "No child found!");
457  return;
458  }
459  const G4String tag = Transcode(child->getTagName());
460 
461  if (tag=="position")
462  { VectorRead(child,position); } else
463  if (tag=="rotation")
464  { VectorRead(child,rotation); } else
465  if (tag=="positionref")
466  { position = GetPosition(GenerateName(RefRead(child))); } else
467  if (tag=="rotationref")
468  { rotation = GetRotation(GenerateName(RefRead(child))); } else
469  if (tag=="direction")
470  { axis=AxisRead(child); } else
471  if (tag=="width")
472  { width=QuantityRead(child); } else
473  if (tag=="offset")
474  { offset=QuantityRead(child); }
475  else
476  {
477  G4String error_msg = "Unknown tag in ReplicaRead: " + tag;
478  G4Exception("G4GDMLReadStructure::ReplicaRead()", "ReadError",
479  FatalException, error_msg);
480  }
481  }
482 
483  G4String pv_name = logvol->GetName() + "_PV";
485  ->Replicate(pv_name,logvol,pMotherLogical,axis,number,width,offset);
486 
487  if (pair.first != 0) { GeneratePhysvolName(name,pair.first); }
488  if (pair.second != 0) { GeneratePhysvolName(name,pair.second); }
489 
490 }
const XML_Char * name
Definition: expat.h:151
G4int first(char) const
Definition: xmlparse.cc:187
G4double QuantityRead(const xercesc::DOMElement *const readElement)
#define width
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
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:81
static G4ReflectionFactory * Instance()
G4String RefRead(const xercesc::DOMElement *const)
#define position
Definition: xmlparse.cc:622
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)
void VectorRead(const xercesc::DOMElement *const, G4ThreeVector &)
std::pair< G4VPhysicalVolume *, G4VPhysicalVolume * > G4PhysicalVolumesPair
G4ThreeVector GetRotation(const G4String &)
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
G4ThreeVector GetPosition(const G4String &)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 400 of file G4GDMLReadStructure.cc.

401 {
402  G4LogicalVolume* logvol = 0;
403  for (xercesc::DOMNode* iter = replicavolElement->getFirstChild();
404  iter != 0; iter = iter->getNextSibling())
405  {
406  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
407 
408  const xercesc::DOMElement* const child
409  = dynamic_cast<xercesc::DOMElement*>(iter);
410  if (!child)
411  {
412  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
413  "InvalidRead", FatalException, "No child found!");
414  return;
415  }
416  const G4String tag = Transcode(child->getTagName());
417 
418  if (tag=="volumeref")
419  {
420  logvol = GetVolume(GenerateName(RefRead(child)));
421  }
422  else if (tag=="replicate_along_axis")
423  {
424  if (logvol) { ReplicaRead(child,logvol,number); }
425  }
426  else
427  {
428  G4String error_msg = "Unknown tag in ReplicavolRead: " + tag;
429  G4Exception("G4GDMLReadStructure::ReplicavolRead()",
430  "ReadError", FatalException, error_msg);
431  }
432  }
433 }
Definition: xmlparse.cc:187
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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 651 of file G4GDMLReadStructure.cc.

652 {
653  G4String name;
654  G4LogicalVolume* logvol = 0;
655  G4SurfaceProperty* prop = 0;
656 
657  const xercesc::DOMNamedNodeMap* const attributes
658  = skinsurfaceElement->getAttributes();
659  XMLSize_t attributeCount = attributes->getLength();
660 
661  for (XMLSize_t attribute_index=0;
662  attribute_index<attributeCount; attribute_index++)
663  {
664  xercesc::DOMNode* attribute_node = attributes->item(attribute_index);
665 
666  if (attribute_node->getNodeType() != xercesc::DOMNode::ATTRIBUTE_NODE)
667  { continue; }
668 
669  const xercesc::DOMAttr* const attribute
670  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
671  if (!attribute)
672  {
673  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()",
674  "InvalidRead", FatalException, "No attribute found!");
675  return;
676  }
677  const G4String attName = Transcode(attribute->getName());
678  const G4String attValue = Transcode(attribute->getValue());
679 
680  if (attName=="name")
681  { name = GenerateName(attValue); } else
682  if (attName=="surfaceproperty")
683  { prop = GetSurfaceProperty(GenerateName(attValue)); }
684  }
685 
686  for (xercesc::DOMNode* iter = skinsurfaceElement->getFirstChild();
687  iter != 0; iter = iter->getNextSibling())
688  {
689  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
690 
691  const xercesc::DOMElement* const child
692  = dynamic_cast<xercesc::DOMElement*>(iter);
693  if (!child)
694  {
695  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()",
696  "InvalidRead", FatalException, "No child found!");
697  return;
698  }
699  const G4String tag = Transcode(child->getTagName());
700 
701  if (tag=="volumeref")
702  {
703  logvol = GetVolume(GenerateName(RefRead(child)));
704  }
705  else
706  {
707  G4String error_msg = "Unknown tag in skinsurface: " + tag;
708  G4Exception("G4GDMLReadStructure::SkinsurfaceRead()", "ReadError",
709  FatalException, error_msg);
710  }
711  }
712 
713  new G4LogicalSkinSurface(Strip(name),logvol,prop);
714 }
const XML_Char * name
Definition: expat.h:151
Definition: xmlparse.cc:187
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
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:97

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4GDMLRead.

Definition at line 795 of file G4GDMLReadStructure.cc.

796 {
797  G4cout << "G4GDML: Reading structure..." << G4endl;
798 
799  for (xercesc::DOMNode* iter = structureElement->getFirstChild();
800  iter != 0; iter = iter->getNextSibling())
801  {
802  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
803 
804  const xercesc::DOMElement* const child
805  = dynamic_cast<xercesc::DOMElement*>(iter);
806  if (!child)
807  {
808  G4Exception("G4GDMLReadStructure::StructureRead()",
809  "InvalidRead", FatalException, "No child found!");
810  return;
811  }
812  const G4String tag = Transcode(child->getTagName());
813 
814  if (tag=="bordersurface") { BorderSurfaceRead(child); } else
815  if (tag=="skinsurface") { SkinSurfaceRead(child); } else
816  if (tag=="volume") { VolumeRead(child); } else
817  if (tag=="assembly") { AssemblyRead(child); } else
818  if (tag=="loop") { LoopRead(child,&G4GDMLRead::StructureRead); }
819  else
820  {
821  G4String error_msg = "Unknown tag in structure: " + tag;
822  G4Exception("G4GDMLReadStructure::StructureRead()",
823  "ReadError", FatalException, error_msg);
824  }
825  }
826 }
Definition: xmlparse.cc:187
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:176
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)

Here is the call graph for this function:

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

Implements G4GDMLRead.

Definition at line 717 of file G4GDMLReadStructure.cc.

718 {
719  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
720  iter != 0; iter = iter->getNextSibling())
721  {
722  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
723 
724  const xercesc::DOMElement* const child
725  = dynamic_cast<xercesc::DOMElement*>(iter);
726  if (!child)
727  {
728  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
729  "InvalidRead", FatalException, "No child found!");
730  return;
731  }
732  const G4String tag = Transcode(child->getTagName());
733 
734  if ((tag=="auxiliary") || (tag=="materialref") || (tag=="solidref"))
735  {
736  // These are already processed in VolumeRead()
737  }
738  else if (tag=="paramvol")
739  {
741  }
742  else if (tag=="physvol")
743  {
744  PhysvolRead(child);
745  }
746  else if (tag=="replicavol")
747  {
748  G4int number = 1;
749  const xercesc::DOMNamedNodeMap* const attributes
750  = child->getAttributes();
751  XMLSize_t attributeCount = attributes->getLength();
752  for (XMLSize_t attribute_index=0;
753  attribute_index<attributeCount; attribute_index++)
754  {
755  xercesc::DOMNode* attribute_node
756  = attributes->item(attribute_index);
757  if (attribute_node->getNodeType()!=xercesc::DOMNode::ATTRIBUTE_NODE)
758  {
759  continue;
760  }
761  const xercesc::DOMAttr* const attribute
762  = dynamic_cast<xercesc::DOMAttr*>(attribute_node);
763  if (!attribute)
764  {
765  G4Exception("G4GDMLReadStructure::Volume_contentRead()",
766  "InvalidRead", FatalException, "No attribute found!");
767  return;
768  }
769  const G4String attName = Transcode(attribute->getName());
770  const G4String attValue = Transcode(attribute->getValue());
771  if (attName=="number")
772  {
773  number = eval.EvaluateInteger(attValue);
774  }
775  }
776  ReplicavolRead(child,number);
777  }
778  else if (tag=="divisionvol")
779  {
780  DivisionvolRead(child);
781  }
782  else if (tag=="loop")
783  {
785  }
786  else
787  {
788  G4cout << "Treating unknown GDML tag in volume '" << tag
789  << "' as GDML extension..." << G4endl;
790  }
791  }
792 }
G4int EvaluateInteger(const G4String &)
G4GDMLEvaluator eval
Definition: G4GDMLRead.hh:155
virtual void ParamvolRead(const xercesc::DOMElement *const, G4LogicalVolume *)
virtual void Volume_contentRead(const xercesc::DOMElement *const)=0
Definition: xmlparse.cc:187
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:176
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)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented in G03ColorReader.

Definition at line 573 of file G4GDMLReadStructure.cc.

574 {
575  G4VSolid* solidPtr = 0;
576  G4Material* materialPtr = 0;
577  G4GDMLAuxListType auxList;
578 
579  XMLCh *name_attr = xercesc::XMLString::transcode("name");
580  const G4String name = Transcode(volumeElement->getAttribute(name_attr));
581  xercesc::XMLString::release(&name_attr);
582 
583  for (xercesc::DOMNode* iter = volumeElement->getFirstChild();
584  iter != 0; iter = iter->getNextSibling())
585  {
586  if (iter->getNodeType() != xercesc::DOMNode::ELEMENT_NODE) { continue; }
587 
588  const xercesc::DOMElement* const child
589  = dynamic_cast<xercesc::DOMElement*>(iter);
590  if (!child)
591  {
592  G4Exception("G4GDMLReadStructure::VolumeRead()",
593  "InvalidRead", FatalException, "No child found!");
594  return;
595  }
596  const G4String tag = Transcode(child->getTagName());
597 
598  if (tag=="auxiliary")
599  { auxList.push_back(AuxiliaryRead(child)); } else
600  if (tag=="materialref")
601  { materialPtr = GetMaterial(GenerateName(RefRead(child),true)); } else
602  if (tag=="solidref")
603  { solidPtr = GetSolid(GenerateName(RefRead(child))); }
604  }
605 
606  pMotherLogical = new G4LogicalVolume(solidPtr,materialPtr,
607  GenerateName(name),0,0,0);
608 
609  if (!auxList.empty()) { auxMap[pMotherLogical] = auxList; }
610 
611  Volume_contentRead(volumeElement);
612 }
const XML_Char * name
Definition: expat.h:151
G4Material * GetMaterial(const G4String &, G4bool verbose=true) const
Definition: xmlparse.cc:187
G4String Transcode(const XMLCh *const)
Definition: G4GDMLRead.cc:55
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
G4String GenerateName(const G4String &name, G4bool strip=false)
Definition: G4GDMLRead.cc:68
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
G4GDMLAuxMapType auxMap
G4GDMLAuxStructType AuxiliaryRead(const xercesc::DOMElement *const auxElem)
Definition: G4GDMLRead.cc:262
virtual void Volume_contentRead(const xercesc::DOMElement *const)

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

G4GDMLAssemblyMapType G4GDMLReadStructure::assemblyMap
protected

Definition at line 93 of file G4GDMLReadStructure.hh.

G4GDMLAuxMapType G4GDMLReadStructure::auxMap
protected

Definition at line 92 of file G4GDMLReadStructure.hh.

G4LogicalVolume* G4GDMLReadStructure::pMotherLogical
protected

Definition at line 94 of file G4GDMLReadStructure.hh.

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

Definition at line 95 of file G4GDMLReadStructure.hh.

G4bool G4GDMLReadStructure::strip
protected

Definition at line 96 of file G4GDMLReadStructure.hh.


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