Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLWriteSolids Class Reference

#include <G4GDMLWriteSolids.hh>

+ Inheritance 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)
 
void SetOutputFileOverwrite (G4bool flag)
 
virtual void DefineWrite (xercesc::DOMElement *)=0
 
virtual void MaterialsWrite (xercesc::DOMElement *)=0
 
virtual void SolidsWrite (xercesc::DOMElement *)=0
 
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)
 
void PropertyWrite (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)
 
void PropertyConstWrite (const G4String &, const G4double, const G4MaterialPropertiesTable *)
 
- 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 = nullptr
 
- Protected Attributes inherited from G4GDMLWriteMaterials
std::vector< const G4Isotope * > isotopeList
 
std::vector< const G4Element * > elementList
 
std::vector< const G4Material * > materialList
 
std::vector< const G4PhysicsOrderedFreeVector * > propertyList
 
xercesc::DOMElement * materialsElement = nullptr
 
- Protected Attributes inherited from G4GDMLWriteDefine
xercesc::DOMElement * defineElement = nullptr
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc = nullptr
 
xercesc::DOMElement * extElement = nullptr
 
xercesc::DOMElement * userinfoElement = nullptr
 
XMLCh tempStr [10000]
 
G4GDMLAuxListType auxList
 
G4bool overwriteOutputFile = false
 

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

Constructor & Destructor Documentation

◆ G4GDMLWriteSolids()

G4GDMLWriteSolids::G4GDMLWriteSolids ( )
protected

Definition at line 71 of file G4GDMLWriteSolids.cc.

◆ ~G4GDMLWriteSolids()

G4GDMLWriteSolids::~G4GDMLWriteSolids ( )
protectedvirtual

Definition at line 77 of file G4GDMLWriteSolids.cc.

78{
79}

Member Function Documentation

◆ AddSolid()

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

Definition at line 1146 of file G4GDMLWriteSolids.cc.

