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

#include <G4GDMLWriteSolids.hh>

Inheritance diagram for G4GDMLWriteSolids:
Collaboration diagram for G4GDMLWriteSolids:

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 ScaledWrite (xercesc::DOMElement *, const G4ScaledSolid *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 ()
 
VolumeMapType & VolumeMap ()
 
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
 
PhysVolumeMapType & PvolumeMap ()
 
DepthMapType & DepthMap ()
 

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 79 of file G4GDMLWriteSolids.hh.

Constructor & Destructor Documentation

G4GDMLWriteSolids::G4GDMLWriteSolids ( )
protected

Definition at line 73 of file G4GDMLWriteSolids.cc.

75 {
76 }
xercesc::DOMElement * solidsElement
G4GDMLWriteSolids::~G4GDMLWriteSolids ( )
protectedvirtual

Definition at line 78 of file G4GDMLWriteSolids.cc.

79 {
80 }

Member Function Documentation

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

Definition at line 1066 of file G4GDMLWriteSolids.cc.

1067 {
1068  for (size_t i=0; i<solidList.size(); i++) // Check if solid is
1069  { // already in the list!
1070  if (solidList[i] == solidPtr) { return; }
1071  }
1072 
1073  solidList.push_back(solidPtr);
1074 
1075  if (const G4BooleanSolid* const booleanPtr
1076  = dynamic_cast<const G4BooleanSolid*>(solidPtr))
1077  { BooleanWrite(solidsElement,booleanPtr); } else
1078  if (const G4ScaledSolid* const scaledPtr
1079  = dynamic_cast<const G4ScaledSolid*>(solidPtr))
1080  { ScaledWrite(solidsElement,scaledPtr); } else
1081  if (solidPtr->GetEntityType()=="G4MultiUnion")
1082  { const G4MultiUnion* const munionPtr
1083  = static_cast<const G4MultiUnion*>(solidPtr);
1084  MultiUnionWrite(solidsElement,munionPtr); } else
1085  if (solidPtr->GetEntityType()=="G4Box")
1086  { const G4Box* const boxPtr
1087  = static_cast<const G4Box*>(solidPtr);
1088  BoxWrite(solidsElement,boxPtr); } else
1089  if (solidPtr->GetEntityType()=="G4Cons")
1090  { const G4Cons* const conePtr
1091  = static_cast<const G4Cons*>(solidPtr);
1092  ConeWrite(solidsElement,conePtr); } else
1093  if (solidPtr->GetEntityType()=="G4EllipticalCone")
1094  { const G4EllipticalCone* const elconePtr
1095  = static_cast<const G4EllipticalCone*>(solidPtr);
1096  ElconeWrite(solidsElement,elconePtr); } else
1097  if (solidPtr->GetEntityType()=="G4Ellipsoid")
1098  { const G4Ellipsoid* const ellipsoidPtr
1099  = static_cast<const G4Ellipsoid*>(solidPtr);
1100  EllipsoidWrite(solidsElement,ellipsoidPtr); } else
1101  if (solidPtr->GetEntityType()=="G4EllipticalTube")
1102  { const G4EllipticalTube* const eltubePtr
1103  = static_cast<const G4EllipticalTube*>(solidPtr);
1104  EltubeWrite(solidsElement,eltubePtr); } else
1105  if (solidPtr->GetEntityType()=="G4ExtrudedSolid")
1106  { const G4ExtrudedSolid* const xtruPtr
1107  = static_cast<const G4ExtrudedSolid*>(solidPtr);
1108  XtruWrite(solidsElement,xtruPtr); } else
1109  if (solidPtr->GetEntityType()=="G4Hype")
1110  { const G4Hype* const hypePtr
1111  = static_cast<const G4Hype*>(solidPtr);
1112  HypeWrite(solidsElement,hypePtr); } else
1113  if (solidPtr->GetEntityType()=="G4Orb")
1114  { const G4Orb* const orbPtr
1115  = static_cast<const G4Orb*>(solidPtr);
1116  OrbWrite(solidsElement,orbPtr); } else
1117  if (solidPtr->GetEntityType()=="G4Para")
1118  { const G4Para* const paraPtr
1119  = static_cast<const G4Para*>(solidPtr);
1120  ParaWrite(solidsElement,paraPtr); } else
1121  if (solidPtr->GetEntityType()=="G4Paraboloid")
1122  { const G4Paraboloid* const paraboloidPtr
1123  = static_cast<const G4Paraboloid*>(solidPtr);
1124  ParaboloidWrite(solidsElement,paraboloidPtr); } else
1125  if (solidPtr->GetEntityType()=="G4Polycone")
1126  { const G4Polycone* const polyconePtr
1127  = static_cast<const G4Polycone*>(solidPtr);
1128  PolyconeWrite(solidsElement,polyconePtr); } else
1129  if (solidPtr->GetEntityType()=="G4GenericPolycone")
1130  { const G4GenericPolycone* const genpolyconePtr
1131  = static_cast<const G4GenericPolycone*>(solidPtr);
1132  GenericPolyconeWrite(solidsElement,genpolyconePtr); } else
1133  if (solidPtr->GetEntityType()=="G4Polyhedra")
1134  { const G4Polyhedra* const polyhedraPtr
1135  = static_cast<const G4Polyhedra*>(solidPtr);
1136  PolyhedraWrite(solidsElement,polyhedraPtr); } else
1137  if (solidPtr->GetEntityType()=="G4Sphere")
1138  { const G4Sphere* const spherePtr
1139  = static_cast<const G4Sphere*>(solidPtr);
1140  SphereWrite(solidsElement,spherePtr); } else
1141  if (solidPtr->GetEntityType()=="G4TessellatedSolid")
1142  { const G4TessellatedSolid* const tessellatedPtr
1143  = static_cast<const G4TessellatedSolid*>(solidPtr);
1144  TessellatedWrite(solidsElement,tessellatedPtr); } else
1145  if (solidPtr->GetEntityType()=="G4Tet")
1146  { const G4Tet* const tetPtr
1147  = static_cast<const G4Tet*>(solidPtr);
1148  TetWrite(solidsElement,tetPtr); } else
1149  if (solidPtr->GetEntityType()=="G4Torus")
1150  { const G4Torus* const torusPtr
1151  = static_cast<const G4Torus*>(solidPtr);
1152  TorusWrite(solidsElement,torusPtr); } else
1153  if (solidPtr->GetEntityType()=="G4GenericTrap")
1154  { const G4GenericTrap* const gtrapPtr
1155  = static_cast<const G4GenericTrap*>(solidPtr);
1156  GenTrapWrite(solidsElement,gtrapPtr); } else
1157  if (solidPtr->GetEntityType()=="G4Trap")
1158  { const G4Trap* const trapPtr
1159  = static_cast<const G4Trap*>(solidPtr);
1160  TrapWrite(solidsElement,trapPtr); } else
1161  if (solidPtr->GetEntityType()=="G4Trd")
1162  { const G4Trd* const trdPtr
1163  = static_cast<const G4Trd*>(solidPtr);
1164  TrdWrite(solidsElement,trdPtr); } else
1165  if (solidPtr->GetEntityType()=="G4Tubs")
1166  { const G4Tubs* const tubePtr
1167  = static_cast<const G4Tubs*>(solidPtr);
1168  TubeWrite(solidsElement,tubePtr); } else
1169  if (solidPtr->GetEntityType()=="G4CutTubs")
1170  { const G4CutTubs* const cuttubePtr
1171  = static_cast<const G4CutTubs*>(solidPtr);
1172  CutTubeWrite(solidsElement,cuttubePtr); } else
1173  if (solidPtr->GetEntityType()=="G4TwistedBox")
1174  { const G4TwistedBox* const twistedboxPtr
1175  = static_cast<const G4TwistedBox*>(solidPtr);
1176  TwistedboxWrite(solidsElement,twistedboxPtr); } else
1177  if (solidPtr->GetEntityType()=="G4TwistedTrap")
1178  { const G4TwistedTrap* const twistedtrapPtr
1179  = static_cast<const G4TwistedTrap*>(solidPtr);
1180  TwistedtrapWrite(solidsElement,twistedtrapPtr); } else
1181  if (solidPtr->GetEntityType()=="G4TwistedTrd")
1182  { const G4TwistedTrd* const twistedtrdPtr
1183  = static_cast<const G4TwistedTrd*>(solidPtr);
1184  TwistedtrdWrite(solidsElement,twistedtrdPtr); } else
1185  if (solidPtr->GetEntityType()=="G4TwistedTubs")
1186  { const G4TwistedTubs* const twistedtubsPtr
1187  = static_cast<const G4TwistedTubs*>(solidPtr);
1188  TwistedtubsWrite(solidsElement,twistedtubsPtr); }
1189  else
1190  {
1191  G4String error_msg = "Unknown solid: " + solidPtr->GetName()
1192  + "; Type: " + solidPtr->GetEntityType();
1193  G4Exception("G4GDMLWriteSolids::AddSolid()", "WriteError",
1194  FatalException, error_msg);
1195  }
1196 }
G4String GetName() const
void ScaledWrite(xercesc::DOMElement *, const G4ScaledSolid *const)
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)
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:66
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:

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

