Geant4  10.02.p03
G4GDMLWriteSolids Class Reference

#include <G4GDMLWriteSolids.hh>

Inheritance diagram for G4GDMLWriteSolids:
Collaboration diagram for G4GDMLWriteSolids:

Classes

class  G4ThreeVectorCompare
 

Public Member Functions

virtual void AddSolid (const G4VSolid *const)
 
virtual void SolidsWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWriteMaterials
void AddIsotope (const G4Isotope *const)
 
void AddElement (const G4Element *const)
 
void AddMaterial (const G4Material *const)
 
virtual void MaterialsWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWriteDefine
G4ThreeVector GetAngles (const G4RotationMatrix &)
 
void ScaleWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
 
void RotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void PositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void FirstrotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void FirstpositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void AddPosition (const G4String &name, const G4ThreeVector &pos)
 
virtual void DefineWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWrite
G4Transform3D Write (const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
 
void AddModule (const G4VPhysicalVolume *const topVol)
 
void AddModule (const G4int depth)
 
void AddAuxiliary (G4GDMLAuxStructType myaux)
 
virtual void StructureWrite (xercesc::DOMElement *)=0
 
virtual G4Transform3D TraverseVolumeTree (const G4LogicalVolume *const, const G4int)=0
 
virtual void SurfacesWrite ()=0
 
virtual void SetupWrite (xercesc::DOMElement *, const G4LogicalVolume *const)=0
 
virtual void ExtensionWrite (xercesc::DOMElement *)
 
virtual void UserinfoWrite (xercesc::DOMElement *)
 
virtual void AddExtension (xercesc::DOMElement *, const G4LogicalVolume *const)
 
G4String GenerateName (const G4String &, const void *const)
 

Protected Member Functions

 G4GDMLWriteSolids ()
 
virtual ~G4GDMLWriteSolids ()
 
void MultiUnionWrite (xercesc::DOMElement *solElement, const G4MultiUnion *const)
 
void BooleanWrite (xercesc::DOMElement *, const G4BooleanSolid *const)
 
void BoxWrite (xercesc::DOMElement *, const G4Box *const)
 
void ConeWrite (xercesc::DOMElement *, const G4Cons *const)
 
void ElconeWrite (xercesc::DOMElement *, const G4EllipticalCone *const)
 
void EllipsoidWrite (xercesc::DOMElement *, const G4Ellipsoid *const)
 
void EltubeWrite (xercesc::DOMElement *, const G4EllipticalTube *const)
 
void XtruWrite (xercesc::DOMElement *, const G4ExtrudedSolid *const)
 
void HypeWrite (xercesc::DOMElement *, const G4Hype *const)
 
void OrbWrite (xercesc::DOMElement *, const G4Orb *const)
 
void ParaWrite (xercesc::DOMElement *, const G4Para *const)
 
void ParaboloidWrite (xercesc::DOMElement *, const G4Paraboloid *const)
 
void PolyconeWrite (xercesc::DOMElement *, const G4Polycone *const)
 
void GenericPolyconeWrite (xercesc::DOMElement *, const G4GenericPolycone *const)
 
void PolyhedraWrite (xercesc::DOMElement *, const G4Polyhedra *const)
 
void SphereWrite (xercesc::DOMElement *, const G4Sphere *const)
 
void TessellatedWrite (xercesc::DOMElement *, const G4TessellatedSolid *const)
 
void TetWrite (xercesc::DOMElement *, const G4Tet *const)
 
void TorusWrite (xercesc::DOMElement *, const G4Torus *const)
 
void GenTrapWrite (xercesc::DOMElement *, const G4GenericTrap *const)
 
void TrapWrite (xercesc::DOMElement *, const G4Trap *const)
 
void TrdWrite (xercesc::DOMElement *, const G4Trd *const)
 
void TubeWrite (xercesc::DOMElement *, const G4Tubs *const)
 
void CutTubeWrite (xercesc::DOMElement *, const G4CutTubs *const)
 
void TwistedboxWrite (xercesc::DOMElement *, const G4TwistedBox *const)
 
void TwistedtrapWrite (xercesc::DOMElement *, const G4TwistedTrap *const)
 
void TwistedtrdWrite (xercesc::DOMElement *, const G4TwistedTrd *const)
 
void TwistedtubsWrite (xercesc::DOMElement *, const G4TwistedTubs *const)
 
void ZplaneWrite (xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
 
void RZPointWrite (xercesc::DOMElement *, const G4double &, const G4double &)
 
void OpticalSurfaceWrite (xercesc::DOMElement *, const G4OpticalSurface *const)
 
- Protected Member Functions inherited from G4GDMLWriteMaterials
 G4GDMLWriteMaterials ()
 
virtual ~G4GDMLWriteMaterials ()
 
void AtomWrite (xercesc::DOMElement *, const G4double &)
 
void DWrite (xercesc::DOMElement *, const G4double &)
 
void PWrite (xercesc::DOMElement *, const G4double &)
 
void TWrite (xercesc::DOMElement *, const G4double &)
 
void MEEWrite (xercesc::DOMElement *, const G4double &)
 
void IsotopeWrite (const G4Isotope *const)
 
void ElementWrite (const G4Element *const)
 
void MaterialWrite (const G4Material *const)
 
void PropertyWrite (xercesc::DOMElement *, const G4Material *const)
 
void PropertyVectorWrite (const G4String &, const G4PhysicsOrderedFreeVector *const)
 
- Protected Member Functions inherited from G4GDMLWriteDefine
 G4GDMLWriteDefine ()
 
virtual ~G4GDMLWriteDefine ()
 
void Scale_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Rotation_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Position_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
- Protected Member Functions inherited from G4GDMLWrite
 G4GDMLWrite ()
 
virtual ~G4GDMLWrite ()
 
VolumeMapTypeVolumeMap ()
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4String &)
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4double &)
 
xercesc::DOMElement * NewElement (const G4String &)
 
G4String Modularize (const G4VPhysicalVolume *const topvol, const G4int depth)
 
void AddAuxInfo (G4GDMLAuxListType *auxInfoList, xercesc::DOMElement *element)
 
G4bool FileExists (const G4String &) const
 
PhysVolumeMapTypePvolumeMap ()
 
DepthMapTypeDepthMap ()
 

Protected Attributes

std::vector< const G4VSolid * > solidList
 
xercesc::DOMElement * solidsElement
 
- Protected Attributes inherited from G4GDMLWriteMaterials
std::vector< const G4Isotope * > isotopeList
 
std::vector< const G4Element * > elementList
 
std::vector< const G4Material * > materialList
 
xercesc::DOMElement * materialsElement
 
- Protected Attributes inherited from G4GDMLWriteDefine
xercesc::DOMElement * defineElement
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc
 
xercesc::DOMElement * extElement
 
xercesc::DOMElement * userinfoElement
 
XMLCh tempStr [10000]
 
G4GDMLAuxListType auxList
 

Static Protected Attributes

static const G4int maxTransforms = 8
 
- Static Protected Attributes inherited from G4GDMLWriteDefine
static const G4double kRelativePrecision = DBL_EPSILON
 
static const G4double kAngularPrecision = DBL_EPSILON
 
static const G4double kLinearPrecision = DBL_EPSILON
 
- Static Protected Attributes inherited from G4GDMLWrite
static G4bool addPointerToName = true
 

Additional Inherited Members

- Static Public Member Functions inherited from G4GDMLWrite
static void SetAddPointerToName (G4bool)
 

Detailed Description

Definition at line 78 of file G4GDMLWriteSolids.hh.

Constructor & Destructor Documentation

◆ G4GDMLWriteSolids()

G4GDMLWriteSolids::G4GDMLWriteSolids ( )
protected

Definition at line 72 of file G4GDMLWriteSolids.cc.

74 {
75 }
xercesc::DOMElement * solidsElement
Here is the caller graph for this function:

◆ ~G4GDMLWriteSolids()

G4GDMLWriteSolids::~G4GDMLWriteSolids ( )
protectedvirtual

Definition at line 77 of file G4GDMLWriteSolids.cc.

78 {
79 }
Here is the call graph for this function:
Here is the caller graph for this function:

Member Function Documentation

◆ AddSolid()

void G4GDMLWriteSolids::AddSolid ( const G4VSolid * const  solidPtr)
virtual

Definition at line 1032 of file G4GDMLWriteSolids.cc.