1147{
1148 for(std::size_t i = 0; i < solidList.size(); ++i) // Check if solid is
1149 { // already in the list!
1150 if(solidList[i] == solidPtr)
1151 {
1152 return;
1153 }
1154 }
1155
1156 solidList.push_back(solidPtr);
1157
1158 if(const G4BooleanSolid* const booleanPtr =
1159 dynamic_cast<const G4BooleanSolid*>(solidPtr))
1160 {
1161 BooleanWrite(solidsElement, booleanPtr);
1162 }
1163 else if(const G4ScaledSolid* const scaledPtr =
1164 dynamic_cast<const G4ScaledSolid*>(solidPtr))
1165 {
1166 ScaledWrite(solidsElement, scaledPtr);
1167 }
1168 else if(solidPtr->GetEntityType() == "G4MultiUnion")
1169 {
1170 const G4MultiUnion* const munionPtr =
1171 static_cast<const G4MultiUnion*>(solidPtr);
1172 MultiUnionWrite(solidsElement, munionPtr);
1173 }
1174 else if(solidPtr->GetEntityType() == "G4Box")
1175 {
1176 const G4Box* const boxPtr = static_cast<const G4Box*>(solidPtr);
1177 BoxWrite(solidsElement, boxPtr);
1178 }
1179 else if(solidPtr->GetEntityType() == "G4Cons")
1180 {
1181 const G4Cons* const conePtr = static_cast<const G4Cons*>(solidPtr);
1182 ConeWrite(solidsElement, conePtr);
1183 }
1184 else if(solidPtr->GetEntityType() == "G4EllipticalCone")
1185 {
1186 const G4EllipticalCone* const elconePtr =
1187 static_cast<const G4EllipticalCone*>(solidPtr);
1188 ElconeWrite(solidsElement, elconePtr);
1189 }
1190 else if(solidPtr->GetEntityType() == "G4Ellipsoid")
1191 {
1192 const G4Ellipsoid* const ellipsoidPtr =
1193 static_cast<const G4Ellipsoid*>(solidPtr);
1194 EllipsoidWrite(solidsElement, ellipsoidPtr);
1195 }
1196 else if(solidPtr->GetEntityType() == "G4EllipticalTube")
1197 {
1198 const G4EllipticalTube* const eltubePtr =
1199 static_cast<const G4EllipticalTube*>(solidPtr);
1200 EltubeWrite(solidsElement, eltubePtr);
1201 }
1202 else if(solidPtr->GetEntityType() == "G4ExtrudedSolid")
1203 {
1204 const G4ExtrudedSolid* const xtruPtr =
1205 static_cast<const G4ExtrudedSolid*>(solidPtr);
1206 XtruWrite(solidsElement, xtruPtr);
1207 }
1208 else if(solidPtr->GetEntityType() == "G4Hype")
1209 {
1210 const G4Hype* const hypePtr = static_cast<const G4Hype*>(solidPtr);
1211 HypeWrite(solidsElement, hypePtr);
1212 }
1213 else if(solidPtr->GetEntityType() == "G4Orb")
1214 {
1215 const G4Orb* const orbPtr = static_cast<const G4Orb*>(solidPtr);
1216 OrbWrite(solidsElement, orbPtr);
1217 }
1218 else if(solidPtr->GetEntityType() == "G4Para")
1219 {
1220 const G4Para* const paraPtr = static_cast<const G4Para*>(solidPtr);
1221 ParaWrite(solidsElement, paraPtr);
1222 }
1223 else if(solidPtr->GetEntityType() == "G4Paraboloid")
1224 {
1225 const G4Paraboloid* const paraboloidPtr =
1226 static_cast<const G4Paraboloid*>(solidPtr);
1227 ParaboloidWrite(solidsElement, paraboloidPtr);
1228 }
1229 else if(solidPtr->GetEntityType() == "G4Polycone")
1230 {
1231 const G4Polycone* const polyconePtr =
1232 static_cast<const G4Polycone*>(solidPtr);
1233 PolyconeWrite(solidsElement, polyconePtr);
1234 }
1235 else if(solidPtr->GetEntityType() == "G4GenericPolycone")
1236 {
1237 const G4GenericPolycone* const genpolyconePtr =
1238 static_cast<const G4GenericPolycone*>(solidPtr);
1239 GenericPolyconeWrite(solidsElement, genpolyconePtr);
1240 }
1241 else if(solidPtr->GetEntityType() == "G4Polyhedra")
1242 {
1243 const G4Polyhedra* const polyhedraPtr =
1244 static_cast<const G4Polyhedra*>(solidPtr);
1245 PolyhedraWrite(solidsElement, polyhedraPtr);
1246 }
1247 else if(solidPtr->GetEntityType() == "G4Sphere")
1248 {
1249 const G4Sphere* const spherePtr = static_cast<const G4Sphere*>(solidPtr);
1250 SphereWrite(solidsElement, spherePtr);
1251 }
1252 else if(solidPtr->GetEntityType() == "G4TessellatedSolid")
1253 {
1254 const G4TessellatedSolid* const tessellatedPtr =
1255 static_cast<const G4TessellatedSolid*>(solidPtr);
1256 TessellatedWrite(solidsElement, tessellatedPtr);
1257 }
1258 else if(solidPtr->GetEntityType() == "G4Tet")
1259 {
1260 const G4Tet* const tetPtr = static_cast<const G4Tet*>(solidPtr);
1261 TetWrite(solidsElement, tetPtr);
1262 }
1263 else if(solidPtr->GetEntityType() == "G4Torus")
1264 {
1265 const G4Torus* const torusPtr = static_cast<const G4Torus*>(solidPtr);
1266 TorusWrite(solidsElement, torusPtr);
1267 }
1268 else if(solidPtr->GetEntityType() == "G4GenericTrap")
1269 {
1270 const G4GenericTrap* const gtrapPtr =
1271 static_cast<const G4GenericTrap*>(solidPtr);
1272 GenTrapWrite(solidsElement, gtrapPtr);
1273 }
1274 else if(solidPtr->GetEntityType() == "G4Trap")
1275 {
1276 const G4Trap* const trapPtr = static_cast<const G4Trap*>(solidPtr);
1277 TrapWrite(solidsElement, trapPtr);
1278 }
1279 else if(solidPtr->GetEntityType() == "G4Trd")
1280 {
1281 const G4Trd* const trdPtr = static_cast<const G4Trd*>(solidPtr);
1282 TrdWrite(solidsElement, trdPtr);
1283 }
1284 else if(solidPtr->GetEntityType() == "G4Tubs")
1285 {
1286 const G4Tubs* const tubePtr = static_cast<const G4Tubs*>(solidPtr);
1287 TubeWrite(solidsElement, tubePtr);
1288 }
1289 else if(solidPtr->GetEntityType() == "G4CutTubs")
1290 {
1291 const G4CutTubs* const cuttubePtr = static_cast<const G4CutTubs*>(solidPtr);
1292 CutTubeWrite(solidsElement, cuttubePtr);
1293 }
1294 else if(solidPtr->GetEntityType() == "G4TwistedBox")
1295 {
1296 const G4TwistedBox* const twistedboxPtr =
1297 static_cast<const G4TwistedBox*>(solidPtr);
1298 TwistedboxWrite(solidsElement, twistedboxPtr);
1299 }
1300 else if(solidPtr->GetEntityType() == "G4TwistedTrap")
1301 {
1302 const G4TwistedTrap* const twistedtrapPtr =
1303 static_cast<const G4TwistedTrap*>(solidPtr);
1304 TwistedtrapWrite(solidsElement, twistedtrapPtr);
1305 }
1306 else if(solidPtr->GetEntityType() == "G4TwistedTrd")
1307 {
1308 const G4TwistedTrd* const twistedtrdPtr =
1309 static_cast<const G4TwistedTrd*>(solidPtr);
1310 TwistedtrdWrite(solidsElement, twistedtrdPtr);
1311 }
1312 else if(solidPtr->GetEntityType() == "G4TwistedTubs")
1313 {
1314 const G4TwistedTubs* const twistedtubsPtr =
1315 static_cast<const G4TwistedTubs*>(solidPtr);
1316 TwistedtubsWrite(solidsElement, twistedtubsPtr);
1317 }
1318 else
1319 {
1320 G4String error_msg = "Unknown solid: " + solidPtr->GetName() +
1321 "; Type: " + solidPtr->GetEntityType();
1322 G4Exception("G4GDMLWriteSolids::AddSolid()", "WriteError", FatalException,
1323 error_msg);
1324 }
1325}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
Definition: G4Box.hh:56
Definition: G4Cons.hh:78
void TwistedtrdWrite(xercesc::DOMElement *, const G4TwistedTrd *const)
void TorusWrite(xercesc::DOMElement *, const G4Torus *const)
void TetWrite(xercesc::DOMElement *, const G4Tet *const)
void TrapWrite(xercesc::DOMElement *, const G4Trap *const)
void HypeWrite(xercesc::DOMElement *, const G4Hype *const)
void TwistedtrapWrite(xercesc::DOMElement *, const G4TwistedTrap *const)
void ScaledWrite(xercesc::DOMElement *, const G4ScaledSolid *const)
void MultiUnionWrite(xercesc::DOMElement *solElement, const G4MultiUnion *const)
void PolyhedraWrite(xercesc::DOMElement *, const G4Polyhedra *const)
void TessellatedWrite(xercesc::DOMElement *, const G4TessellatedSolid *const)
void GenTrapWrite(xercesc::DOMElement *, const G4GenericTrap *const)
void TwistedboxWrite(xercesc::DOMElement *, const G4TwistedBox *const)
void TubeWrite(xercesc::DOMElement *, const G4Tubs *const)
void ParaboloidWrite(xercesc::DOMElement *, const G4Paraboloid *const)
std::vector< const G4VSolid * > solidList
void SphereWrite(xercesc::DOMElement *, const G4Sphere *const)
void BoxWrite(xercesc::DOMElement *, const G4Box *const)
void OrbWrite(xercesc::DOMElement *, const G4Orb *const)
void TwistedtubsWrite(xercesc::DOMElement *, const G4TwistedTubs *const)
void EltubeWrite(xercesc::DOMElement *, const G4EllipticalTube *const)
void BooleanWrite(xercesc::DOMElement *, const G4BooleanSolid *const)
void ElconeWrite(xercesc::DOMElement *, const G4EllipticalCone *const)
void ConeWrite(xercesc::DOMElement *, const G4Cons *const)
void GenericPolyconeWrite(xercesc::DOMElement *, const G4GenericPolycone *const)
xercesc::DOMElement * solidsElement
void CutTubeWrite(xercesc::DOMElement *, const G4CutTubs *const)
void EllipsoidWrite(xercesc::DOMElement *, const G4Ellipsoid *const)
void PolyconeWrite(xercesc::DOMElement *, const G4Polycone *const)
void ParaWrite(xercesc::DOMElement *, const G4Para *const)
void TrdWrite(xercesc::DOMElement *, const G4Trd *const)
void XtruWrite(xercesc::DOMElement *, const G4ExtrudedSolid *const)
Definition: G4Hype.hh:69
Definition: G4Orb.hh:56
Definition: G4Para.hh:79
Definition: G4Tet.hh:56
Definition: G4Trd.hh:63
Definition: G4Tubs.hh:75
G4String GetName() const
virtual G4GeometryType GetEntityType() const =0