Definition at line 154 of file G4GDMLWriteSolids.cc.

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

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

Definition at line 303 of file G4GDMLWriteSolids.cc.

304 {
305  const G4String& name = GenerateName(box->GetName(),box);
306 
307  xercesc::DOMElement* boxElement = NewElement("box");
308  boxElement->setAttributeNode(NewAttribute("name",name));
309  boxElement->setAttributeNode(NewAttribute("x",2.0*box->GetXHalfLength()/mm));
310  boxElement->setAttributeNode(NewAttribute("y",2.0*box->GetYHalfLength()/mm));
311  boxElement->setAttributeNode(NewAttribute("z",2.0*box->GetZHalfLength()/mm));
312  boxElement->setAttributeNode(NewAttribute("lunit","mm"));
313  solElement->appendChild(boxElement);
314 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
G4double GetXHalfLength() const
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetZHalfLength() const
G4double GetYHalfLength() const
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:

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

Definition at line 317 of file G4GDMLWriteSolids.cc.

318 {
319  const G4String& name = GenerateName(cone->GetName(),cone);
320 
321  xercesc::DOMElement* coneElement = NewElement("cone");
322  coneElement->setAttributeNode(NewAttribute("name",name));
323  coneElement->
324  setAttributeNode(NewAttribute("rmin1",cone->GetInnerRadiusMinusZ()/mm));
325  coneElement->
326  setAttributeNode(NewAttribute("rmax1",cone->GetOuterRadiusMinusZ()/mm));
327  coneElement->
328  setAttributeNode(NewAttribute("rmin2",cone->GetInnerRadiusPlusZ()/mm));
329  coneElement->
330  setAttributeNode(NewAttribute("rmax2",cone->GetOuterRadiusPlusZ()/mm));
331  coneElement->
332  setAttributeNode(NewAttribute("z",2.0*cone->GetZHalfLength()/mm));
333  coneElement->
334  setAttributeNode(NewAttribute("startphi",cone->GetStartPhiAngle()/degree));
335  coneElement->
336  setAttributeNode(NewAttribute("deltaphi",cone->GetDeltaPhiAngle()/degree));
337  coneElement->setAttributeNode(NewAttribute("aunit","deg"));
338  coneElement->setAttributeNode(NewAttribute("lunit","mm"));
339  solElement->appendChild(coneElement);
340 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetOuterRadiusMinusZ() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetStartPhiAngle() 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 GetInnerRadiusMinusZ() const
G4double GetOuterRadiusPlusZ() const
G4double GetZHalfLength() const
G4double GetDeltaPhiAngle() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 880 of file G4GDMLWriteSolids.cc.

881 {
882  const G4String& name = GenerateName(cuttube->GetName(),cuttube);
883 
884  xercesc::DOMElement* cuttubeElement = NewElement("cutTube");
885  cuttubeElement->setAttributeNode(NewAttribute("name",name));
886  cuttubeElement->setAttributeNode(NewAttribute("rmin",
887  cuttube->GetInnerRadius()/mm));
888  cuttubeElement->setAttributeNode(NewAttribute("rmax",
889  cuttube->GetOuterRadius()/mm));
890  cuttubeElement->setAttributeNode(NewAttribute("z",
891  2.0*cuttube->GetZHalfLength()/mm));
892  cuttubeElement->setAttributeNode(NewAttribute("startphi",
893  cuttube->GetStartPhiAngle()/degree));
894  cuttubeElement->setAttributeNode(NewAttribute("deltaphi",
895  cuttube->GetDeltaPhiAngle()/degree));
896  cuttubeElement->setAttributeNode(NewAttribute("lowX",
897  cuttube->GetLowNorm().getX()/mm));
898  cuttubeElement->setAttributeNode(NewAttribute("lowY",
899  cuttube->GetLowNorm().getY()/mm));
900  cuttubeElement->setAttributeNode(NewAttribute("lowZ",
901  cuttube->GetLowNorm().getZ()/mm));
902  cuttubeElement->setAttributeNode(NewAttribute("highX",
903  cuttube->GetHighNorm().getX()/mm));
904  cuttubeElement->setAttributeNode(NewAttribute("highY",
905  cuttube->GetHighNorm().getY()/mm));
906  cuttubeElement->setAttributeNode(NewAttribute("highZ",
907  cuttube->GetHighNorm().getZ()/mm));
908  cuttubeElement->setAttributeNode(NewAttribute("aunit","deg"));
909  cuttubeElement->setAttributeNode(NewAttribute("lunit","mm"));
910  solElement->appendChild(cuttubeElement);
911 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4ThreeVector GetLowNorm() const
double getY() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetStartPhiAngle() const
double getX() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetZHalfLength() const
G4double GetOuterRadius() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
double getZ() const
G4ThreeVector GetHighNorm() const
G4double GetInnerRadius() const
G4double GetDeltaPhiAngle() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 343 of file G4GDMLWriteSolids.cc.

345 {
346  const G4String& name = GenerateName(elcone->GetName(),elcone);
347 
348  xercesc::DOMElement* elconeElement = NewElement("elcone");
349  elconeElement->setAttributeNode(NewAttribute("name",name));
350  elconeElement->setAttributeNode(NewAttribute("dx",elcone->GetSemiAxisX()/mm));
351  elconeElement->setAttributeNode(NewAttribute("dy",elcone->GetSemiAxisY()/mm));
352  elconeElement->setAttributeNode(NewAttribute("zmax",elcone->GetZMax()/mm));
353  elconeElement->setAttributeNode(NewAttribute("zcut",elcone->GetZTopCut()/mm));
354  elconeElement->setAttributeNode(NewAttribute("lunit","mm"));
355  solElement->appendChild(elconeElement);
356 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetSemiAxisX() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetZMax() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetSemiAxisY() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetZTopCut() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 359 of file G4GDMLWriteSolids.cc.

361 {
362  const G4String& name = GenerateName(ellipsoid->GetName(),ellipsoid);
363 
364  xercesc::DOMElement* ellipsoidElement = NewElement("ellipsoid");
365  ellipsoidElement->setAttributeNode(NewAttribute("name",name));
366  ellipsoidElement->
367  setAttributeNode(NewAttribute("ax",ellipsoid->GetSemiAxisMax(0)/mm));
368  ellipsoidElement->
369  setAttributeNode(NewAttribute("by",ellipsoid->GetSemiAxisMax(1)/mm));
370  ellipsoidElement->
371  setAttributeNode(NewAttribute("cz",ellipsoid->GetSemiAxisMax(2)/mm));
372  ellipsoidElement->
373  setAttributeNode(NewAttribute("zcut1",ellipsoid->GetZBottomCut()/mm));
374  ellipsoidElement->
375  setAttributeNode(NewAttribute("zcut2",ellipsoid->GetZTopCut()/mm));
376  ellipsoidElement->
377  setAttributeNode(NewAttribute("lunit","mm"));
378  solElement->appendChild(ellipsoidElement);
379 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetZTopCut() const
G4double GetZBottomCut() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetSemiAxisMax(G4int i) const
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:

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

Definition at line 382 of file G4GDMLWriteSolids.cc.

384 {
385  const G4String& name = GenerateName(eltube->GetName(),eltube);
386 
387  xercesc::DOMElement* eltubeElement = NewElement("eltube");
388  eltubeElement->setAttributeNode(NewAttribute("name",name));
389  eltubeElement->setAttributeNode(NewAttribute("dx",eltube->GetDx()/mm));
390  eltubeElement->setAttributeNode(NewAttribute("dy",eltube->GetDy()/mm));
391  eltubeElement->setAttributeNode(NewAttribute("dz",eltube->GetDz()/mm));
392  eltubeElement->setAttributeNode(NewAttribute("lunit","mm"));
393  solElement->appendChild(eltubeElement);
394 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetDy() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetDz() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetDx() const
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:

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

Definition at line 547 of file G4GDMLWriteSolids.cc.

549 {
550  const G4String& name = GenerateName(polycone->GetName(),polycone);
551  xercesc::DOMElement* polyconeElement = NewElement("genericPolycone");
552  const G4double startPhi=polycone->GetStartPhi();
553  polyconeElement->setAttributeNode(NewAttribute("name",name));
554  polyconeElement->setAttributeNode(NewAttribute("startphi",
555  startPhi/degree));
556  polyconeElement->setAttributeNode(NewAttribute("deltaphi",
557  (polycone->GetEndPhi()-startPhi)/degree));
558  polyconeElement->setAttributeNode(NewAttribute("aunit","deg"));
559  polyconeElement->setAttributeNode(NewAttribute("lunit","mm"));
560  solElement->appendChild(polyconeElement);
561 
562  const size_t num_rzpoints = polycone->GetNumRZCorner();
563  for (size_t i=0; i<num_rzpoints; i++)
564  {
565  const G4double r_point=polycone->GetCorner(i).r;
566  const G4double z_point=polycone->GetCorner(i).z;
567  RZPointWrite(polyconeElement,r_point,z_point);
568  }
569 
570 }
G4PolyconeSideRZ GetCorner(G4int index) const
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetEndPhi() const
void RZPointWrite(xercesc::DOMElement *, const G4double &, const G4double &)
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetStartPhi() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
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:

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

Definition at line 769 of file G4GDMLWriteSolids.cc.

771 {
772  const G4String& name = GenerateName(gtrap->GetName(),gtrap);
773 
774  std::vector<G4TwoVector> vertices = gtrap->GetVertices();
775 
776  xercesc::DOMElement* gtrapElement = NewElement("arb8");
777  gtrapElement->setAttributeNode(NewAttribute("name",name));
778  gtrapElement->setAttributeNode(NewAttribute("dz",
779  gtrap->GetZHalfLength()/mm));
780  gtrapElement->setAttributeNode(NewAttribute("v1x", vertices[0].x()));
781  gtrapElement->setAttributeNode(NewAttribute("v1y", vertices[0].y()));
782  gtrapElement->setAttributeNode(NewAttribute("v2x", vertices[1].x()));
783  gtrapElement->setAttributeNode(NewAttribute("v2y", vertices[1].y()));
784  gtrapElement->setAttributeNode(NewAttribute("v3x", vertices[2].x()));
785  gtrapElement->setAttributeNode(NewAttribute("v3y", vertices[2].y()));
786  gtrapElement->setAttributeNode(NewAttribute("v4x", vertices[3].x()));
787  gtrapElement->setAttributeNode(NewAttribute("v4y", vertices[3].y()));
788  gtrapElement->setAttributeNode(NewAttribute("v5x", vertices[4].x()));
789  gtrapElement->setAttributeNode(NewAttribute("v5y", vertices[4].y()));
790  gtrapElement->setAttributeNode(NewAttribute("v6x", vertices[5].x()));
791  gtrapElement->setAttributeNode(NewAttribute("v6y", vertices[5].y()));
792  gtrapElement->setAttributeNode(NewAttribute("v7x", vertices[6].x()));
793  gtrapElement->setAttributeNode(NewAttribute("v7y", vertices[6].y()));
794  gtrapElement->setAttributeNode(NewAttribute("v8x", vertices[7].x()));
795  gtrapElement->setAttributeNode(NewAttribute("v8y", vertices[7].y()));
796  gtrapElement->setAttributeNode(NewAttribute("lunit","mm"));
797  solElement->appendChild(gtrapElement);
798 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
tuple x
Definition: test.py:50
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
const std::vector< G4TwoVector > & GetVertices() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 441 of file G4GDMLWriteSolids.cc.

442 {
443  const G4String& name = GenerateName(hype->GetName(),hype);
444 
445  xercesc::DOMElement* hypeElement = NewElement("hype");
446  hypeElement->setAttributeNode(NewAttribute("name",name));
447  hypeElement->setAttributeNode(NewAttribute("rmin",
448  hype->GetInnerRadius()/mm));
449  hypeElement->setAttributeNode(NewAttribute("rmax",
450  hype->GetOuterRadius()/mm));
451  hypeElement->setAttributeNode(NewAttribute("inst",
452  hype->GetInnerStereo()/degree));
453  hypeElement->setAttributeNode(NewAttribute("outst",
454  hype->GetOuterStereo()/degree));
455  hypeElement->setAttributeNode(NewAttribute("z",
456  2.0*hype->GetZHalfLength()/mm));
457  hypeElement->setAttributeNode(NewAttribute("aunit","deg"));
458  hypeElement->setAttributeNode(NewAttribute("lunit","mm"));
459  solElement->appendChild(hypeElement);
460 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
G4double GetOuterStereo() const
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetInnerStereo() const
G4double GetInnerRadius() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
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
G4double GetZHalfLength() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 84 of file G4GDMLWriteSolids.cc.

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

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

Definition at line 1040 of file G4GDMLWriteSolids.cc.

1042 {
1043  xercesc::DOMElement* optElement = NewElement("opticalsurface");
1044  G4OpticalSurfaceModel smodel = surf->GetModel();
1045  G4double sval = (smodel==glisur) ? surf->GetPolish() : surf->GetSigmaAlpha();
1046 
1047  optElement->setAttributeNode(NewAttribute("name", surf->GetName()));
1048  optElement->setAttributeNode(NewAttribute("model", smodel));
1049  optElement->setAttributeNode(NewAttribute("finish", surf->GetFinish()));
1050  optElement->setAttributeNode(NewAttribute("type", surf->GetType()));
1051  optElement->setAttributeNode(NewAttribute("value", sval));
1052 
1053  solElement->appendChild(optElement);
1054 }
G4OpticalSurfaceModel GetModel() const
const G4String & GetName() const
const G4SurfaceType & GetType() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetPolish() const
G4OpticalSurfaceFinish GetFinish() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4OpticalSurfaceModel
double G4double
Definition: G4Types.hh:76
G4double GetSigmaAlpha() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 463 of file G4GDMLWriteSolids.cc.

464 {
465  const G4String& name = GenerateName(orb->GetName(),orb);
466 
467  xercesc::DOMElement* orbElement = NewElement("orb");
468  orbElement->setAttributeNode(NewAttribute("name",name));
469  orbElement->setAttributeNode(NewAttribute("r",orb->GetRadius()/mm));
470  orbElement->setAttributeNode(NewAttribute("lunit","mm"));
471  solElement->appendChild(orbElement);
472 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 501 of file G4GDMLWriteSolids.cc.

503 {
504  const G4String& name = GenerateName(paraboloid->GetName(),paraboloid);
505 
506  xercesc::DOMElement* paraboloidElement = NewElement("paraboloid");
507  paraboloidElement->setAttributeNode(NewAttribute("name",name));
508  paraboloidElement->setAttributeNode(NewAttribute("rlo",
509  paraboloid->GetRadiusMinusZ()/mm));
510  paraboloidElement->setAttributeNode(NewAttribute("rhi",
511  paraboloid->GetRadiusPlusZ()/mm));
512  paraboloidElement->setAttributeNode(NewAttribute("dz",
513  paraboloid->GetZHalfLength()/mm));
514  paraboloidElement->setAttributeNode(NewAttribute("lunit","mm"));
515  solElement->appendChild(paraboloidElement);
516 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetRadiusMinusZ() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetZHalfLength() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetRadiusPlusZ() const
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:

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

Definition at line 475 of file G4GDMLWriteSolids.cc.

476 {
477  const G4String& name = GenerateName(para->GetName(),para);
478 
479  const G4ThreeVector simaxis = para->GetSymAxis();
480  const G4double alpha = std::atan(para->GetTanAlpha());
481  const G4double phi = simaxis.phi();
482  const G4double theta = simaxis.theta();
483 
484  xercesc::DOMElement* paraElement = NewElement("para");
485  paraElement->setAttributeNode(NewAttribute("name",name));
486  paraElement->setAttributeNode(NewAttribute("x",
487  2.0*para->GetXHalfLength()/mm));
488  paraElement->setAttributeNode(NewAttribute("y",
489  2.0*para->GetYHalfLength()/mm));
490  paraElement->setAttributeNode(NewAttribute("z",
491  2.0*para->GetZHalfLength()/mm));
492  paraElement->setAttributeNode(NewAttribute("alpha",alpha/degree));
493  paraElement->setAttributeNode(NewAttribute("theta",theta/degree));
494  paraElement->setAttributeNode(NewAttribute("phi",phi/degree));
495  paraElement->setAttributeNode(NewAttribute("aunit","deg"));
496  paraElement->setAttributeNode(NewAttribute("lunit","mm"));
497  solElement->appendChild(paraElement);
498 }
G4ThreeVector GetSymAxis() const
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetXHalfLength() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetTanAlpha() const
G4double GetZHalfLength() const
double G4double
Definition: G4Types.hh:76
static const G4double alpha
G4double GetYHalfLength() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 518 of file G4GDMLWriteSolids.cc.

520 {
521  const G4String& name = GenerateName(polycone->GetName(),polycone);
522 
523  xercesc::DOMElement* polyconeElement = NewElement("polycone");
524  polyconeElement->setAttributeNode(NewAttribute("name",name));
525  polyconeElement->setAttributeNode(NewAttribute("startphi",
527  polyconeElement->setAttributeNode(NewAttribute("deltaphi",
529  polyconeElement->setAttributeNode(NewAttribute("aunit","deg"));
530  polyconeElement->setAttributeNode(NewAttribute("lunit","mm"));
531  solElement->appendChild(polyconeElement);
532 
533  const size_t num_zplanes = polycone->GetOriginalParameters()->Num_z_planes;
534  const G4double* z_array = polycone->GetOriginalParameters()->Z_values;
535  const G4double* rmin_array = polycone->GetOriginalParameters()->Rmin;
536  const G4double* rmax_array = polycone->GetOriginalParameters()->Rmax;
537 
538  for (size_t i=0; i<num_zplanes; i++)
539  {
540  ZplaneWrite(polyconeElement,z_array[i],rmin_array[i],rmax_array[i]);
541  }
542 
543 
544 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
static constexpr double degree
Definition: G4SIunits.hh:144
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
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 573 of file G4GDMLWriteSolids.cc.

575 {
576  const G4String& name = GenerateName(polyhedra->GetName(),polyhedra);
577  if(polyhedra->IsGeneric() == false){
578  xercesc::DOMElement* polyhedraElement = NewElement("polyhedra");
579  polyhedraElement->setAttributeNode(NewAttribute("name",name));
580  polyhedraElement->setAttributeNode(NewAttribute("startphi",
581  polyhedra->GetOriginalParameters()->Start_angle/degree));
582  polyhedraElement->setAttributeNode(NewAttribute("deltaphi",
584  polyhedraElement->setAttributeNode(NewAttribute("numsides",
585  polyhedra->GetOriginalParameters()->numSide));
586  polyhedraElement->setAttributeNode(NewAttribute("aunit","deg"));
587  polyhedraElement->setAttributeNode(NewAttribute("lunit","mm"));
588  solElement->appendChild(polyhedraElement);
589 
590  const size_t num_zplanes = polyhedra->GetOriginalParameters()->Num_z_planes;
591  const G4double* z_array = polyhedra->GetOriginalParameters()->Z_values;
592  const G4double* rmin_array = polyhedra->GetOriginalParameters()->Rmin;
593  const G4double* rmax_array = polyhedra->GetOriginalParameters()->Rmax;
594 
595  const G4double convertRad =
596  std::cos(0.5*polyhedra->GetOriginalParameters()->Opening_angle
597  / polyhedra->GetOriginalParameters()->numSide);
598 
599  for (size_t i=0;i<num_zplanes;i++)
600  {
601  ZplaneWrite(polyhedraElement,z_array[i],
602  rmin_array[i]*convertRad, rmax_array[i]*convertRad);
603  }
604  }else{//generic polyhedra
605 
606  xercesc::DOMElement* polyhedraElement = NewElement("genericPolyhedra");
607  polyhedraElement->setAttributeNode(NewAttribute("name",name));
608  polyhedraElement->setAttributeNode(NewAttribute("startphi",
609  polyhedra->GetOriginalParameters()->Start_angle/degree));
610  polyhedraElement->setAttributeNode(NewAttribute("deltaphi",
612  polyhedraElement->setAttributeNode(NewAttribute("numsides",
613  polyhedra->GetOriginalParameters()->numSide));
614  polyhedraElement->setAttributeNode(NewAttribute("aunit","deg"));
615  polyhedraElement->setAttributeNode(NewAttribute("lunit","mm"));
616  solElement->appendChild(polyhedraElement);
617 
618  const size_t num_rzpoints = polyhedra->GetNumRZCorner();
619 
620  for (size_t i=0;i<num_rzpoints;i++)
621  {
622  const G4double r_point = polyhedra->GetCorner(i).r;
623  const G4double z_point = polyhedra->GetCorner(i).z;
624  RZPointWrite(polyhedraElement,r_point,z_point);
625  }
626  }
627 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4int GetNumRZCorner() const
void RZPointWrite(xercesc::DOMElement *, const G4double &, const G4double &)
static constexpr double degree
Definition: G4SIunits.hh:144
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
G4bool IsGeneric() const
G4PolyhedraHistorical * GetOriginalParameters() const
G4PolyhedraSideRZ GetCorner(const G4int index) const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 1030 of file G4GDMLWriteSolids.cc.

1032 {
1033  xercesc::DOMElement* rzpointElement = NewElement("rzpoint");
1034  rzpointElement->setAttributeNode(NewAttribute("r",r/mm));
1035  rzpointElement->setAttributeNode(NewAttribute("z",z/mm));
1036  element->appendChild(rzpointElement);
1037 }
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
tuple z
Definition: test.py:28

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteSolids::ScaledWrite ( xercesc::DOMElement *  solElement,
const G4ScaledSolid * const  scaled 
)
protected

Definition at line 270 of file G4GDMLWriteSolids.cc.

272 {
273  G4String tag("scaledSolid");
274 
275  G4VSolid* solid = const_cast<G4VSolid*>(scaled->GetUnscaledSolid());
276  G4Scale3D scale = scaled->GetScaleTransform();
277  G4ThreeVector sclVector = G4ThreeVector(scale.xx(), scale.yy(), scale.zz());
278 
279  AddSolid(solid); // Add the constituent solid!
280 
281  const G4String& name = GenerateName(scaled->GetName(),scaled);
282  const G4String& solidref = GenerateName(solid->GetName(),solid);
283 
284  xercesc::DOMElement* scaledElement = NewElement(tag);
285  scaledElement->setAttributeNode(NewAttribute("name",name));
286 
287  xercesc::DOMElement* solidElement = NewElement("solidref");
288  solidElement->setAttributeNode(NewAttribute("ref",solidref));
289  scaledElement->appendChild(solidElement);
290 
291  if ( (std::fabs(scale.xx()) > kLinearPrecision)
292  && (std::fabs(scale.yy()) > kLinearPrecision)
293  && (std::fabs(scale.zz()) > kLinearPrecision) )
294  {
295  ScaleWrite(scaledElement, name+"_scl", sclVector);
296  }
297 
298  solElement->appendChild(scaledElement);
299  // Add the scaled solid AFTER its constituent solid!
300 }
G4String GetName() const
double yy() const
Definition: Transform3D.h:264
const XML_Char * name
Definition: expat.h:151
G4Scale3D GetScaleTransform() const
CLHEP::Hep3Vector G4ThreeVector
Definition: xmlparse.cc:187
double zz() const
Definition: Transform3D.h:276
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void ScaleWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const G4double kLinearPrecision
G4VSolid * GetUnscaledSolid() const
double xx() const
Definition: Transform3D.h:252
virtual void AddSolid(const G4VSolid *const)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Implements G4GDMLWrite.

Definition at line 1056 of file G4GDMLWriteSolids.cc.

1057 {
1058  G4cout << "G4GDML: Writing solids..." << G4endl;
1059 
1060  solidsElement = NewElement("solids");
1061  gdmlElement->appendChild(solidsElement);
1062 
1063  solidList.clear();
1064 }
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:

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

Definition at line 630 of file G4GDMLWriteSolids.cc.

631 {
632  const G4String& name = GenerateName(sphere->GetName(),sphere);
633 
634  xercesc::DOMElement* sphereElement = NewElement("sphere");
635  sphereElement->setAttributeNode(NewAttribute("name",name));
636  sphereElement->setAttributeNode(NewAttribute("rmin",
637  sphere->GetInnerRadius()/mm));
638  sphereElement->setAttributeNode(NewAttribute("rmax",
639  sphere->GetOuterRadius()/mm));
640  sphereElement->setAttributeNode(NewAttribute("startphi",
641  sphere->GetStartPhiAngle()/degree));
642  sphereElement->setAttributeNode(NewAttribute("deltaphi",
643  sphere->GetDeltaPhiAngle()/degree));
644  sphereElement->setAttributeNode(NewAttribute("starttheta",
645  sphere->GetStartThetaAngle()/degree));
646  sphereElement->setAttributeNode(NewAttribute("deltatheta",
647  sphere->GetDeltaThetaAngle()/degree));
648  sphereElement->setAttributeNode(NewAttribute("aunit","deg"));
649  sphereElement->setAttributeNode(NewAttribute("lunit","mm"));
650  solElement->appendChild(sphereElement);
651 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetDeltaPhiAngle() const
G4double GetStartThetaAngle() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetInnerRadius() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4double GetDeltaThetaAngle() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 654 of file G4GDMLWriteSolids.cc.

656 {
657  const G4String& solid_name = tessellated->GetName();
658  const G4String& name = GenerateName(solid_name, tessellated);
659 
660  xercesc::DOMElement* tessellatedElement = NewElement("tessellated");
661  tessellatedElement->setAttributeNode(NewAttribute("name",name));
662  tessellatedElement->setAttributeNode(NewAttribute("aunit","deg"));
663  tessellatedElement->setAttributeNode(NewAttribute("lunit","mm"));
664  solElement->appendChild(tessellatedElement);
665 
666  std::map<G4ThreeVector, G4String, G4ThreeVectorCompare> vertexMap;
667 
668  const size_t NumFacets = tessellated->GetNumberOfFacets();
669  size_t NumVertex = 0;
670 
671  for (size_t i=0;i<NumFacets;i++)
672  {
673  const G4VFacet* facet = tessellated->GetFacet(i);
674  const size_t NumVertexPerFacet = facet->GetNumberOfVertices();
675 
676  G4String FacetTag;
677 
678  if (NumVertexPerFacet==3) { FacetTag="triangular"; } else
679  if (NumVertexPerFacet==4) { FacetTag="quadrangular"; }
680  else
681  {
682  G4Exception("G4GDMLWriteSolids::TessellatedWrite()", "InvalidSetup",
683  FatalException, "Facet should contain 3 or 4 vertices!");
684  }
685 
686  xercesc::DOMElement* facetElement = NewElement(FacetTag);
687  tessellatedElement->appendChild(facetElement);
688 
689  for (size_t j=0; j<NumVertexPerFacet; j++)
690  {
691  std::stringstream name_stream;
692  std::stringstream ref_stream;
693 
694  name_stream << "vertex" << (j+1);
695  ref_stream << solid_name << "_v" << NumVertex;
696 
697  const G4String& fname = name_stream.str(); // facet's tag variable
698  G4String ref = ref_stream.str(); // vertex tag to be associated
699 
700  // Now search for the existance of the current vertex in the
701  // map of cached vertices. If existing, do NOT store it as
702  // position in the GDML file, so avoiding duplication; otherwise
703  // cache it in the local map and add it as position in the
704  // "define" section of the GDML file.
705 
706  const G4ThreeVector& vertex = facet->GetVertex(j);
707 
708  if(vertexMap.find(vertex) != vertexMap.end()) // Vertex is cached
709  {
710  ref = vertexMap[vertex]; // Set the proper tag for it
711  }
712  else // Vertex not found
713  {
714  vertexMap.insert(std::make_pair(vertex,ref)); // Cache vertex and ...
715  AddPosition(ref, vertex); // ... add it to define section!
716  NumVertex++;
717  }
718 
719  // Now create association of the vertex with its facet
720  //
721  facetElement->setAttributeNode(NewAttribute(fname,ref));
722  }
723  }
724 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
virtual G4int GetNumberOfVertices() const =0
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4VFacet * GetFacet(G4int i) 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
G4int GetNumberOfFacets() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
string fname
Definition: test.py:308
virtual G4ThreeVector GetVertex(G4int i) const =0

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 727 of file G4GDMLWriteSolids.cc.

728 {
729  const G4String& solid_name = tet->GetName();
730  const G4String& name = GenerateName(solid_name, tet);
731 
732  std::vector<G4ThreeVector> vertexList = tet->GetVertices();
733 
734  xercesc::DOMElement* tetElement = NewElement("tet");
735  tetElement->setAttributeNode(NewAttribute("name",name));
736  tetElement->setAttributeNode(NewAttribute("vertex1",solid_name+"_v1"));
737  tetElement->setAttributeNode(NewAttribute("vertex2",solid_name+"_v2"));
738  tetElement->setAttributeNode(NewAttribute("vertex3",solid_name+"_v3"));
739  tetElement->setAttributeNode(NewAttribute("vertex4",solid_name+"_v4"));
740  tetElement->setAttributeNode(NewAttribute("lunit","mm"));
741  solElement->appendChild(tetElement);
742 
743  AddPosition(solid_name+"_v1",vertexList[0]);
744  AddPosition(solid_name+"_v2",vertexList[1]);
745  AddPosition(solid_name+"_v3",vertexList[2]);
746  AddPosition(solid_name+"_v4",vertexList[3]);
747 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:733
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:

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

Definition at line 750 of file G4GDMLWriteSolids.cc.

751 {
752  const G4String& name = GenerateName(torus->GetName(),torus);
753 
754  xercesc::DOMElement* torusElement = NewElement("torus");
755  torusElement->setAttributeNode(NewAttribute("name",name));
756  torusElement->setAttributeNode(NewAttribute("rmin",torus->GetRmin()/mm));
757  torusElement->setAttributeNode(NewAttribute("rmax",torus->GetRmax()/mm));
758  torusElement->setAttributeNode(NewAttribute("rtor",torus->GetRtor()/mm));
759  torusElement->
760  setAttributeNode(NewAttribute("startphi",torus->GetSPhi()/degree));
761  torusElement->
762  setAttributeNode(NewAttribute("deltaphi",torus->GetDPhi()/degree));
763  torusElement->setAttributeNode(NewAttribute("aunit","deg"));
764  torusElement->setAttributeNode(NewAttribute("lunit","mm"));
765  solElement->appendChild(torusElement);
766 }
G4double GetSPhi() const
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetRmax() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetRtor() const
G4double GetRmin() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetDPhi() const
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:

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

Definition at line 801 of file G4GDMLWriteSolids.cc.

802 {
803  const G4String& name = GenerateName(trap->GetName(),trap);
804 
805  const G4ThreeVector& simaxis = trap->GetSymAxis();
806  const G4double phi = simaxis.phi();
807  const G4double theta = simaxis.theta();
808  const G4double alpha1 = std::atan(trap->GetTanAlpha1());
809  const G4double alpha2 = std::atan(trap->GetTanAlpha2());
810 
811  xercesc::DOMElement* trapElement = NewElement("trap");
812  trapElement->setAttributeNode(NewAttribute("name",name));
813  trapElement->setAttributeNode(NewAttribute("z",
814  2.0*trap->GetZHalfLength()/mm));
815  trapElement->setAttributeNode(NewAttribute("theta",theta/degree));
816  trapElement->setAttributeNode(NewAttribute("phi",phi/degree));
817  trapElement->setAttributeNode(NewAttribute("y1",
818  2.0*trap->GetYHalfLength1()/mm));
819  trapElement->setAttributeNode(NewAttribute("x1",
820  2.0*trap->GetXHalfLength1()/mm));
821  trapElement->setAttributeNode(NewAttribute("x2",
822  2.0*trap->GetXHalfLength2()/mm));
823  trapElement->setAttributeNode(NewAttribute("alpha1",alpha1/degree));
824  trapElement->setAttributeNode(NewAttribute("y2",
825  2.0*trap->GetYHalfLength2()/mm));
826  trapElement->setAttributeNode(NewAttribute("x3",
827  2.0*trap->GetXHalfLength3()/mm));
828  trapElement->setAttributeNode(NewAttribute("x4",
829  2.0*trap->GetXHalfLength4()/mm));
830  trapElement->setAttributeNode(NewAttribute("alpha2",alpha2/degree));
831  trapElement->setAttributeNode(NewAttribute("aunit","deg"));
832  trapElement->setAttributeNode(NewAttribute("lunit","mm"));
833  solElement->appendChild(trapElement);
834 }
G4String GetName() const
G4double GetXHalfLength4() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetYHalfLength2() const
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetXHalfLength2() const
G4double GetTanAlpha2() const
G4double GetXHalfLength1() const
G4double GetXHalfLength3() const
static constexpr double degree
Definition: G4SIunits.hh:144
double phi() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4ThreeVector GetSymAxis() const
G4double GetYHalfLength1() const
double G4double
Definition: G4Types.hh:76
G4double GetTanAlpha1() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 837 of file G4GDMLWriteSolids.cc.

838 {
839  const G4String& name = GenerateName(trd->GetName(),trd);
840 
841  xercesc::DOMElement* trdElement = NewElement("trd");
842  trdElement->setAttributeNode(NewAttribute("name",name));
843  trdElement->setAttributeNode(NewAttribute("x1",
844  2.0*trd->GetXHalfLength1()/mm));
845  trdElement->setAttributeNode(NewAttribute("x2",
846  2.0*trd->GetXHalfLength2()/mm));
847  trdElement->setAttributeNode(NewAttribute("y1",
848  2.0*trd->GetYHalfLength1()/mm));
849  trdElement->setAttributeNode(NewAttribute("y2",
850  2.0*trd->GetYHalfLength2()/mm));
851  trdElement->setAttributeNode(NewAttribute("z",
852  2.0*trd->GetZHalfLength()/mm));
853  trdElement->setAttributeNode(NewAttribute("lunit","mm"));
854  solElement->appendChild(trdElement);
855 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetYHalfLength1() const
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetXHalfLength2() const
G4double GetYHalfLength2() 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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 858 of file G4GDMLWriteSolids.cc.

859 {
860  const G4String& name = GenerateName(tube->GetName(),tube);
861 
862  xercesc::DOMElement* tubeElement = NewElement("tube");
863  tubeElement->setAttributeNode(NewAttribute("name",name));
864  tubeElement->setAttributeNode(NewAttribute("rmin",
865  tube->GetInnerRadius()/mm));
866  tubeElement->setAttributeNode(NewAttribute("rmax",
867  tube->GetOuterRadius()/mm));
868  tubeElement->setAttributeNode(NewAttribute("z",
869  2.0*tube->GetZHalfLength()/mm));
870  tubeElement->setAttributeNode(NewAttribute("startphi",
871  tube->GetStartPhiAngle()/degree));
872  tubeElement->setAttributeNode(NewAttribute("deltaphi",
873  tube->GetDeltaPhiAngle()/degree));
874  tubeElement->setAttributeNode(NewAttribute("aunit","deg"));
875  tubeElement->setAttributeNode(NewAttribute("lunit","mm"));
876  solElement->appendChild(tubeElement);
877 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetDeltaPhiAngle() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetStartPhiAngle() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetInnerRadius() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetZHalfLength() const
G4double GetOuterRadius() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 914 of file G4GDMLWriteSolids.cc.

916 {
917  const G4String& name = GenerateName(twistedbox->GetName(),twistedbox);
918 
919  xercesc::DOMElement* twistedboxElement = NewElement("twistedbox");
920  twistedboxElement->setAttributeNode(NewAttribute("name",name));
921  twistedboxElement->setAttributeNode(NewAttribute("x",
922  2.0*twistedbox->GetXHalfLength()/mm));
923  twistedboxElement->setAttributeNode(NewAttribute("y",
924  2.0*twistedbox->GetYHalfLength()/mm));
925  twistedboxElement->setAttributeNode(NewAttribute("z",
926  2.0*twistedbox->GetZHalfLength()/mm));
927  twistedboxElement->setAttributeNode(NewAttribute("PhiTwist",
928  twistedbox->GetPhiTwist()/degree));
929  twistedboxElement->setAttributeNode(NewAttribute("aunit","deg"));
930  twistedboxElement->setAttributeNode(NewAttribute("lunit","mm"));
931  solElement->appendChild(twistedboxElement);
932 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetZHalfLength() const
Definition: G4TwistedBox.hh:75
static constexpr double degree
Definition: G4SIunits.hh:144
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetYHalfLength() const
Definition: G4TwistedBox.hh:74
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetPhiTwist() const
Definition: G4TwistedBox.hh:76
G4double GetXHalfLength() const
Definition: G4TwistedBox.hh:73

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 935 of file G4GDMLWriteSolids.cc.

937 {
938  const G4String& name = GenerateName(twistedtrap->GetName(),twistedtrap);
939 
940  xercesc::DOMElement* twistedtrapElement = NewElement("twistedtrap");
941  twistedtrapElement->setAttributeNode(NewAttribute("name",name));
942  twistedtrapElement->setAttributeNode(NewAttribute("y1",
943  2.0*twistedtrap->GetY1HalfLength()/mm));
944  twistedtrapElement->setAttributeNode(NewAttribute("x1",
945  2.0*twistedtrap->GetX1HalfLength()/mm));
946  twistedtrapElement->setAttributeNode(NewAttribute("x2",
947  2.0*twistedtrap->GetX2HalfLength()/mm));
948  twistedtrapElement->setAttributeNode(NewAttribute("y2",
949  2.0*twistedtrap->GetY2HalfLength()/mm));
950  twistedtrapElement->setAttributeNode(NewAttribute("x3",
951  2.0*twistedtrap->GetX3HalfLength()/mm));
952  twistedtrapElement->setAttributeNode(NewAttribute("x4",
953  2.0*twistedtrap->GetX4HalfLength()/mm));
954  twistedtrapElement->setAttributeNode(NewAttribute("z",
955  2.0*twistedtrap->GetZHalfLength()/mm));
956  twistedtrapElement->setAttributeNode(NewAttribute("Alph",
957  twistedtrap->GetTiltAngleAlpha()/degree));
958  twistedtrapElement->setAttributeNode(NewAttribute("Theta",
959  twistedtrap->GetPolarAngleTheta()/degree));
960  twistedtrapElement->setAttributeNode(NewAttribute("Phi",
961  twistedtrap->GetAzimuthalAnglePhi()/degree));
962  twistedtrapElement->setAttributeNode(NewAttribute("PhiTwist",
963  twistedtrap->GetPhiTwist()/degree));
964  twistedtrapElement->setAttributeNode(NewAttribute("aunit","deg"));
965  twistedtrapElement->setAttributeNode(NewAttribute("lunit","mm"));
966 
967  solElement->appendChild(twistedtrapElement);
968 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetPolarAngleTheta() const
G4double GetZHalfLength() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetY1HalfLength() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetX3HalfLength() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetX4HalfLength() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetAzimuthalAnglePhi() const
G4double GetX1HalfLength() const
G4double GetTiltAngleAlpha() const
G4double GetX2HalfLength() const
G4double GetPhiTwist() const
G4double GetY2HalfLength() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 971 of file G4GDMLWriteSolids.cc.

973 {
974  const G4String& name = GenerateName(twistedtrd->GetName(),twistedtrd);
975 
976  xercesc::DOMElement* twistedtrdElement = NewElement("twistedtrd");
977  twistedtrdElement->setAttributeNode(NewAttribute("name",name));
978  twistedtrdElement->setAttributeNode(NewAttribute("x1",
979  2.0*twistedtrd->GetX1HalfLength()/mm));
980  twistedtrdElement->setAttributeNode(NewAttribute("x2",
981  2.0*twistedtrd->GetX2HalfLength()/mm));
982  twistedtrdElement->setAttributeNode(NewAttribute("y1",
983  2.0*twistedtrd->GetY1HalfLength()/mm));
984  twistedtrdElement->setAttributeNode(NewAttribute("y2",
985  2.0*twistedtrd->GetY2HalfLength()/mm));
986  twistedtrdElement->setAttributeNode(NewAttribute("z",
987  2.0*twistedtrd->GetZHalfLength()/mm));
988  twistedtrdElement->setAttributeNode(NewAttribute("PhiTwist",
989  twistedtrd->GetPhiTwist()/degree));
990  twistedtrdElement->setAttributeNode(NewAttribute("aunit","deg"));
991  twistedtrdElement->setAttributeNode(NewAttribute("lunit","mm"));
992  solElement->appendChild(twistedtrdElement);
993 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
G4double GetX1HalfLength() const
Definition: G4TwistedTrd.hh:77
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetY1HalfLength() const
Definition: G4TwistedTrd.hh:79
G4double GetX2HalfLength() const
Definition: G4TwistedTrd.hh:78
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetZHalfLength() const
Definition: G4TwistedTrd.hh:81
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetY2HalfLength() const
Definition: G4TwistedTrd.hh:80
G4double GetPhiTwist() const
Definition: G4TwistedTrd.hh:82
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:

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

Definition at line 996 of file G4GDMLWriteSolids.cc.

998 {
999  const G4String& name = GenerateName(twistedtubs->GetName(),twistedtubs);
1000 
1001  xercesc::DOMElement* twistedtubsElement = NewElement("twistedtubs");
1002  twistedtubsElement->setAttributeNode(NewAttribute("name",name));
1003  twistedtubsElement->setAttributeNode(NewAttribute("twistedangle",
1004  twistedtubs->GetPhiTwist()/degree));
1005  twistedtubsElement->setAttributeNode(NewAttribute("endinnerrad",
1006  twistedtubs->GetInnerRadius()/mm));
1007  twistedtubsElement->setAttributeNode(NewAttribute("endouterrad",
1008  twistedtubs->GetOuterRadius()/mm));
1009  twistedtubsElement->setAttributeNode(NewAttribute("zlen",
1010  2.0*twistedtubs->GetZHalfLength()/mm));
1011  twistedtubsElement->setAttributeNode(NewAttribute("phi",
1012  twistedtubs->GetDPhi()/degree));
1013  twistedtubsElement->setAttributeNode(NewAttribute("aunit","deg"));
1014  twistedtubsElement->setAttributeNode(NewAttribute("lunit","mm"));
1015  solElement->appendChild(twistedtubsElement);
1016 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetInnerRadius() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4double GetPhiTwist() const
static constexpr double degree
Definition: G4SIunits.hh:144
G4double GetDPhi() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4double GetZHalfLength() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
G4double GetOuterRadius() const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 397 of file G4GDMLWriteSolids.cc.

399 {
400  const G4String& name = GenerateName(xtru->GetName(),xtru);
401 
402  xercesc::DOMElement* xtruElement = NewElement("xtru");
403  xtruElement->setAttributeNode(NewAttribute("name",name));
404  xtruElement->setAttributeNode(NewAttribute("lunit","mm"));
405  solElement->appendChild(xtruElement);
406 
407  const G4int NumVertex = xtru->GetNofVertices();
408 
409  for (G4int i=0;i<NumVertex;i++)
410  {
411  xercesc::DOMElement* twoDimVertexElement = NewElement("twoDimVertex");
412  xtruElement->appendChild(twoDimVertexElement);
413 
414  const G4TwoVector& vertex = xtru->GetVertex(i);
415 
416  twoDimVertexElement->setAttributeNode(NewAttribute("x",vertex.x()/mm));
417  twoDimVertexElement->setAttributeNode(NewAttribute("y",vertex.y()/mm));
418  }
419 
420  const G4int NumSection = xtru->GetNofZSections();
421 
422  for (G4int i=0;i<NumSection;i++)
423  {
424  xercesc::DOMElement* sectionElement = NewElement("section");
425  xtruElement->appendChild(sectionElement);
426 
427  const G4ExtrudedSolid::ZSection section = xtru->GetZSection(i);
428 
429  sectionElement->setAttributeNode(NewAttribute("zOrder",i));
430  sectionElement->setAttributeNode(NewAttribute("zPosition",section.fZ/mm));
431  sectionElement->
432  setAttributeNode(NewAttribute("xOffset",section.fOffset.x()/mm));
433  sectionElement->
434  setAttributeNode(NewAttribute("yOffset",section.fOffset.y()/mm));
435  sectionElement->
436  setAttributeNode(NewAttribute("scalingFactor",section.fScale));
437  }
438 }
G4String GetName() const
const XML_Char * name
Definition: expat.h:151
double y() const
double x() const
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
int G4int
Definition: G4Types.hh:78
G4TwoVector GetVertex(G4int index) const
G4int GetNofVertices() const
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
G4int GetNofZSections() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
ZSection GetZSection(G4int index) const

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 1019 of file G4GDMLWriteSolids.cc.

1021 {
1022  xercesc::DOMElement* zplaneElement = NewElement("zplane");
1023  zplaneElement->setAttributeNode(NewAttribute("z",z/mm));
1024  zplaneElement->setAttributeNode(NewAttribute("rmin",rmin/mm));
1025  zplaneElement->setAttributeNode(NewAttribute("rmax",rmax/mm));
1026  element->appendChild(zplaneElement);
1027 }
static constexpr double mm
Definition: G4SIunits.hh:115
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
tuple z
Definition: test.py:28

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

const G4int G4GDMLWriteSolids::maxTransforms = 8
staticprotected

Definition at line 149 of file G4GDMLWriteSolids.hh.

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

Definition at line 147 of file G4GDMLWriteSolids.hh.

xercesc::DOMElement* G4GDMLWriteSolids::solidsElement
protected

Definition at line 148 of file G4GDMLWriteSolids.hh.


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