1033 {
1034  for (size_t i=0; i<solidList.size(); i++) // Check if solid is
1035  { // already in the list!
1036  if (solidList[i] == solidPtr) { return; }
1037  }
1038 
1039  solidList.push_back(solidPtr);
1040 
1041  if (const G4BooleanSolid* const booleanPtr
1042  = dynamic_cast<const G4BooleanSolid*>(solidPtr))
1043  { BooleanWrite(solidsElement,booleanPtr); } else
1044  if (solidPtr->GetEntityType()=="G4MultiUnion")
1045  { const G4MultiUnion* const munionPtr
1046  = static_cast<const G4MultiUnion*>(solidPtr);
1047  MultiUnionWrite(solidsElement,munionPtr); } else
1048  if (solidPtr->GetEntityType()=="G4Box")
1049  { const G4Box* const boxPtr
1050  = static_cast<const G4Box*>(solidPtr);
1051  BoxWrite(solidsElement,boxPtr); } else
1052  if (solidPtr->GetEntityType()=="G4Cons")
1053  { const G4Cons* const conePtr
1054  = static_cast<const G4Cons*>(solidPtr);
1055  ConeWrite(solidsElement,conePtr); } else
1056  if (solidPtr->GetEntityType()=="G4EllipticalCone")
1057  { const G4EllipticalCone* const elconePtr
1058  = static_cast<const G4EllipticalCone*>(solidPtr);
1059  ElconeWrite(solidsElement,elconePtr); } else
1060  if (solidPtr->GetEntityType()=="G4Ellipsoid")
1061  { const G4Ellipsoid* const ellipsoidPtr
1062  = static_cast<const G4Ellipsoid*>(solidPtr);
1063  EllipsoidWrite(solidsElement,ellipsoidPtr); } else
1064  if (solidPtr->GetEntityType()=="G4EllipticalTube")
1065  { const G4EllipticalTube* const eltubePtr
1066  = static_cast<const G4EllipticalTube*>(solidPtr);
1067  EltubeWrite(solidsElement,eltubePtr); } else
1068  if (solidPtr->GetEntityType()=="G4ExtrudedSolid")
1069  { const G4ExtrudedSolid* const xtruPtr
1070  = static_cast<const G4ExtrudedSolid*>(solidPtr);
1071  XtruWrite(solidsElement,xtruPtr); } else
1072  if (solidPtr->GetEntityType()=="G4Hype")
1073  { const G4Hype* const hypePtr
1074  = static_cast<const G4Hype*>(solidPtr);
1075  HypeWrite(solidsElement,hypePtr); } else
1076  if (solidPtr->GetEntityType()=="G4Orb")
1077  { const G4Orb* const orbPtr
1078  = static_cast<const G4Orb*>(solidPtr);
1079  OrbWrite(solidsElement,orbPtr); } else
1080  if (solidPtr->GetEntityType()=="G4Para")
1081  { const G4Para* const paraPtr
1082  = static_cast<const G4Para*>(solidPtr);
1083  ParaWrite(solidsElement,paraPtr); } else
1084  if (solidPtr->GetEntityType()=="G4Paraboloid")
1085  { const G4Paraboloid* const paraboloidPtr
1086  = static_cast<const G4Paraboloid*>(solidPtr);
1087  ParaboloidWrite(solidsElement,paraboloidPtr); } else
1088  if (solidPtr->GetEntityType()=="G4Polycone")
1089  { const G4Polycone* const polyconePtr
1090  = static_cast<const G4Polycone*>(solidPtr);
1091  PolyconeWrite(solidsElement,polyconePtr); } else
1092  if (solidPtr->GetEntityType()=="G4GenericPolycone")
1093  { const G4GenericPolycone* const genpolyconePtr
1094  = static_cast<const G4GenericPolycone*>(solidPtr);
1095  GenericPolyconeWrite(solidsElement,genpolyconePtr); } else
1096  if (solidPtr->GetEntityType()=="G4Polyhedra")
1097  { const G4Polyhedra* const polyhedraPtr
1098  = static_cast<const G4Polyhedra*>(solidPtr);
1099  PolyhedraWrite(solidsElement,polyhedraPtr); } else
1100  if (solidPtr->GetEntityType()=="G4Sphere")
1101  { const G4Sphere* const spherePtr
1102  = static_cast<const G4Sphere*>(solidPtr);
1103  SphereWrite(solidsElement,spherePtr); } else
1104  if (solidPtr->GetEntityType()=="G4TessellatedSolid")
1105  { const G4TessellatedSolid* const tessellatedPtr
1106  = static_cast<const G4TessellatedSolid*>(solidPtr);
1107  TessellatedWrite(solidsElement,tessellatedPtr); } else
1108  if (solidPtr->GetEntityType()=="G4Tet")
1109  { const G4Tet* const tetPtr
1110  = static_cast<const G4Tet*>(solidPtr);
1111  TetWrite(solidsElement,tetPtr); } else
1112  if (solidPtr->GetEntityType()=="G4Torus")
1113  { const G4Torus* const torusPtr
1114  = static_cast<const G4Torus*>(solidPtr);
1115  TorusWrite(solidsElement,torusPtr); } else
1116  if (solidPtr->GetEntityType()=="G4GenericTrap")
1117  { const G4GenericTrap* const gtrapPtr
1118  = static_cast<const G4GenericTrap*>(solidPtr);
1119  GenTrapWrite(solidsElement,gtrapPtr); } else
1120  if (solidPtr->GetEntityType()=="G4Trap")
1121  { const G4Trap* const trapPtr
1122  = static_cast<const G4Trap*>(solidPtr);
1123  TrapWrite(solidsElement,trapPtr); } else
1124  if (solidPtr->GetEntityType()=="G4Trd")
1125  { const G4Trd* const trdPtr
1126  = static_cast<const G4Trd*>(solidPtr);
1127  TrdWrite(solidsElement,trdPtr); } else
1128  if (solidPtr->GetEntityType()=="G4Tubs")
1129  { const G4Tubs* const tubePtr
1130  = static_cast<const G4Tubs*>(solidPtr);
1131  TubeWrite(solidsElement,tubePtr); } else
1132  if (solidPtr->GetEntityType()=="G4CutTubs")
1133  { const G4CutTubs* const cuttubePtr
1134  = static_cast<const G4CutTubs*>(solidPtr);
1135  CutTubeWrite(solidsElement,cuttubePtr); } else
1136  if (solidPtr->GetEntityType()=="G4TwistedBox")
1137  { const G4TwistedBox* const twistedboxPtr
1138  = static_cast<const G4TwistedBox*>(solidPtr);
1139  TwistedboxWrite(solidsElement,twistedboxPtr); } else
1140  if (solidPtr->GetEntityType()=="G4TwistedTrap")
1141  { const G4TwistedTrap* const twistedtrapPtr
1142  = static_cast<const G4TwistedTrap*>(solidPtr);
1143  TwistedtrapWrite(solidsElement,twistedtrapPtr); } else
1144  if (solidPtr->GetEntityType()=="G4TwistedTrd")
1145  { const G4TwistedTrd* const twistedtrdPtr
1146  = static_cast<const G4TwistedTrd*>(solidPtr);
1147  TwistedtrdWrite(solidsElement,twistedtrdPtr); } else
1148  if (solidPtr->GetEntityType()=="G4TwistedTubs")
1149  { const G4TwistedTubs* const twistedtubsPtr
1150  = static_cast<const G4TwistedTubs*>(solidPtr);
1151  TwistedtubsWrite(solidsElement,twistedtubsPtr); }
1152  else
1153  {
1154  G4String error_msg = "Unknown solid: " + solidPtr->GetName()
1155  + "; Type: " + solidPtr->GetEntityType();
1156  G4Exception("G4GDMLWriteSolids::AddSolid()", "WriteError",
1157  FatalException, error_msg);
1158  }
1159 }
Definition: G4Para.hh:77
void BoxWrite(xercesc::DOMElement *, const G4Box *const)
void PolyhedraWrite(xercesc::DOMElement *, const G4Polyhedra *const)
void TessellatedWrite(xercesc::DOMElement *, const G4TessellatedSolid *const)
Definition: G4Box.hh:64
Definition: G4Tubs.hh:85
void TrdWrite(xercesc::DOMElement *, const G4Trd *const)
void XtruWrite(xercesc::DOMElement *, const G4ExtrudedSolid *const)
Definition: G4Trd.hh:72
virtual G4GeometryType GetEntityType() const =0
void TwistedtrdWrite(xercesc::DOMElement *, const G4TwistedTrd *const)
G4String GetName() const
void TwistedboxWrite(xercesc::DOMElement *, const G4TwistedBox *const)
void MultiUnionWrite(xercesc::DOMElement *solElement, const G4MultiUnion *const)
void OrbWrite(xercesc::DOMElement *, const G4Orb *const)
void TwistedtubsWrite(xercesc::DOMElement *, const G4TwistedTubs *const)
void TrapWrite(xercesc::DOMElement *, const G4Trap *const)
void ConeWrite(xercesc::DOMElement *, const G4Cons *const)
Definition: G4Tet.hh:65
void ElconeWrite(xercesc::DOMElement *, const G4EllipticalCone *const)
Definition: G4Hype.hh:67
void ParaboloidWrite(xercesc::DOMElement *, const G4Paraboloid *const)
void PolyconeWrite(xercesc::DOMElement *, const G4Polycone *const)
Definition: G4Cons.hh:83
void GenericPolyconeWrite(xercesc::DOMElement *, const G4GenericPolycone *const)
void CutTubeWrite(xercesc::DOMElement *, const G4CutTubs *const)
Definition: G4Orb.hh:61
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void G4MultiUnion
Definition: G4MultiUnion.hh:53
void HypeWrite(xercesc::DOMElement *, const G4Hype *const)
void EltubeWrite(xercesc::DOMElement *, const G4EllipticalTube *const)
void ParaWrite(xercesc::DOMElement *, const G4Para *const)
void EllipsoidWrite(xercesc::DOMElement *, const G4Ellipsoid *const)
void TetWrite(xercesc::DOMElement *, const G4Tet *const)
void TubeWrite(xercesc::DOMElement *, const G4Tubs *const)
void BooleanWrite(xercesc::DOMElement *, const G4BooleanSolid *const)
void TorusWrite(xercesc::DOMElement *, const G4Torus *const)
void TwistedtrapWrite(xercesc::DOMElement *, const G4TwistedTrap *const)
void GenTrapWrite(xercesc::DOMElement *, const G4GenericTrap *const)
xercesc::DOMElement * solidsElement
void SphereWrite(xercesc::DOMElement *, const G4Sphere *const)
std::vector< const G4VSolid * > solidList
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BooleanWrite()

void G4GDMLWriteSolids::BooleanWrite ( xercesc::DOMElement *  solElement,
const G4BooleanSolid * const  boolean 
)
protected