Referenced by BooleanWrite(), MultiUnionWrite(), ScaledWrite(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ BooleanWrite()

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

Definition at line 141 of file G4GDMLWriteSolids.cc.

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

Referenced by AddSolid().

◆ BoxWrite()

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

Definition at line 292 of file G4GDMLWriteSolids.cc.

294{
295 const G4String& name = GenerateName(box->GetName(), box);
296
297 xercesc::DOMElement* boxElement = NewElement("box");
298 boxElement->setAttributeNode(NewAttribute("name", name));
299 boxElement->setAttributeNode(
300 NewAttribute("x", 2.0 * box->GetXHalfLength() / mm));
301 boxElement->setAttributeNode(
302 NewAttribute("y", 2.0 * box->GetYHalfLength() / mm));
303 boxElement->setAttributeNode(
304 NewAttribute("z", 2.0 * box->GetZHalfLength() / mm));
305 boxElement->setAttributeNode(NewAttribute("lunit", "mm"));
306 solElement->appendChild(boxElement);
307}
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const

Referenced by AddSolid().

◆ ConeWrite()

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

Definition at line 310 of file G4GDMLWriteSolids.cc.

312{
313 const G4String& name = GenerateName(cone->GetName(), cone);
314
315 xercesc::DOMElement* coneElement = NewElement("cone");
316 coneElement->setAttributeNode(NewAttribute("name", name));
317 coneElement->setAttributeNode(
318 NewAttribute("rmin1", cone->GetInnerRadiusMinusZ() / mm));
319 coneElement->setAttributeNode(
320 NewAttribute("rmax1", cone->GetOuterRadiusMinusZ() / mm));
321 coneElement->setAttributeNode(
322 NewAttribute("rmin2", cone->GetInnerRadiusPlusZ() / mm));
323 coneElement->setAttributeNode(
324 NewAttribute("rmax2", cone->GetOuterRadiusPlusZ() / mm));
325 coneElement->setAttributeNode(
326 NewAttribute("z", 2.0 * cone->GetZHalfLength() / mm));
327 coneElement->setAttributeNode(
328 NewAttribute("startphi", cone->GetStartPhiAngle() / degree));
329 coneElement->setAttributeNode(
330 NewAttribute("deltaphi", cone->GetDeltaPhiAngle() / degree));
331 coneElement->setAttributeNode(NewAttribute("aunit", "deg"));
332 coneElement->setAttributeNode(NewAttribute("lunit", "mm"));
333 solElement->appendChild(coneElement);
334}
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetOuterRadiusMinusZ() const
G4double GetZHalfLength() const

Referenced by AddSolid().

◆ CutTubeWrite()

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

Definition at line 894 of file G4GDMLWriteSolids.cc.

896{
897 const G4String& name = GenerateName(cuttube->GetName(), cuttube);
898
899 xercesc::DOMElement* cuttubeElement = NewElement("cutTube");
900 cuttubeElement->setAttributeNode(NewAttribute("name", name));
901 cuttubeElement->setAttributeNode(
902 NewAttribute("rmin", cuttube->GetInnerRadius() / mm));
903 cuttubeElement->setAttributeNode(
904 NewAttribute("rmax", cuttube->GetOuterRadius() / mm));
905 cuttubeElement->setAttributeNode(
906 NewAttribute("z", 2.0 * cuttube->GetZHalfLength() / mm));
907 cuttubeElement->setAttributeNode(
908 NewAttribute("startphi", cuttube->GetStartPhiAngle() / degree));
909 cuttubeElement->setAttributeNode(
910 NewAttribute("deltaphi", cuttube->GetDeltaPhiAngle() / degree));
911 cuttubeElement->setAttributeNode(
912 NewAttribute("lowX", cuttube->GetLowNorm().getX() / mm));
913 cuttubeElement->setAttributeNode(
914 NewAttribute("lowY", cuttube->GetLowNorm().getY() / mm));
915 cuttubeElement->setAttributeNode(
916 NewAttribute("lowZ", cuttube->GetLowNorm().getZ() / mm));
917 cuttubeElement->setAttributeNode(
918 NewAttribute("highX", cuttube->GetHighNorm().getX() / mm));
919 cuttubeElement->setAttributeNode(
920 NewAttribute("highY", cuttube->GetHighNorm().getY() / mm));
921 cuttubeElement->setAttributeNode(
922 NewAttribute("highZ", cuttube->GetHighNorm().getZ() / mm));
923 cuttubeElement->setAttributeNode(NewAttribute("aunit", "deg"));
924 cuttubeElement->setAttributeNode(NewAttribute("lunit", "mm"));
925 solElement->appendChild(cuttubeElement);
926}
double getZ() const
double getX() const
double getY() const
G4ThreeVector GetHighNorm() const
G4double GetStartPhiAngle() const
G4double GetZHalfLength() const
G4double GetInnerRadius() const
G4ThreeVector GetLowNorm() const
G4double GetDeltaPhiAngle() const
G4double GetOuterRadius() const

Referenced by AddSolid().

◆ ElconeWrite()

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

Definition at line 337 of file G4GDMLWriteSolids.cc.

339{
340 const G4String& name = GenerateName(elcone->GetName(), elcone);
341
342 xercesc::DOMElement* elconeElement = NewElement("elcone");
343 elconeElement->setAttributeNode(NewAttribute("name", name));
344 elconeElement->setAttributeNode(
345 NewAttribute("dx", elcone->GetSemiAxisX() / mm));
346 elconeElement->setAttributeNode(
347 NewAttribute("dy", elcone->GetSemiAxisY() / mm));
348 elconeElement->setAttributeNode(NewAttribute("zmax", elcone->GetZMax() / mm));
349 elconeElement->setAttributeNode(
350 NewAttribute("zcut", elcone->GetZTopCut() / mm));
351 elconeElement->setAttributeNode(NewAttribute("lunit", "mm"));
352 solElement->appendChild(elconeElement);
353}
G4double GetSemiAxisX() const
G4double GetSemiAxisY() const
G4double GetZMax() const
G4double GetZTopCut() const

Referenced by AddSolid().

◆ EllipsoidWrite()

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

Definition at line 356 of file G4GDMLWriteSolids.cc.

358{
359 const G4String& name = GenerateName(ellipsoid->GetName(), ellipsoid);
360
361 xercesc::DOMElement* ellipsoidElement = NewElement("ellipsoid");
362 ellipsoidElement->setAttributeNode(NewAttribute("name", name));
363 ellipsoidElement->setAttributeNode(
364 NewAttribute("ax", ellipsoid->GetSemiAxisMax(0) / mm));
365 ellipsoidElement->setAttributeNode(
366 NewAttribute("by", ellipsoid->GetSemiAxisMax(1) / mm));
367 ellipsoidElement->setAttributeNode(
368 NewAttribute("cz", ellipsoid->GetSemiAxisMax(2) / mm));
369 ellipsoidElement->setAttributeNode(
370 NewAttribute("zcut1", ellipsoid->GetZBottomCut() / mm));
371 ellipsoidElement->setAttributeNode(
372 NewAttribute("zcut2", ellipsoid->GetZTopCut() / mm));
373 ellipsoidElement->setAttributeNode(NewAttribute("lunit", "mm"));
374 solElement->appendChild(ellipsoidElement);
375}
G4double GetSemiAxisMax(G4int i) const
G4double GetZTopCut() const
G4double GetZBottomCut() const

Referenced by AddSolid().

◆ EltubeWrite()

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

Definition at line 378 of file G4GDMLWriteSolids.cc.

380{
381 const G4String& name = GenerateName(eltube->GetName(), eltube);
382
383 xercesc::DOMElement* eltubeElement = NewElement("eltube");
384 eltubeElement->setAttributeNode(NewAttribute("name", name));
385 eltubeElement->setAttributeNode(NewAttribute("dx", eltube->GetDx() / mm));
386 eltubeElement->setAttributeNode(NewAttribute("dy", eltube->GetDy() / mm));
387 eltubeElement->setAttributeNode(NewAttribute("dz", eltube->GetDz() / mm));
388 eltubeElement->setAttributeNode(NewAttribute("lunit", "mm"));
389 solElement->appendChild(eltubeElement);
390}
G4double GetDy() const
G4double GetDx() const
G4double GetDz() const

Referenced by AddSolid().

◆ GenericPolyconeWrite()

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(
555 NewAttribute("startphi", startPhi / degree));
556 polyconeElement->setAttributeNode(
557 NewAttribute("deltaphi", (polycone->GetEndPhi() - startPhi) / degree));
558 polyconeElement->setAttributeNode(NewAttribute("aunit", "deg"));
559 polyconeElement->setAttributeNode(NewAttribute("lunit", "mm"));
560 solElement->appendChild(polyconeElement);
561
562 const std::size_t num_rzpoints = polycone->GetNumRZCorner();
563 for(std::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}
double G4double
Definition: G4Types.hh:83
void RZPointWrite(xercesc::DOMElement *, const G4double &, const G4double &)
G4double GetStartPhi() const
G4double GetEndPhi() const
G4int GetNumRZCorner() const
G4PolyconeSideRZ GetCorner(G4int index) const

Referenced by AddSolid().

◆ GenTrapWrite()

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

Definition at line 780 of file G4GDMLWriteSolids.cc.

782{
783 const G4String& name = GenerateName(gtrap->GetName(), gtrap);
784
785 std::vector<G4TwoVector> vertices = gtrap->GetVertices();
786
787 xercesc::DOMElement* gtrapElement = NewElement("arb8");
788 gtrapElement->setAttributeNode(NewAttribute("name", name));
789 gtrapElement->setAttributeNode(
790 NewAttribute("dz", gtrap->GetZHalfLength() / mm));
791 gtrapElement->setAttributeNode(NewAttribute("v1x", vertices[0].x()));
792 gtrapElement->setAttributeNode(NewAttribute("v1y", vertices[0].y()));
793 gtrapElement->setAttributeNode(NewAttribute("v2x", vertices[1].x()));
794 gtrapElement->setAttributeNode(NewAttribute("v2y", vertices[1].y()));
795 gtrapElement->setAttributeNode(NewAttribute("v3x", vertices[2].x()));
796 gtrapElement->setAttributeNode(NewAttribute("v3y", vertices[2].y()));
797 gtrapElement->setAttributeNode(NewAttribute("v4x", vertices[3].x()));
798 gtrapElement->setAttributeNode(NewAttribute("v4y", vertices[3].y()));
799 gtrapElement->setAttributeNode(NewAttribute("v5x", vertices[4].x()));
800 gtrapElement->setAttributeNode(NewAttribute("v5y", vertices[4].y()));
801 gtrapElement->setAttributeNode(NewAttribute("v6x", vertices[5].x()));
802 gtrapElement->setAttributeNode(NewAttribute("v6y", vertices[5].y()));
803 gtrapElement->setAttributeNode(NewAttribute("v7x", vertices[6].x()));
804 gtrapElement->setAttributeNode(NewAttribute("v7y", vertices[6].y()));
805 gtrapElement->setAttributeNode(NewAttribute("v8x", vertices[7].x()));
806 gtrapElement->setAttributeNode(NewAttribute("v8y", vertices[7].y()));
807 gtrapElement->setAttributeNode(NewAttribute("lunit", "mm"));
808 solElement->appendChild(gtrapElement);
809}
G4double GetZHalfLength() const
const std::vector< G4TwoVector > & GetVertices() const

Referenced by AddSolid().

◆ HypeWrite()

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

Definition at line 438 of file G4GDMLWriteSolids.cc.

440{
441 const G4String& name = GenerateName(hype->GetName(), hype);
442
443 xercesc::DOMElement* hypeElement = NewElement("hype");
444 hypeElement->setAttributeNode(NewAttribute("name", name));
445 hypeElement->setAttributeNode(
446 NewAttribute("rmin", hype->GetInnerRadius() / mm));
447 hypeElement->setAttributeNode(
448 NewAttribute("rmax", hype->GetOuterRadius() / mm));
449 hypeElement->setAttributeNode(
450 NewAttribute("inst", hype->GetInnerStereo() / degree));
451 hypeElement->setAttributeNode(
452 NewAttribute("outst", hype->GetOuterStereo() / degree));
453 hypeElement->setAttributeNode(
454 NewAttribute("z", 2.0 * hype->GetZHalfLength() / mm));
455 hypeElement->setAttributeNode(NewAttribute("aunit", "deg"));
456 hypeElement->setAttributeNode(NewAttribute("lunit", "mm"));
457 solElement->appendChild(hypeElement);
458}
G4double GetInnerStereo() const
G4double GetZHalfLength() const
G4double GetOuterStereo() const
G4double GetOuterRadius() const
G4double GetInnerRadius() const

Referenced by AddSolid().

◆ MultiUnionWrite()

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

Definition at line 82 of file G4GDMLWriteSolids.cc.

84{
85 G4int numSolids = munionSolid->GetNumberOfSolids();
86 G4String tag("multiUnion");
87
88 G4VSolid* solid;
89 G4Transform3D transform;
90
91 const G4String& name = GenerateName(munionSolid->GetName(), munionSolid);
92 xercesc::DOMElement* multiUnionElement = NewElement(tag);
93 multiUnionElement->setAttributeNode(NewAttribute("name", name));
94
95 for(G4int i = 0; i < numSolids; ++i)
96 {
97 solid = munionSolid->GetSolid(i);
98 transform = munionSolid->GetTransformation(i);
99
100 HepGeom::Rotate3D rot3d;
102 HepGeom::Scale3D scale;
103 transform.getDecomposition(scale, rot3d, transl);
104
106 G4RotationMatrix rotm(CLHEP::HepRep3x3(rot3d.xx(), rot3d.xy(), rot3d.xz(),
107 rot3d.yx(), rot3d.yy(), rot3d.yz(),
108 rot3d.zx(), rot3d.zy(), rot3d.zz()));
109 G4ThreeVector rot = GetAngles(rotm);
110
111 AddSolid(solid);
112 const G4String& solidref = GenerateName(solid->GetName(), solid);
113 std::ostringstream os;
114 os << i + 1;
115 const G4String& nodeName = "Node-" + G4String(os.str());
116 xercesc::DOMElement* solidElement = NewElement("solid");
117 solidElement->setAttributeNode(NewAttribute("ref", solidref));
118 xercesc::DOMElement* multiUnionNodeElement = NewElement("multiUnionNode");
119 multiUnionNodeElement->setAttributeNode(NewAttribute("name", nodeName));
120 multiUnionNodeElement->appendChild(solidElement); // Append solid to node
121 if((std::fabs(pos.x()) > kLinearPrecision) ||
122 (std::fabs(pos.y()) > kLinearPrecision) ||
123 (std::fabs(pos.z()) > kLinearPrecision))
124 {
125 PositionWrite(multiUnionNodeElement, name + "_pos", pos);
126 }
127 if((std::fabs(rot.x()) > kAngularPrecision) ||
128 (std::fabs(rot.y()) > kAngularPrecision) ||
129 (std::fabs(rot.z()) > kAngularPrecision))
130 {
131 RotationWrite(multiUnionNodeElement, name + "_rot", rot);
132 }
133 multiUnionElement->appendChild(multiUnionNodeElement); // Append node
134 }
135
136 solElement->appendChild(multiUnionElement);
137 // Add the multiUnion solid AFTER the constituent nodes!
138}
const G4Transform3D & GetTransformation(G4int index) const
G4int GetNumberOfSolids() const
G4VSolid * GetSolid(G4int index) const
double zz() const
Definition: Transform3D.h:281
double yz() const
Definition: Transform3D.h:272
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:173
double xy() const
Definition: Transform3D.h:260
CLHEP::Hep3Vector getTranslation() const
double zx() const
Definition: Transform3D.h:275
double yx() const
Definition: Transform3D.h:266
double zy() const
Definition: Transform3D.h:278
double xx() const
Definition: Transform3D.h:257
double yy() const
Definition: Transform3D.h:269
double xz() const
Definition: Transform3D.h:263

Referenced by AddSolid().

◆ OpticalSurfaceWrite()

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

Definition at line 1058 of file G4GDMLWriteSolids.cc.

1060{
1061 xercesc::DOMElement* optElement = NewElement("opticalsurface");
1062 G4OpticalSurfaceModel smodel = surf->GetModel();
1063 G4double sval =
1064 (smodel == glisur) ? surf->GetPolish() : surf->GetSigmaAlpha();
1065 const G4String& name = GenerateName(surf->GetName(), surf);
1066
1067 optElement->setAttributeNode(NewAttribute("name", name));
1068 optElement->setAttributeNode(NewAttribute("model", smodel));
1069 optElement->setAttributeNode(NewAttribute("finish", surf->GetFinish()));
1070 optElement->setAttributeNode(NewAttribute("type", surf->GetType()));
1071 optElement->setAttributeNode(NewAttribute("value", sval));
1072
1073 // Write any property attached to the optical surface...
1074 //
1075 if(surf->GetMaterialPropertiesTable())
1076 {
1077 PropertyWrite(optElement, surf);
1078 }
1079
1080 solElement->appendChild(optElement);
1081}
G4OpticalSurfaceModel
@ glisur
void PropertyWrite(xercesc::DOMElement *, const G4OpticalSurface *const)
G4OpticalSurfaceModel GetModel() const
G4double GetSigmaAlpha() const
G4OpticalSurfaceFinish GetFinish() const
G4double GetPolish() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
const G4String & GetName() const
const G4SurfaceType & GetType() const

Referenced by G4GDMLWriteStructure::BorderSurfaceCache(), and G4GDMLWriteStructure::SkinSurfaceCache().

◆ OrbWrite()

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

Definition at line 461 of file G4GDMLWriteSolids.cc.

463{
464 const G4String& name = GenerateName(orb->GetName(), orb);
465
466 xercesc::DOMElement* orbElement = NewElement("orb");
467 orbElement->setAttributeNode(NewAttribute("name", name));
468 orbElement->setAttributeNode(NewAttribute("r", orb->GetRadius() / mm));
469 orbElement->setAttributeNode(NewAttribute("lunit", "mm"));
470 solElement->appendChild(orbElement);
471}
G4double GetRadius() const

Referenced by AddSolid().

◆ ParaboloidWrite()

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(
509 NewAttribute("rlo", paraboloid->GetRadiusMinusZ() / mm));
510 paraboloidElement->setAttributeNode(
511 NewAttribute("rhi", paraboloid->GetRadiusPlusZ() / mm));
512 paraboloidElement->setAttributeNode(
513 NewAttribute("dz", paraboloid->GetZHalfLength() / mm));
514 paraboloidElement->setAttributeNode(NewAttribute("lunit", "mm"));
515 solElement->appendChild(paraboloidElement);
516}
G4double GetRadiusPlusZ() const
G4double GetRadiusMinusZ() const
G4double GetZHalfLength() const

Referenced by AddSolid().

◆ ParaWrite()

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

Definition at line 474 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(
487 NewAttribute("x", 2.0 * para->GetXHalfLength() / mm));
488 paraElement->setAttributeNode(
489 NewAttribute("y", 2.0 * para->GetYHalfLength() / mm));
490 paraElement->setAttributeNode(
491 NewAttribute("z", 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}
double phi() const
double theta() const
G4double GetTanAlpha() const
G4ThreeVector GetSymAxis() const
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const

Referenced by AddSolid().

◆ PolyconeWrite()

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

Definition at line 519 of file G4GDMLWriteSolids.cc.

521{
522 const G4String& name = GenerateName(polycone->GetName(), polycone);
523
524 xercesc::DOMElement* polyconeElement = NewElement("polycone");
525 polyconeElement->setAttributeNode(NewAttribute("name", name));
526 polyconeElement->setAttributeNode(NewAttribute(
527 "startphi", polycone->GetOriginalParameters()->Start_angle / degree));
528 polyconeElement->setAttributeNode(NewAttribute(
529 "deltaphi", polycone->GetOriginalParameters()->Opening_angle / degree));
530 polyconeElement->setAttributeNode(NewAttribute("aunit", "deg"));
531 polyconeElement->setAttributeNode(NewAttribute("lunit", "mm"));
532 solElement->appendChild(polyconeElement);
533
534 const std::size_t num_zplanes
535 = polycone->GetOriginalParameters()->Num_z_planes;
536 const G4double* z_array = polycone->GetOriginalParameters()->Z_values;
537 const G4double* rmin_array = polycone->GetOriginalParameters()->Rmin;
538 const G4double* rmax_array = polycone->GetOriginalParameters()->Rmax;
539
540 for(std::size_t i = 0; i < num_zplanes; ++i)
541 {
542 ZplaneWrite(polyconeElement, z_array[i], rmin_array[i], rmax_array[i]);
543 }
544}
void ZplaneWrite(xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
G4PolyconeHistorical * GetOriginalParameters() const

Referenced by AddSolid().

◆ PolyhedraWrite()

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

Definition at line 572 of file G4GDMLWriteSolids.cc.

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

Referenced by AddSolid().

◆ PropertyWrite()

void G4GDMLWriteSolids::PropertyWrite ( xercesc::DOMElement *  optElement,
const G4OpticalSurface * const  surf 
)
protected

Definition at line 1084 of file G4GDMLWriteSolids.cc.

1086{
1087 xercesc::DOMElement* propElement;
1089 auto pmap = ptable->GetPropertyMap();
1090 auto cmap = ptable->GetConstPropertyMap();
1091
1092 for(auto mpos = pmap->cbegin(); mpos != pmap->cend(); ++mpos)
1093 {
1094 propElement = NewElement("property");
1095 propElement->setAttributeNode(
1096 NewAttribute("name", ptable->GetMaterialPropertyNames()[mpos->first]));
1097 propElement->setAttributeNode(NewAttribute(
1098 "ref", GenerateName(ptable->GetMaterialPropertyNames()[mpos->first],
1099 mpos->second)));
1100 if(mpos->second)
1101 {
1102 PropertyVectorWrite(ptable->GetMaterialPropertyNames()[mpos->first],
1103 mpos->second);
1104 optElement->appendChild(propElement);
1105 }
1106 else
1107 {
1108 G4String warn_message = "Null pointer for material property -" +
1109 ptable->GetMaterialPropertyNames()[mpos->first] +
1110 "- of optical surface -" + surf->GetName() +
1111 "- !";
1112 G4Exception("G4GDMLWriteSolids::PropertyWrite()", "NullPointer",
1113 JustWarning, warn_message);
1114 continue;
1115 }
1116 }
1117 for(auto cpos = cmap->cbegin(); cpos != cmap->cend(); ++cpos)
1118 {
1119 propElement = NewElement("property");
1120 propElement->setAttributeNode(NewAttribute(
1121 "name", ptable->GetMaterialConstPropertyNames()[cpos->first]));
1122 propElement->setAttributeNode(NewAttribute(
1123 "ref", ptable->GetMaterialConstPropertyNames()[cpos->first]));
1124 xercesc::DOMElement* constElement = NewElement("constant");
1125 constElement->setAttributeNode(NewAttribute(
1126 "name", ptable->GetMaterialConstPropertyNames()[cpos->first]));
1127 constElement->setAttributeNode(NewAttribute("value", cpos->second));
1128 defineElement->appendChild(constElement);
1129 optElement->appendChild(propElement);
1130 }
1131}
@ JustWarning
xercesc::DOMElement * defineElement
void PropertyVectorWrite(const G4String &, const G4PhysicsOrderedFreeVector *const)
std::vector< G4String > GetMaterialPropertyNames() const
std::vector< G4String > GetMaterialConstPropertyNames() const
const std::map< G4int, G4double, std::less< G4int > > * GetConstPropertyMap() const
const std::map< G4int, G4MaterialPropertyVector *, std::less< G4int > > * GetPropertyMap() const

Referenced by OpticalSurfaceWrite().

◆ RZPointWrite()

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

Definition at line 1048 of file G4GDMLWriteSolids.cc.

1050{
1051 xercesc::DOMElement* rzpointElement = NewElement("rzpoint");
1052 rzpointElement->setAttributeNode(NewAttribute("r", r / mm));
1053 rzpointElement->setAttributeNode(NewAttribute("z", z / mm));
1054 element->appendChild(rzpointElement);
1055}

Referenced by GenericPolyconeWrite(), and PolyhedraWrite().

◆ ScaledWrite()

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

Definition at line 259 of file G4GDMLWriteSolids.cc.

261{
262 G4String tag("scaledSolid");
263
264 G4VSolid* solid = const_cast<G4VSolid*>(scaled->GetUnscaledSolid());
265 G4Scale3D scale = scaled->GetScaleTransform();
266 G4ThreeVector sclVector = G4ThreeVector(scale.xx(), scale.yy(), scale.zz());
267
268 AddSolid(solid); // Add the constituent solid!
269
270 const G4String& name = GenerateName(scaled->GetName(), scaled);
271 const G4String& solidref = GenerateName(solid->GetName(), solid);
272
273 xercesc::DOMElement* scaledElement = NewElement(tag);
274 scaledElement->setAttributeNode(NewAttribute("name", name));
275
276 xercesc::DOMElement* solidElement = NewElement("solidref");
277 solidElement->setAttributeNode(NewAttribute("ref", solidref));
278 scaledElement->appendChild(solidElement);
279
280 if((std::fabs(scale.xx()) > kLinearPrecision) &&
281 (std::fabs(scale.yy()) > kLinearPrecision) &&
282 (std::fabs(scale.zz()) > kLinearPrecision))
283 {
284 ScaleWrite(scaledElement, name + "_scl", sclVector);
285 }
286
287 solElement->appendChild(scaledElement);
288 // Add the scaled solid AFTER its constituent solid!
289}
CLHEP::Hep3Vector G4ThreeVector
void ScaleWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
G4VSolid * GetUnscaledSolid() const
G4Scale3D GetScaleTransform() const

Referenced by AddSolid().

◆ SolidsWrite()

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

Implements G4GDMLWrite.

Definition at line 1134 of file G4GDMLWriteSolids.cc.

1135{
1136#ifdef G4VERBOSE
1137 G4cout << "G4GDML: Writing solids..." << G4endl;
1138#endif
1139 solidsElement = NewElement("solids");
1140 gdmlElement->appendChild(solidsElement);
1141
1142 solidList.clear();
1143}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ SphereWrite()

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

Definition at line 632 of file G4GDMLWriteSolids.cc.

634{
635 const G4String& name = GenerateName(sphere->GetName(), sphere);
636
637 xercesc::DOMElement* sphereElement = NewElement("sphere");
638 sphereElement->setAttributeNode(NewAttribute("name", name));
639 sphereElement->setAttributeNode(
640 NewAttribute("rmin", sphere->GetInnerRadius() / mm));
641 sphereElement->setAttributeNode(
642 NewAttribute("rmax", sphere->GetOuterRadius() / mm));
643 sphereElement->setAttributeNode(
644 NewAttribute("startphi", sphere->GetStartPhiAngle() / degree));
645 sphereElement->setAttributeNode(
646 NewAttribute("deltaphi", sphere->GetDeltaPhiAngle() / degree));
647 sphereElement->setAttributeNode(
648 NewAttribute("starttheta", sphere->GetStartThetaAngle() / degree));
649 sphereElement->setAttributeNode(
650 NewAttribute("deltatheta", sphere->GetDeltaThetaAngle() / degree));
651 sphereElement->setAttributeNode(NewAttribute("aunit", "deg"));
652 sphereElement->setAttributeNode(NewAttribute("lunit", "mm"));
653 solElement->appendChild(sphereElement);
654}
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetDeltaThetaAngle() const
G4double GetStartThetaAngle() const

Referenced by AddSolid().

◆ TessellatedWrite()

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

Definition at line 657 of file G4GDMLWriteSolids.cc.

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

Referenced by AddSolid().

◆ TetWrite()

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

Definition at line 736 of file G4GDMLWriteSolids.cc.

738{
739 const G4String& solid_name = tet->GetName();
740 const G4String& name = GenerateName(solid_name, tet);
741
742 std::vector<G4ThreeVector> vertexList = tet->GetVertices();
743
744 xercesc::DOMElement* tetElement = NewElement("tet");
745 tetElement->setAttributeNode(NewAttribute("name", name));
746 tetElement->setAttributeNode(NewAttribute("vertex1", solid_name + "_v1"));
747 tetElement->setAttributeNode(NewAttribute("vertex2", solid_name + "_v2"));
748 tetElement->setAttributeNode(NewAttribute("vertex3", solid_name + "_v3"));
749 tetElement->setAttributeNode(NewAttribute("vertex4", solid_name + "_v4"));
750 tetElement->setAttributeNode(NewAttribute("lunit", "mm"));
751 solElement->appendChild(tetElement);
752
753 AddPosition(solid_name + "_v1", vertexList[0]);
754 AddPosition(solid_name + "_v2", vertexList[1]);
755 AddPosition(solid_name + "_v3", vertexList[2]);
756 AddPosition(solid_name + "_v4", vertexList[3]);
757}
void GetVertices(G4ThreeVector &anchor, G4ThreeVector &p1, G4ThreeVector &p2, G4ThreeVector &p3) const
Definition: G4Tet.cc:282

Referenced by AddSolid().

◆ TorusWrite()

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

Definition at line 760 of file G4GDMLWriteSolids.cc.

762{
763 const G4String& name = GenerateName(torus->GetName(), torus);
764
765 xercesc::DOMElement* torusElement = NewElement("torus");
766 torusElement->setAttributeNode(NewAttribute("name", name));
767 torusElement->setAttributeNode(NewAttribute("rmin", torus->GetRmin() / mm));
768 torusElement->setAttributeNode(NewAttribute("rmax", torus->GetRmax() / mm));
769 torusElement->setAttributeNode(NewAttribute("rtor", torus->GetRtor() / mm));
770 torusElement->setAttributeNode(
771 NewAttribute("startphi", torus->GetSPhi() / degree));
772 torusElement->setAttributeNode(
773 NewAttribute("deltaphi", torus->GetDPhi() / degree));
774 torusElement->setAttributeNode(NewAttribute("aunit", "deg"));
775 torusElement->setAttributeNode(NewAttribute("lunit", "mm"));
776 solElement->appendChild(torusElement);
777}
G4double GetDPhi() const
G4double GetRmin() const
G4double GetRtor() const
G4double GetRmax() const
G4double GetSPhi() const

Referenced by AddSolid().

◆ TrapWrite()

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

Definition at line 812 of file G4GDMLWriteSolids.cc.

814{
815 const G4String& name = GenerateName(trap->GetName(), trap);
816
817 const G4ThreeVector& simaxis = trap->GetSymAxis();
818 const G4double phi = simaxis.phi();
819 const G4double theta = simaxis.theta();
820 const G4double alpha1 = std::atan(trap->GetTanAlpha1());
821 const G4double alpha2 = std::atan(trap->GetTanAlpha2());
822
823 xercesc::DOMElement* trapElement = NewElement("trap");
824 trapElement->setAttributeNode(NewAttribute("name", name));
825 trapElement->setAttributeNode(
826 NewAttribute("z", 2.0 * trap->GetZHalfLength() / mm));
827 trapElement->setAttributeNode(NewAttribute("theta", theta / degree));
828 trapElement->setAttributeNode(NewAttribute("phi", phi / degree));
829 trapElement->setAttributeNode(
830 NewAttribute("y1", 2.0 * trap->GetYHalfLength1() / mm));
831 trapElement->setAttributeNode(
832 NewAttribute("x1", 2.0 * trap->GetXHalfLength1() / mm));
833 trapElement->setAttributeNode(
834 NewAttribute("x2", 2.0 * trap->GetXHalfLength2() / mm));
835 trapElement->setAttributeNode(NewAttribute("alpha1", alpha1 / degree));
836 trapElement->setAttributeNode(
837 NewAttribute("y2", 2.0 * trap->GetYHalfLength2() / mm));
838 trapElement->setAttributeNode(
839 NewAttribute("x3", 2.0 * trap->GetXHalfLength3() / mm));
840 trapElement->setAttributeNode(
841 NewAttribute("x4", 2.0 * trap->GetXHalfLength4() / mm));
842 trapElement->setAttributeNode(NewAttribute("alpha2", alpha2 / degree));
843 trapElement->setAttributeNode(NewAttribute("aunit", "deg"));
844 trapElement->setAttributeNode(NewAttribute("lunit", "mm"));
845 solElement->appendChild(trapElement);
846}
G4double GetYHalfLength1() const
G4double GetTanAlpha2() const
G4double GetXHalfLength2() const
G4ThreeVector GetSymAxis() const
G4double GetXHalfLength4() const
G4double GetZHalfLength() const
G4double GetYHalfLength2() const
G4double GetTanAlpha1() const
G4double GetXHalfLength3() const
G4double GetXHalfLength1() const

Referenced by AddSolid().

◆ TrdWrite()

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

Definition at line 849 of file G4GDMLWriteSolids.cc.

851{
852 const G4String& name = GenerateName(trd->GetName(), trd);
853
854 xercesc::DOMElement* trdElement = NewElement("trd");
855 trdElement->setAttributeNode(NewAttribute("name", name));
856 trdElement->setAttributeNode(
857 NewAttribute("x1", 2.0 * trd->GetXHalfLength1() / mm));
858 trdElement->setAttributeNode(
859 NewAttribute("x2", 2.0 * trd->GetXHalfLength2() / mm));
860 trdElement->setAttributeNode(
861 NewAttribute("y1", 2.0 * trd->GetYHalfLength1() / mm));
862 trdElement->setAttributeNode(
863 NewAttribute("y2", 2.0 * trd->GetYHalfLength2() / mm));
864 trdElement->setAttributeNode(
865 NewAttribute("z", 2.0 * trd->GetZHalfLength() / mm));
866 trdElement->setAttributeNode(NewAttribute("lunit", "mm"));
867 solElement->appendChild(trdElement);
868}
G4double GetXHalfLength2() const
G4double GetYHalfLength2() const
G4double GetXHalfLength1() const
G4double GetYHalfLength1() const
G4double GetZHalfLength() const

Referenced by AddSolid().

◆ TubeWrite()

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

Definition at line 871 of file G4GDMLWriteSolids.cc.

873{
874 const G4String& name = GenerateName(tube->GetName(), tube);
875
876 xercesc::DOMElement* tubeElement = NewElement("tube");
877 tubeElement->setAttributeNode(NewAttribute("name", name));
878 tubeElement->setAttributeNode(
879 NewAttribute("rmin", tube->GetInnerRadius() / mm));
880 tubeElement->setAttributeNode(
881 NewAttribute("rmax", tube->GetOuterRadius() / mm));
882 tubeElement->setAttributeNode(
883 NewAttribute("z", 2.0 * tube->GetZHalfLength() / mm));
884 tubeElement->setAttributeNode(
885 NewAttribute("startphi", tube->GetStartPhiAngle() / degree));
886 tubeElement->setAttributeNode(
887 NewAttribute("deltaphi", tube->GetDeltaPhiAngle() / degree));
888 tubeElement->setAttributeNode(NewAttribute("aunit", "deg"));
889 tubeElement->setAttributeNode(NewAttribute("lunit", "mm"));
890 solElement->appendChild(tubeElement);
891}
G4double GetZHalfLength() const
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const

Referenced by AddSolid().

◆ TwistedboxWrite()

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

Definition at line 929 of file G4GDMLWriteSolids.cc.

931{
932 const G4String& name = GenerateName(twistedbox->GetName(), twistedbox);
933
934 xercesc::DOMElement* twistedboxElement = NewElement("twistedbox");
935 twistedboxElement->setAttributeNode(NewAttribute("name", name));
936 twistedboxElement->setAttributeNode(
937 NewAttribute("x", 2.0 * twistedbox->GetXHalfLength() / mm));
938 twistedboxElement->setAttributeNode(
939 NewAttribute("y", 2.0 * twistedbox->GetYHalfLength() / mm));
940 twistedboxElement->setAttributeNode(
941 NewAttribute("z", 2.0 * twistedbox->GetZHalfLength() / mm));
942 twistedboxElement->setAttributeNode(
943 NewAttribute("PhiTwist", twistedbox->GetPhiTwist() / degree));
944 twistedboxElement->setAttributeNode(NewAttribute("aunit", "deg"));
945 twistedboxElement->setAttributeNode(NewAttribute("lunit", "mm"));
946 solElement->appendChild(twistedboxElement);
947}
G4double GetPhiTwist() const
Definition: G4TwistedBox.hh:65
G4double GetXHalfLength() const
Definition: G4TwistedBox.hh:62
G4double GetZHalfLength() const
Definition: G4TwistedBox.hh:64
G4double GetYHalfLength() const
Definition: G4TwistedBox.hh:63

Referenced by AddSolid().

◆ TwistedtrapWrite()

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

Definition at line 950 of file G4GDMLWriteSolids.cc.

952{
953 const G4String& name = GenerateName(twistedtrap->GetName(), twistedtrap);
954
955 xercesc::DOMElement* twistedtrapElement = NewElement("twistedtrap");
956 twistedtrapElement->setAttributeNode(NewAttribute("name", name));
957 twistedtrapElement->setAttributeNode(
958 NewAttribute("y1", 2.0 * twistedtrap->GetY1HalfLength() / mm));
959 twistedtrapElement->setAttributeNode(
960 NewAttribute("x1", 2.0 * twistedtrap->GetX1HalfLength() / mm));
961 twistedtrapElement->setAttributeNode(
962 NewAttribute("x2", 2.0 * twistedtrap->GetX2HalfLength() / mm));
963 twistedtrapElement->setAttributeNode(
964 NewAttribute("y2", 2.0 * twistedtrap->GetY2HalfLength() / mm));
965 twistedtrapElement->setAttributeNode(
966 NewAttribute("x3", 2.0 * twistedtrap->GetX3HalfLength() / mm));
967 twistedtrapElement->setAttributeNode(
968 NewAttribute("x4", 2.0 * twistedtrap->GetX4HalfLength() / mm));
969 twistedtrapElement->setAttributeNode(
970 NewAttribute("z", 2.0 * twistedtrap->GetZHalfLength() / mm));
971 twistedtrapElement->setAttributeNode(
972 NewAttribute("Alph", twistedtrap->GetTiltAngleAlpha() / degree));
973 twistedtrapElement->setAttributeNode(
974 NewAttribute("Theta", twistedtrap->GetPolarAngleTheta() / degree));
975 twistedtrapElement->setAttributeNode(
976 NewAttribute("Phi", twistedtrap->GetAzimuthalAnglePhi() / degree));
977 twistedtrapElement->setAttributeNode(
978 NewAttribute("PhiTwist", twistedtrap->GetPhiTwist() / degree));
979 twistedtrapElement->setAttributeNode(NewAttribute("aunit", "deg"));
980 twistedtrapElement->setAttributeNode(NewAttribute("lunit", "mm"));
981
982 solElement->appendChild(twistedtrapElement);
983}
G4double GetPolarAngleTheta() const
G4double GetAzimuthalAnglePhi() const
G4double GetTiltAngleAlpha() const
G4double GetZHalfLength() const
G4double GetX1HalfLength() const
G4double GetX2HalfLength() const
G4double GetX3HalfLength() const
G4double GetX4HalfLength() const
G4double GetY2HalfLength() const
G4double GetPhiTwist() const
G4double GetY1HalfLength() const

Referenced by AddSolid().

◆ TwistedtrdWrite()

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

Definition at line 986 of file G4GDMLWriteSolids.cc.

988{
989 const G4String& name = GenerateName(twistedtrd->GetName(), twistedtrd);
990
991 xercesc::DOMElement* twistedtrdElement = NewElement("twistedtrd");
992 twistedtrdElement->setAttributeNode(NewAttribute("name", name));
993 twistedtrdElement->setAttributeNode(
994 NewAttribute("x1", 2.0 * twistedtrd->GetX1HalfLength() / mm));
995 twistedtrdElement->setAttributeNode(
996 NewAttribute("x2", 2.0 * twistedtrd->GetX2HalfLength() / mm));
997 twistedtrdElement->setAttributeNode(
998 NewAttribute("y1", 2.0 * twistedtrd->GetY1HalfLength() / mm));
999 twistedtrdElement->setAttributeNode(
1000 NewAttribute("y2", 2.0 * twistedtrd->GetY2HalfLength() / mm));
1001 twistedtrdElement->setAttributeNode(
1002 NewAttribute("z", 2.0 * twistedtrd->GetZHalfLength() / mm));
1003 twistedtrdElement->setAttributeNode(
1004 NewAttribute("PhiTwist", twistedtrd->GetPhiTwist() / degree));
1005 twistedtrdElement->setAttributeNode(NewAttribute("aunit", "deg"));
1006 twistedtrdElement->setAttributeNode(NewAttribute("lunit", "mm"));
1007 solElement->appendChild(twistedtrdElement);
1008}
G4double GetX2HalfLength() const
Definition: G4TwistedTrd.hh:67
G4double GetY2HalfLength() const
Definition: G4TwistedTrd.hh:69
G4double GetY1HalfLength() const
Definition: G4TwistedTrd.hh:68
G4double GetZHalfLength() const
Definition: G4TwistedTrd.hh:70
G4double GetPhiTwist() const
Definition: G4TwistedTrd.hh:71
G4double GetX1HalfLength() const
Definition: G4TwistedTrd.hh:66

Referenced by AddSolid().

◆ TwistedtubsWrite()

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

Definition at line 1011 of file G4GDMLWriteSolids.cc.

1013{
1014 const G4String& name = GenerateName(twistedtubs->GetName(), twistedtubs);
1015
1016 xercesc::DOMElement* twistedtubsElement = NewElement("twistedtubs");
1017 twistedtubsElement->setAttributeNode(NewAttribute("name", name));
1018 twistedtubsElement->setAttributeNode(
1019 NewAttribute("twistedangle", twistedtubs->GetPhiTwist() / degree));
1020 twistedtubsElement->setAttributeNode(
1021 NewAttribute("midinnerrad", twistedtubs->GetInnerRadius() / mm));
1022 twistedtubsElement->setAttributeNode(
1023 NewAttribute("midouterrad", twistedtubs->GetOuterRadius() / mm));
1024 twistedtubsElement->setAttributeNode(
1025 NewAttribute("negativeEndz", twistedtubs->GetEndZ(0) / mm));
1026 twistedtubsElement->setAttributeNode(
1027 NewAttribute("positiveEndz", twistedtubs->GetEndZ(1) / mm));
1028 twistedtubsElement->setAttributeNode(
1029 NewAttribute("phi", twistedtubs->GetDPhi() / degree));
1030 twistedtubsElement->setAttributeNode(NewAttribute("aunit", "deg"));
1031 twistedtubsElement->setAttributeNode(NewAttribute("lunit", "mm"));
1032 solElement->appendChild(twistedtubsElement);
1033}
G4double GetOuterRadius() const
G4double GetPhiTwist() const
G4double GetEndZ(G4int i) const
G4double GetInnerRadius() const
G4double GetDPhi() const

Referenced by AddSolid().

◆ XtruWrite()

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

Definition at line 393 of file G4GDMLWriteSolids.cc.

395{
396 const G4String& name = GenerateName(xtru->GetName(), xtru);
397
398 xercesc::DOMElement* xtruElement = NewElement("xtru");
399 xtruElement->setAttributeNode(NewAttribute("name", name));
400 xtruElement->setAttributeNode(NewAttribute("lunit", "mm"));
401 solElement->appendChild(xtruElement);
402
403 const G4int NumVertex = xtru->GetNofVertices();
404
405 for(G4int i = 0; i < NumVertex; ++i)
406 {
407 xercesc::DOMElement* twoDimVertexElement = NewElement("twoDimVertex");
408 xtruElement->appendChild(twoDimVertexElement);
409
410 const G4TwoVector& vertex = xtru->GetVertex(i);
411
412 twoDimVertexElement->setAttributeNode(NewAttribute("x", vertex.x() / mm));
413 twoDimVertexElement->setAttributeNode(NewAttribute("y", vertex.y() / mm));
414 }
415
416 const G4int NumSection = xtru->GetNofZSections();
417
418 for(G4int i = 0; i < NumSection; ++i)
419 {
420 xercesc::DOMElement* sectionElement = NewElement("section");
421 xtruElement->appendChild(sectionElement);
422
423 const G4ExtrudedSolid::ZSection section = xtru->GetZSection(i);
424
425 sectionElement->setAttributeNode(NewAttribute("zOrder", i));
426 sectionElement->setAttributeNode(
427 NewAttribute("zPosition", section.fZ / mm));
428 sectionElement->setAttributeNode(
429 NewAttribute("xOffset", section.fOffset.x() / mm));
430 sectionElement->setAttributeNode(
431 NewAttribute("yOffset", section.fOffset.y() / mm));
432 sectionElement->setAttributeNode(
433 NewAttribute("scalingFactor", section.fScale));
434 }
435}
double x() const
double y() const
ZSection GetZSection(G4int index) const
G4int GetNofZSections() const
G4int GetNofVertices() const
G4TwoVector GetVertex(G4int index) const

Referenced by AddSolid().

◆ ZplaneWrite()

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

Definition at line 1036 of file G4GDMLWriteSolids.cc.

1039{
1040 xercesc::DOMElement* zplaneElement = NewElement("zplane");
1041 zplaneElement->setAttributeNode(NewAttribute("z", z / mm));
1042 zplaneElement->setAttributeNode(NewAttribute("rmin", rmin / mm));
1043 zplaneElement->setAttributeNode(NewAttribute("rmax", rmax / mm));
1044 element->appendChild(zplaneElement);
1045}

Referenced by G4GDMLWriteParamvol::Polycone_dimensionsWrite(), PolyconeWrite(), G4GDMLWriteParamvol::Polyhedra_dimensionsWrite(), and PolyhedraWrite().

Member Data Documentation

◆ maxTransforms

const G4int G4GDMLWriteSolids::maxTransforms = 8
staticprotected

Definition at line 147 of file G4GDMLWriteSolids.hh.

Referenced by BooleanWrite(), and G4GDMLWriteStructure::TraverseVolumeTree().

◆ solidList

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

Definition at line 145 of file G4GDMLWriteSolids.hh.

Referenced by AddSolid(), and SolidsWrite().

◆ solidsElement

xercesc::DOMElement* G4GDMLWriteSolids::solidsElement = nullptr
protected

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