Definition at line 153 of file G4GDMLWriteSolids.cc.

155 {
156  G4int displaced=0;
157 
158  G4String tag("undefined");
159  if (dynamic_cast<const G4IntersectionSolid*>(boolean))
160  { tag = "intersection"; } else
161  if (dynamic_cast<const G4SubtractionSolid*>(boolean))
162  { tag = "subtraction"; } else
163  if (dynamic_cast<const G4UnionSolid*>(boolean))
164  { tag = "union"; }
165 
166  G4VSolid* firstPtr = const_cast<G4VSolid*>(boolean->GetConstituentSolid(0));
167  G4VSolid* secondPtr = const_cast<G4VSolid*>(boolean->GetConstituentSolid(1));
168 
169  G4ThreeVector firstpos,firstrot,pos,rot;
170 
171  // Solve possible displacement of referenced solids!
172  //
173  while (true)
174  {
175  if ( displaced>8 )
176  {
177  G4String ErrorMessage = "The referenced solid '"
178  + firstPtr->GetName() +
179  + "in the Boolean shape '" +
180  + boolean->GetName() +
181  + "' was displaced too many times!";
182  G4Exception("G4GDMLWriteSolids::BooleanWrite()",
183  "InvalidSetup", FatalException, ErrorMessage);
184  }
185 
186  if (G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(firstPtr))
187  {
188  firstpos += disp->GetObjectTranslation();
189  firstrot += GetAngles(disp->GetObjectRotation());
190  firstPtr = disp->GetConstituentMovedSolid();
191  displaced++;
192  continue;
193  }
194  break;
195  }
196  displaced = 0;
197  while (true)
198  {
199  if ( displaced>maxTransforms )
200  {
201  G4String ErrorMessage = "The referenced solid '"
202  + secondPtr->GetName() +
203  + "in the Boolean shape '" +
204  + boolean->GetName() +
205  + "' was displaced too many times!";
206  G4Exception("G4GDMLWriteSolids::BooleanWrite()",
207  "InvalidSetup", FatalException, ErrorMessage);
208  }
209 
210  if (G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(secondPtr))
211  {
212  pos += disp->GetObjectTranslation();
213  rot += GetAngles(disp->GetObjectRotation());
214  secondPtr = disp->GetConstituentMovedSolid();
215  displaced++;
216  continue;
217  }
218  break;
219  }
220 
221  AddSolid(firstPtr); // At first add the constituent solids!
222  AddSolid(secondPtr);
223 
224  const G4String& name = GenerateName(boolean->GetName(),boolean);
225  const G4String& firstref = GenerateName(firstPtr->GetName(),firstPtr);
226  const G4String& secondref = GenerateName(secondPtr->GetName(),secondPtr);
227 
228  xercesc::DOMElement* booleanElement = NewElement(tag);
229  booleanElement->setAttributeNode(NewAttribute("name",name));
230  xercesc::DOMElement* firstElement = NewElement("first");
231  firstElement->setAttributeNode(NewAttribute("ref",firstref));
232  booleanElement->appendChild(firstElement);
233  xercesc::DOMElement* secondElement = NewElement("second");
234  secondElement->setAttributeNode(NewAttribute("ref",secondref));
235  booleanElement->appendChild(secondElement);
236  solElement->appendChild(booleanElement);
237  // Add the boolean solid AFTER the constituent solids!
238 
239  if ( (std::fabs(pos.x()) > kLinearPrecision)
240  || (std::fabs(pos.y()) > kLinearPrecision)
241  || (std::fabs(pos.z()) > kLinearPrecision) )
242  {
243  PositionWrite(booleanElement,name+"_pos",pos);
244  }
245 
246  if ( (std::fabs(rot.x()) > kAngularPrecision)
247  || (std::fabs(rot.y()) > kAngularPrecision)
248  || (std::fabs(rot.z()) > kAngularPrecision) )
249  {
250  RotationWrite(booleanElement,name+"_rot",rot);
251  }
252 
253  if ( (std::fabs(firstpos.x()) > kLinearPrecision)
254  || (std::fabs(firstpos.y()) > kLinearPrecision)
255  || (std::fabs(firstpos.z()) > kLinearPrecision) )
256  {
257  FirstpositionWrite(booleanElement,name+"_fpos",firstpos);
258  }
259 
260  if ( (std::fabs(firstrot.x()) > kAngularPrecision)
261  || (std::fabs(firstrot.y()) > kAngularPrecision)
262  || (std::fabs(firstrot.z()) > kAngularPrecision) )
263  {
264  FirstrotationWrite(booleanElement,name+"_frot",firstrot);
265  }
266 }
static const G4double kAngularPrecision
G4String name
Definition: TRTMaterials.hh:40
Definition: xmlparse.cc:187
void FirstpositionWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
int G4int
Definition: G4Types.hh:78
G4String GetName() const
static const G4int maxTransforms
void FirstrotationWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
char boolean
Definition: csz_inflate.cc:208
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void RotationWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double x() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const G4double kLinearPrecision
double y() const
G4ThreeVector GetAngles(const G4RotationMatrix &)
double z() const
void PositionWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
static const G4double pos
virtual void AddSolid(const G4VSolid *const)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BoxWrite()

void G4GDMLWriteSolids::BoxWrite ( xercesc::DOMElement *  solElement,
const G4Box * const  box 
)
protected

Definition at line 269 of file G4GDMLWriteSolids.cc.

270 {
271  const G4String& name = GenerateName(box->GetName(),box);
272 
273  xercesc::DOMElement* boxElement = NewElement("box");
274  boxElement->setAttributeNode(NewAttribute("name",name));
275  boxElement->setAttributeNode(NewAttribute("x",2.0*box->GetXHalfLength()/mm));
276  boxElement->setAttributeNode(NewAttribute("y",2.0*box->GetYHalfLength()/mm));
277  boxElement->setAttributeNode(NewAttribute("z",2.0*box->GetZHalfLength()/mm));
278  boxElement->setAttributeNode(NewAttribute("lunit","mm"));
279  solElement->appendChild(boxElement);
280 }
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetXHalfLength() const
G4double GetZHalfLength() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetYHalfLength() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConeWrite()

void G4GDMLWriteSolids::ConeWrite ( xercesc::DOMElement *  solElement,
const G4Cons * const  cone 
)
protected

Definition at line 283 of file G4GDMLWriteSolids.cc.

284 {
285  const G4String& name = GenerateName(cone->GetName(),cone);
286 
287  xercesc::DOMElement* coneElement = NewElement("cone");
288  coneElement->setAttributeNode(NewAttribute("name",name));
289  coneElement->
290  setAttributeNode(NewAttribute("rmin1",cone->GetInnerRadiusMinusZ()/mm));
291  coneElement->
292  setAttributeNode(NewAttribute("rmax1",cone->GetOuterRadiusMinusZ()/mm));
293  coneElement->
294  setAttributeNode(NewAttribute("rmin2",cone->GetInnerRadiusPlusZ()/mm));
295  coneElement->
296  setAttributeNode(NewAttribute("rmax2",cone->GetOuterRadiusPlusZ()/mm));
297  coneElement->
298  setAttributeNode(NewAttribute("z",2.0*cone->GetZHalfLength()/mm));
299  coneElement->
300  setAttributeNode(NewAttribute("startphi",cone->GetStartPhiAngle()/degree));
301  coneElement->
302  setAttributeNode(NewAttribute("deltaphi",cone->GetDeltaPhiAngle()/degree));
303  coneElement->setAttributeNode(NewAttribute("aunit","deg"));
304  coneElement->setAttributeNode(NewAttribute("lunit","mm"));
305  solElement->appendChild(coneElement);
306 }
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4String name
Definition: TRTMaterials.hh:40
G4double GetInnerRadiusMinusZ() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetInnerRadiusPlusZ() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetOuterRadiusMinusZ() const
static const double degree
Definition: G4SIunits.hh:143
G4double GetZHalfLength() const
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CutTubeWrite()

void G4GDMLWriteSolids::CutTubeWrite ( xercesc::DOMElement *  solElement,
const G4CutTubs * const  cuttube 
)
protected

Definition at line 846 of file G4GDMLWriteSolids.cc.

847 {
848  const G4String& name = GenerateName(cuttube->GetName(),cuttube);
849 
850  xercesc::DOMElement* cuttubeElement = NewElement("cutTube");
851  cuttubeElement->setAttributeNode(NewAttribute("name",name));
852  cuttubeElement->setAttributeNode(NewAttribute("rmin",
853  cuttube->GetInnerRadius()/mm));
854  cuttubeElement->setAttributeNode(NewAttribute("rmax",
855  cuttube->GetOuterRadius()/mm));
856  cuttubeElement->setAttributeNode(NewAttribute("z",
857  2.0*cuttube->GetZHalfLength()/mm));
858  cuttubeElement->setAttributeNode(NewAttribute("startphi",
859  cuttube->GetStartPhiAngle()/degree));
860  cuttubeElement->setAttributeNode(NewAttribute("deltaphi",
861  cuttube->GetDeltaPhiAngle()/degree));
862  cuttubeElement->setAttributeNode(NewAttribute("lowX",
863  cuttube->GetLowNorm().getX()/mm));
864  cuttubeElement->setAttributeNode(NewAttribute("lowY",
865  cuttube->GetLowNorm().getY()/mm));
866  cuttubeElement->setAttributeNode(NewAttribute("lowZ",
867  cuttube->GetLowNorm().getZ()/mm));
868  cuttubeElement->setAttributeNode(NewAttribute("highX",
869  cuttube->GetHighNorm().getX()/mm));
870  cuttubeElement->setAttributeNode(NewAttribute("highY",
871  cuttube->GetHighNorm().getY()/mm));
872  cuttubeElement->setAttributeNode(NewAttribute("highZ",
873  cuttube->GetHighNorm().getZ()/mm));
874  cuttubeElement->setAttributeNode(NewAttribute("aunit","deg"));
875  cuttubeElement->setAttributeNode(NewAttribute("lunit","mm"));
876  solElement->appendChild(cuttubeElement);
877 }
G4double GetStartPhiAngle() const
G4String name
Definition: TRTMaterials.hh:40
G4double GetInnerRadius() const
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
double getY() const
double getX() const
double getZ() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetOuterRadius() const
G4double GetDeltaPhiAngle() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double degree
Definition: G4SIunits.hh:143
G4ThreeVector GetLowNorm() const
static const double mm
Definition: G4SIunits.hh:114
G4ThreeVector GetHighNorm() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ElconeWrite()

void G4GDMLWriteSolids::ElconeWrite ( xercesc::DOMElement *  solElement,
const G4EllipticalCone * const  elcone 
)
protected

Definition at line 309 of file G4GDMLWriteSolids.cc.

311 {
312  const G4String& name = GenerateName(elcone->GetName(),elcone);
313 
314  xercesc::DOMElement* elconeElement = NewElement("elcone");
315  elconeElement->setAttributeNode(NewAttribute("name",name));
316  elconeElement->setAttributeNode(NewAttribute("dx",elcone->GetSemiAxisX()/mm));
317  elconeElement->setAttributeNode(NewAttribute("dy",elcone->GetSemiAxisY()/mm));
318  elconeElement->setAttributeNode(NewAttribute("zmax",elcone->GetZMax()/mm));
319  elconeElement->setAttributeNode(NewAttribute("zcut",elcone->GetZTopCut()/mm));
320  elconeElement->setAttributeNode(NewAttribute("lunit","mm"));
321  solElement->appendChild(elconeElement);
322 }
G4double GetZMax() const
G4double GetZTopCut() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetSemiAxisY() const
static const double mm
Definition: G4SIunits.hh:114
G4double GetSemiAxisX() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EllipsoidWrite()

void G4GDMLWriteSolids::EllipsoidWrite ( xercesc::DOMElement *  solElement,
const G4Ellipsoid * const  ellipsoid 
)
protected

Definition at line 325 of file G4GDMLWriteSolids.cc.

327 {
328  const G4String& name = GenerateName(ellipsoid->GetName(),ellipsoid);
329 
330  xercesc::DOMElement* ellipsoidElement = NewElement("ellipsoid");
331  ellipsoidElement->setAttributeNode(NewAttribute("name",name));
332  ellipsoidElement->
333  setAttributeNode(NewAttribute("ax",ellipsoid->GetSemiAxisMax(0)/mm));
334  ellipsoidElement->
335  setAttributeNode(NewAttribute("by",ellipsoid->GetSemiAxisMax(1)/mm));
336  ellipsoidElement->
337  setAttributeNode(NewAttribute("cz",ellipsoid->GetSemiAxisMax(2)/mm));
338  ellipsoidElement->
339  setAttributeNode(NewAttribute("zcut1",ellipsoid->GetZBottomCut()/mm));
340  ellipsoidElement->
341  setAttributeNode(NewAttribute("zcut2",ellipsoid->GetZTopCut()/mm));
342  ellipsoidElement->
343  setAttributeNode(NewAttribute("lunit","mm"));
344  solElement->appendChild(ellipsoidElement);
345 }
G4double GetSemiAxisMax(G4int i) const
G4String name
Definition: TRTMaterials.hh:40
G4double GetZBottomCut() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetZTopCut() const
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ EltubeWrite()

void G4GDMLWriteSolids::EltubeWrite ( xercesc::DOMElement *  solElement,
const G4EllipticalTube * const  eltube 
)
protected

Definition at line 348 of file G4GDMLWriteSolids.cc.

350 {
351  const G4String& name = GenerateName(eltube->GetName(),eltube);
352 
353  xercesc::DOMElement* eltubeElement = NewElement("eltube");
354  eltubeElement->setAttributeNode(NewAttribute("name",name));
355  eltubeElement->setAttributeNode(NewAttribute("dx",eltube->GetDx()/mm));
356  eltubeElement->setAttributeNode(NewAttribute("dy",eltube->GetDy()/mm));
357  eltubeElement->setAttributeNode(NewAttribute("dz",eltube->GetDz()/mm));
358  eltubeElement->setAttributeNode(NewAttribute("lunit","mm"));
359  solElement->appendChild(eltubeElement);
360 }
G4double GetDx() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetDz() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetDy() const
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GenericPolyconeWrite()

void G4GDMLWriteSolids::GenericPolyconeWrite ( xercesc::DOMElement *  solElement,
const G4GenericPolycone * const  polycone 
)
protected

Definition at line 513 of file G4GDMLWriteSolids.cc.

515 {
516  const G4String& name = GenerateName(polycone->GetName(),polycone);
517  xercesc::DOMElement* polyconeElement = NewElement("genericPolycone");
518  const G4double startPhi=polycone->GetStartPhi();
519  polyconeElement->setAttributeNode(NewAttribute("name",name));
520  polyconeElement->setAttributeNode(NewAttribute("startphi",
521  startPhi/degree));
522  polyconeElement->setAttributeNode(NewAttribute("deltaphi",
523  (polycone->GetEndPhi()-startPhi)/degree));
524  polyconeElement->setAttributeNode(NewAttribute("aunit","deg"));
525  polyconeElement->setAttributeNode(NewAttribute("lunit","mm"));
526  solElement->appendChild(polyconeElement);
527 
528  const size_t num_rzpoints = polycone->GetNumRZCorner();
529  for (size_t i=0; i<num_rzpoints; i++)
530  {
531  const G4double r_point=polycone->GetCorner(i).r;
532  const G4double z_point=polycone->GetCorner(i).z;
533  RZPointWrite(polyconeElement,r_point,z_point);
534  }
535 
536 }
G4double GetStartPhi() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
void RZPointWrite(xercesc::DOMElement *, const G4double &, const G4double &)
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double degree
Definition: G4SIunits.hh:143
G4int GetNumRZCorner() const
G4PolyconeSideRZ GetCorner(G4int index) const
double G4double
Definition: G4Types.hh:76
G4double GetEndPhi() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GenTrapWrite()

void G4GDMLWriteSolids::GenTrapWrite ( xercesc::DOMElement *  solElement,
const G4GenericTrap * const  gtrap 
)
protected

Definition at line 735 of file G4GDMLWriteSolids.cc.

737 {
738  const G4String& name = GenerateName(gtrap->GetName(),gtrap);
739 
740  std::vector<G4TwoVector> vertices = gtrap->GetVertices();
741 
742  xercesc::DOMElement* gtrapElement = NewElement("arb8");
743  gtrapElement->setAttributeNode(NewAttribute("name",name));
744  gtrapElement->setAttributeNode(NewAttribute("dz",
745  gtrap->GetZHalfLength()/mm));
746  gtrapElement->setAttributeNode(NewAttribute("v1x", vertices[0].x()));
747  gtrapElement->setAttributeNode(NewAttribute("v1y", vertices[0].y()));
748  gtrapElement->setAttributeNode(NewAttribute("v2x", vertices[1].x()));
749  gtrapElement->setAttributeNode(NewAttribute("v2y", vertices[1].y()));
750  gtrapElement->setAttributeNode(NewAttribute("v3x", vertices[2].x()));
751  gtrapElement->setAttributeNode(NewAttribute("v3y", vertices[2].y()));
752  gtrapElement->setAttributeNode(NewAttribute("v4x", vertices[3].x()));
753  gtrapElement->setAttributeNode(NewAttribute("v4y", vertices[3].y()));
754  gtrapElement->setAttributeNode(NewAttribute("v5x", vertices[4].x()));
755  gtrapElement->setAttributeNode(NewAttribute("v5y", vertices[4].y()));
756  gtrapElement->setAttributeNode(NewAttribute("v6x", vertices[5].x()));
757  gtrapElement->setAttributeNode(NewAttribute("v6y", vertices[5].y()));
758  gtrapElement->setAttributeNode(NewAttribute("v7x", vertices[6].x()));
759  gtrapElement->setAttributeNode(NewAttribute("v7y", vertices[6].y()));
760  gtrapElement->setAttributeNode(NewAttribute("v8x", vertices[7].x()));
761  gtrapElement->setAttributeNode(NewAttribute("v8y", vertices[7].y()));
762  gtrapElement->setAttributeNode(NewAttribute("lunit","mm"));
763  solElement->appendChild(gtrapElement);
764 }
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
const std::vector< G4TwoVector > & GetVertices() const
Double_t y
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double mm
Definition: G4SIunits.hh:114
G4double GetZHalfLength() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ HypeWrite()

void G4GDMLWriteSolids::HypeWrite ( xercesc::DOMElement *  solElement,
const G4Hype * const  hype 
)
protected

Definition at line 407 of file G4GDMLWriteSolids.cc.

408 {
409  const G4String& name = GenerateName(hype->GetName(),hype);
410 
411  xercesc::DOMElement* hypeElement = NewElement("hype");
412  hypeElement->setAttributeNode(NewAttribute("name",name));
413  hypeElement->setAttributeNode(NewAttribute("rmin",
414  hype->GetInnerRadius()/mm));
415  hypeElement->setAttributeNode(NewAttribute("rmax",
416  hype->GetOuterRadius()/mm));
417  hypeElement->setAttributeNode(NewAttribute("inst",
418  hype->GetInnerStereo()/degree));
419  hypeElement->setAttributeNode(NewAttribute("outst",
420  hype->GetOuterStereo()/degree));
421  hypeElement->setAttributeNode(NewAttribute("z",
422  2.0*hype->GetZHalfLength()/mm));
423  hypeElement->setAttributeNode(NewAttribute("aunit","deg"));
424  hypeElement->setAttributeNode(NewAttribute("lunit","mm"));
425  solElement->appendChild(hypeElement);
426 }
G4double GetZHalfLength() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetInnerRadius() const
G4double GetInnerStereo() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetOuterStereo() const
static const double degree
Definition: G4SIunits.hh:143
static const double mm
Definition: G4SIunits.hh:114
G4double GetOuterRadius() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MultiUnionWrite()

void G4GDMLWriteSolids::MultiUnionWrite ( xercesc::DOMElement *  solElement,
const G4MultiUnion * const   
)
protected

Definition at line 83 of file G4GDMLWriteSolids.cc.

85 {
86  G4Exception("G4GDMLWriteSolids::MultiUnionWrite()",
87  "InvalidSetup", FatalException,
88  "Installation with USolids primitives required!");
89  return;
90 }
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:

◆ OpticalSurfaceWrite()

void G4GDMLWriteSolids::OpticalSurfaceWrite ( xercesc::DOMElement *  solElement,
const G4OpticalSurface * const  surf 
)
protected

Definition at line 1006 of file G4GDMLWriteSolids.cc.

1008 {
1009  xercesc::DOMElement* optElement = NewElement("opticalsurface");
1010  G4OpticalSurfaceModel smodel = surf->GetModel();
1011  G4double sval = (smodel==glisur) ? surf->GetPolish() : surf->GetSigmaAlpha();
1012 
1013  optElement->setAttributeNode(NewAttribute("name", surf->GetName()));
1014  optElement->setAttributeNode(NewAttribute("model", smodel));
1015  optElement->setAttributeNode(NewAttribute("finish", surf->GetFinish()));
1016  optElement->setAttributeNode(NewAttribute("type", surf->GetType()));
1017  optElement->setAttributeNode(NewAttribute("value", sval));
1018 
1019  solElement->appendChild(optElement);
1020 }
const G4SurfaceType & GetType() const
const G4String & GetName() const
G4OpticalSurfaceModel GetModel() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetPolish() const
G4double GetSigmaAlpha() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4OpticalSurfaceFinish GetFinish() const
G4OpticalSurfaceModel
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ OrbWrite()

void G4GDMLWriteSolids::OrbWrite ( xercesc::DOMElement *  solElement,
const G4Orb * const  orb 
)
protected

Definition at line 429 of file G4GDMLWriteSolids.cc.

430 {
431  const G4String& name = GenerateName(orb->GetName(),orb);
432 
433  xercesc::DOMElement* orbElement = NewElement("orb");
434  orbElement->setAttributeNode(NewAttribute("name",name));
435  orbElement->setAttributeNode(NewAttribute("r",orb->GetRadius()/mm));
436  orbElement->setAttributeNode(NewAttribute("lunit","mm"));
437  solElement->appendChild(orbElement);
438 }
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetRadius() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ParaboloidWrite()

void G4GDMLWriteSolids::ParaboloidWrite ( xercesc::DOMElement *  solElement,
const G4Paraboloid * const  paraboloid 
)
protected

Definition at line 467 of file G4GDMLWriteSolids.cc.

469 {
470  const G4String& name = GenerateName(paraboloid->GetName(),paraboloid);
471 
472  xercesc::DOMElement* paraboloidElement = NewElement("paraboloid");
473  paraboloidElement->setAttributeNode(NewAttribute("name",name));
474  paraboloidElement->setAttributeNode(NewAttribute("rlo",
475  paraboloid->GetRadiusMinusZ()/mm));
476  paraboloidElement->setAttributeNode(NewAttribute("rhi",
477  paraboloid->GetRadiusPlusZ()/mm));
478  paraboloidElement->setAttributeNode(NewAttribute("dz",
479  paraboloid->GetZHalfLength()/mm));
480  paraboloidElement->setAttributeNode(NewAttribute("lunit","mm"));
481  solElement->appendChild(paraboloidElement);
482 }
G4String name
Definition: TRTMaterials.hh:40
G4double GetRadiusPlusZ() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetRadiusMinusZ() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double mm
Definition: G4SIunits.hh:114
G4double GetZHalfLength() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ParaWrite()

void G4GDMLWriteSolids::ParaWrite ( xercesc::DOMElement *  solElement,
const G4Para * const  para 
)
protected

Definition at line 441 of file G4GDMLWriteSolids.cc.

442 {
443  const G4String& name = GenerateName(para->GetName(),para);
444 
445  const G4ThreeVector simaxis = para->GetSymAxis();
446  const G4double alpha = std::atan(para->GetTanAlpha());
447  const G4double phi = simaxis.phi();
448  const G4double theta = simaxis.theta();
449 
450  xercesc::DOMElement* paraElement = NewElement("para");
451  paraElement->setAttributeNode(NewAttribute("name",name));
452  paraElement->setAttributeNode(NewAttribute("x",
453  2.0*para->GetXHalfLength()/mm));
454  paraElement->setAttributeNode(NewAttribute("y",
455  2.0*para->GetYHalfLength()/mm));
456  paraElement->setAttributeNode(NewAttribute("z",
457  2.0*para->GetZHalfLength()/mm));
458  paraElement->setAttributeNode(NewAttribute("alpha",alpha/degree));
459  paraElement->setAttributeNode(NewAttribute("theta",theta/degree));
460  paraElement->setAttributeNode(NewAttribute("phi",phi/degree));
461  paraElement->setAttributeNode(NewAttribute("aunit","deg"));
462  paraElement->setAttributeNode(NewAttribute("lunit","mm"));
463  solElement->appendChild(paraElement);
464 }
G4double GetZHalfLength() const
G4String name
Definition: TRTMaterials.hh:40
G4ThreeVector GetSymAxis() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetTanAlpha() const
G4double GetXHalfLength() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double degree
Definition: G4SIunits.hh:143
G4double GetYHalfLength() const
double G4double
Definition: G4Types.hh:76
static const G4double alpha
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PolyconeWrite()

void G4GDMLWriteSolids::PolyconeWrite ( xercesc::DOMElement *  solElement,
const G4Polycone * const  polycone 
)
protected

Definition at line 484 of file G4GDMLWriteSolids.cc.

486 {
487  const G4String& name = GenerateName(polycone->GetName(),polycone);
488 
489  xercesc::DOMElement* polyconeElement = NewElement("polycone");
490  polyconeElement->setAttributeNode(NewAttribute("name",name));
491  polyconeElement->setAttributeNode(NewAttribute("startphi",
493  polyconeElement->setAttributeNode(NewAttribute("deltaphi",
495  polyconeElement->setAttributeNode(NewAttribute("aunit","deg"));
496  polyconeElement->setAttributeNode(NewAttribute("lunit","mm"));
497  solElement->appendChild(polyconeElement);
498 
499  const size_t num_zplanes = polycone->GetOriginalParameters()->Num_z_planes;
500  const G4double* z_array = polycone->GetOriginalParameters()->Z_values;
501  const G4double* rmin_array = polycone->GetOriginalParameters()->Rmin;
502  const G4double* rmax_array = polycone->GetOriginalParameters()->Rmax;
503 
504  for (size_t i=0; i<num_zplanes; i++)
505  {
506  ZplaneWrite(polyconeElement,z_array[i],rmin_array[i],rmax_array[i]);
507  }
508 
509 
510 }
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4PolyconeHistorical * GetOriginalParameters() const
void ZplaneWrite(xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double degree
Definition: G4SIunits.hh:143
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PolyhedraWrite()

void G4GDMLWriteSolids::PolyhedraWrite ( xercesc::DOMElement *  solElement,
const G4Polyhedra * const  polyhedra 
)
protected

Definition at line 539 of file G4GDMLWriteSolids.cc.

541 {
542  const G4String& name = GenerateName(polyhedra->GetName(),polyhedra);
543  if(polyhedra->IsGeneric() == false){
544  xercesc::DOMElement* polyhedraElement = NewElement("polyhedra");
545  polyhedraElement->setAttributeNode(NewAttribute("name",name));
546  polyhedraElement->setAttributeNode(NewAttribute("startphi",
547  polyhedra->GetOriginalParameters()->Start_angle/degree));
548  polyhedraElement->setAttributeNode(NewAttribute("deltaphi",
550  polyhedraElement->setAttributeNode(NewAttribute("numsides",
551  polyhedra->GetOriginalParameters()->numSide));
552  polyhedraElement->setAttributeNode(NewAttribute("aunit","deg"));
553  polyhedraElement->setAttributeNode(NewAttribute("lunit","mm"));
554  solElement->appendChild(polyhedraElement);
555 
556  const size_t num_zplanes = polyhedra->GetOriginalParameters()->Num_z_planes;
557  const G4double* z_array = polyhedra->GetOriginalParameters()->Z_values;
558  const G4double* rmin_array = polyhedra->GetOriginalParameters()->Rmin;
559  const G4double* rmax_array = polyhedra->GetOriginalParameters()->Rmax;
560 
561  const G4double convertRad =
562  std::cos(0.5*polyhedra->GetOriginalParameters()->Opening_angle
563  / polyhedra->GetOriginalParameters()->numSide);
564 
565  for (size_t i=0;i<num_zplanes;i++)
566  {
567  ZplaneWrite(polyhedraElement,z_array[i],
568  rmin_array[i]*convertRad, rmax_array[i]*convertRad);
569  }
570  }else{//generic polyhedra
571 
572  xercesc::DOMElement* polyhedraElement = NewElement("genericPolyhedra");
573  polyhedraElement->setAttributeNode(NewAttribute("name",name));
574  polyhedraElement->setAttributeNode(NewAttribute("startphi",
575  polyhedra->GetOriginalParameters()->Start_angle/degree));
576  polyhedraElement->setAttributeNode(NewAttribute("deltaphi",
578  polyhedraElement->setAttributeNode(NewAttribute("numsides",
579  polyhedra->GetOriginalParameters()->numSide));
580  polyhedraElement->setAttributeNode(NewAttribute("aunit","deg"));
581  polyhedraElement->setAttributeNode(NewAttribute("lunit","mm"));
582  solElement->appendChild(polyhedraElement);
583 
584  const size_t num_rzpoints = polyhedra->GetNumRZCorner();
585 
586  for (size_t i=0;i<num_rzpoints;i++)
587  {
588  const G4double r_point = polyhedra->GetCorner(i).r;
589  const G4double z_point = polyhedra->GetCorner(i).z;
590  RZPointWrite(polyhedraElement,r_point,z_point);
591  }
592  }
593 }
G4bool IsGeneric() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
void RZPointWrite(xercesc::DOMElement *, const G4double &, const G4double &)
G4PolyhedraHistorical * GetOriginalParameters() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void ZplaneWrite(xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4PolyhedraSideRZ GetCorner(const G4int index) const
static const double degree
Definition: G4SIunits.hh:143
double G4double
Definition: G4Types.hh:76
G4int GetNumRZCorner() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ RZPointWrite()

void G4GDMLWriteSolids::RZPointWrite ( xercesc::DOMElement *  element,
const G4double r,
const G4double z 
)
protected

Definition at line 996 of file G4GDMLWriteSolids.cc.

998 {
999  xercesc::DOMElement* rzpointElement = NewElement("rzpoint");
1000  rzpointElement->setAttributeNode(NewAttribute("r",r/mm));
1001  rzpointElement->setAttributeNode(NewAttribute("z",z/mm));
1002  element->appendChild(rzpointElement);
1003 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SolidsWrite()

void G4GDMLWriteSolids::SolidsWrite ( xercesc::DOMElement *  gdmlElement)
virtual

Implements G4GDMLWrite.

Definition at line 1022 of file G4GDMLWriteSolids.cc.

1023 {
1024  G4cout << "G4GDML: Writing solids..." << G4endl;
1025 
1026  solidsElement = NewElement("solids");
1027  gdmlElement->appendChild(solidsElement);
1028 
1029  solidList.clear();
1030 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
xercesc::DOMElement * solidsElement
std::vector< const G4VSolid * > solidList
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SphereWrite()

void G4GDMLWriteSolids::SphereWrite ( xercesc::DOMElement *  solElement,
const G4Sphere * const  sphere 
)
protected

Definition at line 596 of file G4GDMLWriteSolids.cc.

597 {
598  const G4String& name = GenerateName(sphere->GetName(),sphere);
599 
600  xercesc::DOMElement* sphereElement = NewElement("sphere");
601  sphereElement->setAttributeNode(NewAttribute("name",name));
602  sphereElement->setAttributeNode(NewAttribute("rmin",
603  sphere->GetInnerRadius()/mm));
604  sphereElement->setAttributeNode(NewAttribute("rmax",
605  sphere->GetOuterRadius()/mm));
606  sphereElement->setAttributeNode(NewAttribute("startphi",
607  sphere->GetStartPhiAngle()/degree));
608  sphereElement->setAttributeNode(NewAttribute("deltaphi",
609  sphere->GetDeltaPhiAngle()/degree));
610  sphereElement->setAttributeNode(NewAttribute("starttheta",
611  sphere->GetStartThetaAngle()/degree));
612  sphereElement->setAttributeNode(NewAttribute("deltatheta",
613  sphere->GetDeltaThetaAngle()/degree));
614  sphereElement->setAttributeNode(NewAttribute("aunit","deg"));
615  sphereElement->setAttributeNode(NewAttribute("lunit","mm"));
616  solElement->appendChild(sphereElement);
617 }
G4double GetStartThetaAngle() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetStartPhiAngle() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetOuterRadius() const
G4double GetDeltaPhiAngle() const
G4double GetDeltaThetaAngle() const
static const double degree
Definition: G4SIunits.hh:143
G4double GetInnerRadius() const
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TessellatedWrite()

void G4GDMLWriteSolids::TessellatedWrite ( xercesc::DOMElement *  solElement,
const G4TessellatedSolid * const  tessellated 
)
protected

Definition at line 620 of file G4GDMLWriteSolids.cc.

622 {
623  const G4String& solid_name = tessellated->GetName();
624  const G4String& name = GenerateName(solid_name, tessellated);
625 
626  xercesc::DOMElement* tessellatedElement = NewElement("tessellated");
627  tessellatedElement->setAttributeNode(NewAttribute("name",name));
628  tessellatedElement->setAttributeNode(NewAttribute("aunit","deg"));
629  tessellatedElement->setAttributeNode(NewAttribute("lunit","mm"));
630  solElement->appendChild(tessellatedElement);
631 
632  std::map<G4ThreeVector, G4String, G4ThreeVectorCompare> vertexMap;
633 
634  const size_t NumFacets = tessellated->GetNumberOfFacets();
635  size_t NumVertex = 0;
636 
637  for (size_t i=0;i<NumFacets;i++)
638  {
639  const G4VFacet* facet = tessellated->GetFacet(i);
640  const size_t NumVertexPerFacet = facet->GetNumberOfVertices();
641 
642  G4String FacetTag;
643 
644  if (NumVertexPerFacet==3) { FacetTag="triangular"; } else
645  if (NumVertexPerFacet==4) { FacetTag="quadrangular"; }
646  else
647  {
648  G4Exception("G4GDMLWriteSolids::TessellatedWrite()", "InvalidSetup",
649  FatalException, "Facet should contain 3 or 4 vertices!");
650  }
651 
652  xercesc::DOMElement* facetElement = NewElement(FacetTag);
653  tessellatedElement->appendChild(facetElement);
654 
655  for (size_t j=0; j<NumVertexPerFacet; j++)
656  {
657  std::stringstream name_stream;
658  std::stringstream ref_stream;
659 
660  name_stream << "vertex" << (j+1);
661  ref_stream << solid_name << "_v" << NumVertex;
662 
663  const G4String& fname = name_stream.str(); // facet's tag variable
664  G4String ref = ref_stream.str(); // vertex tag to be associated
665 
666  // Now search for the existance of the current vertex in the
667  // map of cached vertices. If existing, do NOT store it as
668  // position in the GDML file, so avoiding duplication; otherwise
669  // cache it in the local map and add it as position in the
670  // "define" section of the GDML file.
671 
672  const G4ThreeVector& vertex = facet->GetVertex(j);
673 
674  if(vertexMap.find(vertex) != vertexMap.end()) // Vertex is cached
675  {
676  ref = vertexMap[vertex]; // Set the proper tag for it
677  }
678  else // Vertex not found
679  {
680  vertexMap.insert(std::make_pair(vertex,ref)); // Cache vertex and ...
681  AddPosition(ref, vertex); // ... add it to define section!
682  NumVertex++;
683  }
684 
685  // Now create association of the vertex with its facet
686  //
687  facetElement->setAttributeNode(NewAttribute(fname,ref));
688  }
689  }
690 }
G4String name
Definition: TRTMaterials.hh:40
virtual G4int GetNumberOfVertices() const =0
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4VFacet * GetFacet(G4int i) const
G4String GetName() const
G4int GetNumberOfFacets() const
void AddPosition(const G4String &name, const G4ThreeVector &pos)
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
virtual G4ThreeVector GetVertex(G4int i) const =0
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TetWrite()

void G4GDMLWriteSolids::TetWrite ( xercesc::DOMElement *  solElement,
const G4Tet * const  tet 
)
protected

Definition at line 693 of file G4GDMLWriteSolids.cc.

694 {
695  const G4String& solid_name = tet->GetName();
696  const G4String& name = GenerateName(solid_name, tet);
697 
698  std::vector<G4ThreeVector> vertexList = tet->GetVertices();
699 
700  xercesc::DOMElement* tetElement = NewElement("tet");
701  tetElement->setAttributeNode(NewAttribute("name",name));
702  tetElement->setAttributeNode(NewAttribute("vertex1",solid_name+"_v1"));
703  tetElement->setAttributeNode(NewAttribute("vertex2",solid_name+"_v2"));
704  tetElement->setAttributeNode(NewAttribute("vertex3",solid_name+"_v3"));
705  tetElement->setAttributeNode(NewAttribute("vertex4",solid_name+"_v4"));
706  tetElement->setAttributeNode(NewAttribute("lunit","mm"));
707  solElement->appendChild(tetElement);
708 
709  AddPosition(solid_name+"_v1",vertexList[0]);
710  AddPosition(solid_name+"_v2",vertexList[1]);
711  AddPosition(solid_name+"_v3",vertexList[2]);
712  AddPosition(solid_name+"_v4",vertexList[3]);
713 }
G4String name
Definition: TRTMaterials.hh:40
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:790
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
void AddPosition(const G4String &name, const G4ThreeVector &pos)
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TorusWrite()

void G4GDMLWriteSolids::TorusWrite ( xercesc::DOMElement *  solElement,
const G4Torus * const  torus 
)
protected

Definition at line 716 of file G4GDMLWriteSolids.cc.

717 {
718  const G4String& name = GenerateName(torus->GetName(),torus);
719 
720  xercesc::DOMElement* torusElement = NewElement("torus");
721  torusElement->setAttributeNode(NewAttribute("name",name));
722  torusElement->setAttributeNode(NewAttribute("rmin",torus->GetRmin()/mm));
723  torusElement->setAttributeNode(NewAttribute("rmax",torus->GetRmax()/mm));
724  torusElement->setAttributeNode(NewAttribute("rtor",torus->GetRtor()/mm));
725  torusElement->
726  setAttributeNode(NewAttribute("startphi",torus->GetSPhi()/degree));
727  torusElement->
728  setAttributeNode(NewAttribute("deltaphi",torus->GetDPhi()/degree));
729  torusElement->setAttributeNode(NewAttribute("aunit","deg"));
730  torusElement->setAttributeNode(NewAttribute("lunit","mm"));
731  solElement->appendChild(torusElement);
732 }
G4String name
Definition: TRTMaterials.hh:40
G4double GetRmax() const
G4double GetRtor() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetDPhi() const
G4double GetRmin() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double degree
Definition: G4SIunits.hh:143
G4double GetSPhi() const
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TrapWrite()

void G4GDMLWriteSolids::TrapWrite ( xercesc::DOMElement *  solElement,
const G4Trap * const  trap 
)
protected

Definition at line 767 of file G4GDMLWriteSolids.cc.

768 {
769  const G4String& name = GenerateName(trap->GetName(),trap);
770 
771  const G4ThreeVector& simaxis = trap->GetSymAxis();
772  const G4double phi = simaxis.phi();
773  const G4double theta = simaxis.theta();
774  const G4double alpha1 = std::atan(trap->GetTanAlpha1());
775  const G4double alpha2 = std::atan(trap->GetTanAlpha2());
776 
777  xercesc::DOMElement* trapElement = NewElement("trap");
778  trapElement->setAttributeNode(NewAttribute("name",name));
779  trapElement->setAttributeNode(NewAttribute("z",
780  2.0*trap->GetZHalfLength()/mm));
781  trapElement->setAttributeNode(NewAttribute("theta",theta/degree));
782  trapElement->setAttributeNode(NewAttribute("phi",phi/degree));
783  trapElement->setAttributeNode(NewAttribute("y1",
784  2.0*trap->GetYHalfLength1()/mm));
785  trapElement->setAttributeNode(NewAttribute("x1",
786  2.0*trap->GetXHalfLength1()/mm));
787  trapElement->setAttributeNode(NewAttribute("x2",
788  2.0*trap->GetXHalfLength2()/mm));
789  trapElement->setAttributeNode(NewAttribute("alpha1",alpha1/degree));
790  trapElement->setAttributeNode(NewAttribute("y2",
791  2.0*trap->GetYHalfLength2()/mm));
792  trapElement->setAttributeNode(NewAttribute("x3",
793  2.0*trap->GetXHalfLength3()/mm));
794  trapElement->setAttributeNode(NewAttribute("x4",
795  2.0*trap->GetXHalfLength4()/mm));
796  trapElement->setAttributeNode(NewAttribute("alpha2",alpha2/degree));
797  trapElement->setAttributeNode(NewAttribute("aunit","deg"));
798  trapElement->setAttributeNode(NewAttribute("lunit","mm"));
799  solElement->appendChild(trapElement);
800 }
double phi() const
G4ThreeVector GetSymAxis() const
G4String name
Definition: TRTMaterials.hh:40
G4double GetXHalfLength3() const
G4double GetXHalfLength4() const
G4double GetYHalfLength2() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetYHalfLength1() const
G4double GetXHalfLength2() const
G4double GetZHalfLength() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetXHalfLength1() const
G4double GetTanAlpha2() const
static const double degree
Definition: G4SIunits.hh:143
double G4double
Definition: G4Types.hh:76
static const double mm
Definition: G4SIunits.hh:114
G4double GetTanAlpha1() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TrdWrite()

void G4GDMLWriteSolids::TrdWrite ( xercesc::DOMElement *  solElement,
const G4Trd * const  trd 
)
protected

Definition at line 803 of file G4GDMLWriteSolids.cc.

804 {
805  const G4String& name = GenerateName(trd->GetName(),trd);
806 
807  xercesc::DOMElement* trdElement = NewElement("trd");
808  trdElement->setAttributeNode(NewAttribute("name",name));
809  trdElement->setAttributeNode(NewAttribute("x1",
810  2.0*trd->GetXHalfLength1()/mm));
811  trdElement->setAttributeNode(NewAttribute("x2",
812  2.0*trd->GetXHalfLength2()/mm));
813  trdElement->setAttributeNode(NewAttribute("y1",
814  2.0*trd->GetYHalfLength1()/mm));
815  trdElement->setAttributeNode(NewAttribute("y2",
816  2.0*trd->GetYHalfLength2()/mm));
817  trdElement->setAttributeNode(NewAttribute("z",
818  2.0*trd->GetZHalfLength()/mm));
819  trdElement->setAttributeNode(NewAttribute("lunit","mm"));
820  solElement->appendChild(trdElement);
821 }
G4double GetXHalfLength1() const
G4double GetYHalfLength1() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetXHalfLength2() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetYHalfLength2() const
static const double mm
Definition: G4SIunits.hh:114
G4double GetZHalfLength() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TubeWrite()

void G4GDMLWriteSolids::TubeWrite ( xercesc::DOMElement *  solElement,
const G4Tubs * const  tube 
)
protected

Definition at line 824 of file G4GDMLWriteSolids.cc.

825 {
826  const G4String& name = GenerateName(tube->GetName(),tube);
827 
828  xercesc::DOMElement* tubeElement = NewElement("tube");
829  tubeElement->setAttributeNode(NewAttribute("name",name));
830  tubeElement->setAttributeNode(NewAttribute("rmin",
831  tube->GetInnerRadius()/mm));
832  tubeElement->setAttributeNode(NewAttribute("rmax",
833  tube->GetOuterRadius()/mm));
834  tubeElement->setAttributeNode(NewAttribute("z",
835  2.0*tube->GetZHalfLength()/mm));
836  tubeElement->setAttributeNode(NewAttribute("startphi",
837  tube->GetStartPhiAngle()/degree));
838  tubeElement->setAttributeNode(NewAttribute("deltaphi",
839  tube->GetDeltaPhiAngle()/degree));
840  tubeElement->setAttributeNode(NewAttribute("aunit","deg"));
841  tubeElement->setAttributeNode(NewAttribute("lunit","mm"));
842  solElement->appendChild(tubeElement);
843 }
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetZHalfLength() const
G4double GetOuterRadius() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double degree
Definition: G4SIunits.hh:143
G4double GetInnerRadius() const
G4double GetStartPhiAngle() const
static const double mm
Definition: G4SIunits.hh:114
G4double GetDeltaPhiAngle() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TwistedboxWrite()

void G4GDMLWriteSolids::TwistedboxWrite ( xercesc::DOMElement *  solElement,
const G4TwistedBox * const  twistedbox 
)
protected

Definition at line 880 of file G4GDMLWriteSolids.cc.

882 {
883  const G4String& name = GenerateName(twistedbox->GetName(),twistedbox);
884 
885  xercesc::DOMElement* twistedboxElement = NewElement("twistedbox");
886  twistedboxElement->setAttributeNode(NewAttribute("name",name));
887  twistedboxElement->setAttributeNode(NewAttribute("x",
888  2.0*twistedbox->GetXHalfLength()/mm));
889  twistedboxElement->setAttributeNode(NewAttribute("y",
890  2.0*twistedbox->GetYHalfLength()/mm));
891  twistedboxElement->setAttributeNode(NewAttribute("z",
892  2.0*twistedbox->GetZHalfLength()/mm));
893  twistedboxElement->setAttributeNode(NewAttribute("PhiTwist",
894  twistedbox->GetPhiTwist()/degree));
895  twistedboxElement->setAttributeNode(NewAttribute("aunit","deg"));
896  twistedboxElement->setAttributeNode(NewAttribute("lunit","mm"));
897  solElement->appendChild(twistedboxElement);
898 }
G4String name
Definition: TRTMaterials.hh:40
G4double GetXHalfLength() const
Definition: G4TwistedBox.hh:73
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetYHalfLength() const
Definition: G4TwistedBox.hh:74
G4double GetPhiTwist() const
Definition: G4TwistedBox.hh:76
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetZHalfLength() const
Definition: G4TwistedBox.hh:75
static const double degree
Definition: G4SIunits.hh:143
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TwistedtrapWrite()

void G4GDMLWriteSolids::TwistedtrapWrite ( xercesc::DOMElement *  solElement,
const G4TwistedTrap * const  twistedtrap 
)
protected

Definition at line 901 of file G4GDMLWriteSolids.cc.

903 {
904  const G4String& name = GenerateName(twistedtrap->GetName(),twistedtrap);
905 
906  xercesc::DOMElement* twistedtrapElement = NewElement("twistedtrap");
907  twistedtrapElement->setAttributeNode(NewAttribute("name",name));
908  twistedtrapElement->setAttributeNode(NewAttribute("y1",
909  2.0*twistedtrap->GetY1HalfLength()/mm));
910  twistedtrapElement->setAttributeNode(NewAttribute("x1",
911  2.0*twistedtrap->GetX1HalfLength()/mm));
912  twistedtrapElement->setAttributeNode(NewAttribute("x2",
913  2.0*twistedtrap->GetX2HalfLength()/mm));
914  twistedtrapElement->setAttributeNode(NewAttribute("y2",
915  2.0*twistedtrap->GetY2HalfLength()/mm));
916  twistedtrapElement->setAttributeNode(NewAttribute("x3",
917  2.0*twistedtrap->GetX3HalfLength()/mm));
918  twistedtrapElement->setAttributeNode(NewAttribute("x4",
919  2.0*twistedtrap->GetX4HalfLength()/mm));
920  twistedtrapElement->setAttributeNode(NewAttribute("z",
921  2.0*twistedtrap->GetZHalfLength()/mm));
922  twistedtrapElement->setAttributeNode(NewAttribute("Alph",
923  twistedtrap->GetTiltAngleAlpha()/degree));
924  twistedtrapElement->setAttributeNode(NewAttribute("Theta",
925  twistedtrap->GetPolarAngleTheta()/degree));
926  twistedtrapElement->setAttributeNode(NewAttribute("Phi",
927  twistedtrap->GetAzimuthalAnglePhi()/degree));
928  twistedtrapElement->setAttributeNode(NewAttribute("PhiTwist",
929  twistedtrap->GetPhiTwist()/degree));
930  twistedtrapElement->setAttributeNode(NewAttribute("aunit","deg"));
931  twistedtrapElement->setAttributeNode(NewAttribute("lunit","mm"));
932 
933  solElement->appendChild(twistedtrapElement);
934 }
G4double GetPolarAngleTheta() const
G4double GetX2HalfLength() const
G4double GetY2HalfLength() const
G4String name
Definition: TRTMaterials.hh:40
G4double GetX1HalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetZHalfLength() const
G4String GetName() const
G4double GetY1HalfLength() const
G4double GetTiltAngleAlpha() const
G4double GetPhiTwist() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetX3HalfLength() const
static const double degree
Definition: G4SIunits.hh:143
G4double GetX4HalfLength() const
static const double mm
Definition: G4SIunits.hh:114
G4double GetAzimuthalAnglePhi() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TwistedtrdWrite()

void G4GDMLWriteSolids::TwistedtrdWrite ( xercesc::DOMElement *  solElement,
const G4TwistedTrd * const  twistedtrd 
)
protected

Definition at line 937 of file G4GDMLWriteSolids.cc.

939 {
940  const G4String& name = GenerateName(twistedtrd->GetName(),twistedtrd);
941 
942  xercesc::DOMElement* twistedtrdElement = NewElement("twistedtrd");
943  twistedtrdElement->setAttributeNode(NewAttribute("name",name));
944  twistedtrdElement->setAttributeNode(NewAttribute("x1",
945  2.0*twistedtrd->GetX1HalfLength()/mm));
946  twistedtrdElement->setAttributeNode(NewAttribute("x2",
947  2.0*twistedtrd->GetX2HalfLength()/mm));
948  twistedtrdElement->setAttributeNode(NewAttribute("y1",
949  2.0*twistedtrd->GetY1HalfLength()/mm));
950  twistedtrdElement->setAttributeNode(NewAttribute("y2",
951  2.0*twistedtrd->GetY2HalfLength()/mm));
952  twistedtrdElement->setAttributeNode(NewAttribute("z",
953  2.0*twistedtrd->GetZHalfLength()/mm));
954  twistedtrdElement->setAttributeNode(NewAttribute("PhiTwist",
955  twistedtrd->GetPhiTwist()/degree));
956  twistedtrdElement->setAttributeNode(NewAttribute("aunit","deg"));
957  twistedtrdElement->setAttributeNode(NewAttribute("lunit","mm"));
958  solElement->appendChild(twistedtrdElement);
959 }
G4double GetX2HalfLength() const
Definition: G4TwistedTrd.hh:78
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetY2HalfLength() const
Definition: G4TwistedTrd.hh:80
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetZHalfLength() const
Definition: G4TwistedTrd.hh:81
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetX1HalfLength() const
Definition: G4TwistedTrd.hh:77
static const double degree
Definition: G4SIunits.hh:143
G4double GetY1HalfLength() const
Definition: G4TwistedTrd.hh:79
static const double mm
Definition: G4SIunits.hh:114
G4double GetPhiTwist() const
Definition: G4TwistedTrd.hh:82
Here is the call graph for this function:
Here is the caller graph for this function:

◆ TwistedtubsWrite()

void G4GDMLWriteSolids::TwistedtubsWrite ( xercesc::DOMElement *  solElement,
const G4TwistedTubs * const  twistedtubs 
)
protected

Definition at line 962 of file G4GDMLWriteSolids.cc.

964 {
965  const G4String& name = GenerateName(twistedtubs->GetName(),twistedtubs);
966 
967  xercesc::DOMElement* twistedtubsElement = NewElement("twistedtubs");
968  twistedtubsElement->setAttributeNode(NewAttribute("name",name));
969  twistedtubsElement->setAttributeNode(NewAttribute("twistedangle",
970  twistedtubs->GetPhiTwist()/degree));
971  twistedtubsElement->setAttributeNode(NewAttribute("endinnerrad",
972  twistedtubs->GetInnerRadius()/mm));
973  twistedtubsElement->setAttributeNode(NewAttribute("endouterrad",
974  twistedtubs->GetOuterRadius()/mm));
975  twistedtubsElement->setAttributeNode(NewAttribute("zlen",
976  2.0*twistedtubs->GetZHalfLength()/mm));
977  twistedtubsElement->setAttributeNode(NewAttribute("phi",
978  twistedtubs->GetDPhi()/degree));
979  twistedtubsElement->setAttributeNode(NewAttribute("aunit","deg"));
980  twistedtubsElement->setAttributeNode(NewAttribute("lunit","mm"));
981  solElement->appendChild(twistedtubsElement);
982 }
G4double GetOuterRadius() const
G4String name
Definition: TRTMaterials.hh:40
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GetName() const
G4double GetZHalfLength() const
G4double GetPhiTwist() const
G4double GetDPhi() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetInnerRadius() const
static const double degree
Definition: G4SIunits.hh:143
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ XtruWrite()

void G4GDMLWriteSolids::XtruWrite ( xercesc::DOMElement *  solElement,
const G4ExtrudedSolid * const  xtru 
)
protected

Definition at line 363 of file G4GDMLWriteSolids.cc.

365 {
366  const G4String& name = GenerateName(xtru->GetName(),xtru);
367 
368  xercesc::DOMElement* xtruElement = NewElement("xtru");
369  xtruElement->setAttributeNode(NewAttribute("name",name));
370  xtruElement->setAttributeNode(NewAttribute("lunit","mm"));
371  solElement->appendChild(xtruElement);
372 
373  const G4int NumVertex = xtru->GetNofVertices();
374 
375  for (G4int i=0;i<NumVertex;i++)
376  {
377  xercesc::DOMElement* twoDimVertexElement = NewElement("twoDimVertex");
378  xtruElement->appendChild(twoDimVertexElement);
379 
380  const G4TwoVector& vertex = xtru->GetVertex(i);
381 
382  twoDimVertexElement->setAttributeNode(NewAttribute("x",vertex.x()/mm));
383  twoDimVertexElement->setAttributeNode(NewAttribute("y",vertex.y()/mm));
384  }
385 
386  const G4int NumSection = xtru->GetNofZSections();
387 
388  for (G4int i=0;i<NumSection;i++)
389  {
390  xercesc::DOMElement* sectionElement = NewElement("section");
391  xtruElement->appendChild(sectionElement);
392 
393  const G4ExtrudedSolid::ZSection section = xtru->GetZSection(i);
394 
395  sectionElement->setAttributeNode(NewAttribute("zOrder",i));
396  sectionElement->setAttributeNode(NewAttribute("zPosition",section.fZ/mm));
397  sectionElement->
398  setAttributeNode(NewAttribute("xOffset",section.fOffset.x()/mm));
399  sectionElement->
400  setAttributeNode(NewAttribute("yOffset",section.fOffset.y()/mm));
401  sectionElement->
402  setAttributeNode(NewAttribute("scalingFactor",section.fScale));
403  }
404 }
G4int GetNofZSections() const
double y() const
G4String name
Definition: TRTMaterials.hh:40
ZSection GetZSection(G4int index) const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
int G4int
Definition: G4Types.hh:78
G4String GetName() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4int GetNofVertices() const
G4TwoVector GetVertex(G4int index) const
double x() const
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ZplaneWrite()

void G4GDMLWriteSolids::ZplaneWrite ( xercesc::DOMElement *  element,
const G4double z,
const G4double rmin,
const G4double rmax 
)
protected

Definition at line 985 of file G4GDMLWriteSolids.cc.

987 {
988  xercesc::DOMElement* zplaneElement = NewElement("zplane");
989  zplaneElement->setAttributeNode(NewAttribute("z",z/mm));
990  zplaneElement->setAttributeNode(NewAttribute("rmin",rmin/mm));
991  zplaneElement->setAttributeNode(NewAttribute("rmax",rmax/mm));
992  element->appendChild(zplaneElement);
993 }
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ maxTransforms

const G4int G4GDMLWriteSolids::maxTransforms = 8
staticprotected

Definition at line 147 of file G4GDMLWriteSolids.hh.

◆ solidList

std::vector<const G4VSolid*> G4GDMLWriteSolids::solidList
protected

Definition at line 145 of file G4GDMLWriteSolids.hh.

◆ solidsElement

xercesc::DOMElement* G4GDMLWriteSolids::solidsElement
protected

Definition at line 146 of file G4GDMLWriteSolids.hh.


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