Geant4 11.2.2
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 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 G4PhysicsFreeVector *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 G4PhysicsFreeVector * > 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
 
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 1144 of file G4GDMLWriteSolids.cc.

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

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 901 of file G4GDMLWriteSolids.cc.

903{
904 const G4String& name = GenerateName(cuttube->GetName(), cuttube);
905
906 xercesc::DOMElement* cuttubeElement = NewElement("cutTube");
907 cuttubeElement->setAttributeNode(NewAttribute("name", name));
908 cuttubeElement->setAttributeNode(
909 NewAttribute("rmin", cuttube->GetInnerRadius() / mm));
910 cuttubeElement->setAttributeNode(
911 NewAttribute("rmax", cuttube->GetOuterRadius() / mm));
912 cuttubeElement->setAttributeNode(
913 NewAttribute("z", 2.0 * cuttube->GetZHalfLength() / mm));
914 cuttubeElement->setAttributeNode(
915 NewAttribute("startphi", cuttube->GetStartPhiAngle() / degree));
916 cuttubeElement->setAttributeNode(
917 NewAttribute("deltaphi", cuttube->GetDeltaPhiAngle() / degree));
918 cuttubeElement->setAttributeNode(
919 NewAttribute("lowX", cuttube->GetLowNorm().getX() / mm));
920 cuttubeElement->setAttributeNode(
921 NewAttribute("lowY", cuttube->GetLowNorm().getY() / mm));
922 cuttubeElement->setAttributeNode(
923 NewAttribute("lowZ", cuttube->GetLowNorm().getZ() / mm));
924 cuttubeElement->setAttributeNode(
925 NewAttribute("highX", cuttube->GetHighNorm().getX() / mm));
926 cuttubeElement->setAttributeNode(
927 NewAttribute("highY", cuttube->GetHighNorm().getY() / mm));
928 cuttubeElement->setAttributeNode(
929 NewAttribute("highZ", cuttube->GetHighNorm().getZ() / mm));
930 cuttubeElement->setAttributeNode(NewAttribute("aunit", "deg"));
931 cuttubeElement->setAttributeNode(NewAttribute("lunit", "mm"));
932 solElement->appendChild(cuttubeElement);
933}
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 G4int num_rzpoints = (G4int)polycone->GetNumRZCorner();
563 for(G4int 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 787 of file G4GDMLWriteSolids.cc.

789{
790 const G4String& name = GenerateName(gtrap->GetName(), gtrap);
791
792 std::vector<G4TwoVector> vertices = gtrap->GetVertices();
793
794 xercesc::DOMElement* gtrapElement = NewElement("arb8");
795 gtrapElement->setAttributeNode(NewAttribute("name", name));
796 gtrapElement->setAttributeNode(
797 NewAttribute("dz", gtrap->GetZHalfLength() / mm));
798 gtrapElement->setAttributeNode(NewAttribute("v1x", vertices[0].x()));
799 gtrapElement->setAttributeNode(NewAttribute("v1y", vertices[0].y()));
800 gtrapElement->setAttributeNode(NewAttribute("v2x", vertices[1].x()));
801 gtrapElement->setAttributeNode(NewAttribute("v2y", vertices[1].y()));
802 gtrapElement->setAttributeNode(NewAttribute("v3x", vertices[2].x()));
803 gtrapElement->setAttributeNode(NewAttribute("v3y", vertices[2].y()));
804 gtrapElement->setAttributeNode(NewAttribute("v4x", vertices[3].x()));
805 gtrapElement->setAttributeNode(NewAttribute("v4y", vertices[3].y()));
806 gtrapElement->setAttributeNode(NewAttribute("v5x", vertices[4].x()));
807 gtrapElement->setAttributeNode(NewAttribute("v5y", vertices[4].y()));
808 gtrapElement->setAttributeNode(NewAttribute("v6x", vertices[5].x()));
809 gtrapElement->setAttributeNode(NewAttribute("v6y", vertices[5].y()));
810 gtrapElement->setAttributeNode(NewAttribute("v7x", vertices[6].x()));
811 gtrapElement->setAttributeNode(NewAttribute("v7y", vertices[6].y()));
812 gtrapElement->setAttributeNode(NewAttribute("v8x", vertices[7].x()));
813 gtrapElement->setAttributeNode(NewAttribute("v8y", vertices[7].y()));
814 gtrapElement->setAttributeNode(NewAttribute("lunit", "mm"));
815 solElement->appendChild(gtrapElement);
816}
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", 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+"_"+nodeName+"_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+"_"+nodeName+"_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
double yz() const
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
double xy() const
CLHEP::Hep3Vector getTranslation() const
double zx() const
double yx() const
double zy() const
double xx() const
double yy() const
double xz() const

Referenced by AddSolid().

◆ OpticalSurfaceWrite()

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

Definition at line 1065 of file G4GDMLWriteSolids.cc.

1067{
1068 xercesc::DOMElement* optElement = NewElement("opticalsurface");
1069 G4OpticalSurfaceModel smodel = surf->GetModel();
1070 G4double sval =
1071 (smodel == glisur) ? surf->GetPolish() : surf->GetSigmaAlpha();
1072 const G4String& name = GenerateName(surf->GetName(), surf);
1073
1074 optElement->setAttributeNode(NewAttribute("name", name));
1075 optElement->setAttributeNode(NewAttribute("model", smodel));
1076 optElement->setAttributeNode(NewAttribute("finish", surf->GetFinish()));
1077 optElement->setAttributeNode(NewAttribute("type", surf->GetType()));
1078 optElement->setAttributeNode(NewAttribute("value", sval));
1079
1080 // Write any property attached to the optical surface...
1081 //
1082 if(surf->GetMaterialPropertiesTable())
1083 {
1084 PropertyWrite(optElement, surf);
1085 }
1086
1087 solElement->appendChild(optElement);
1088}
G4OpticalSurfaceModel
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 G4int num_rzpoints = (G4int)polyhedra->GetNumRZCorner();
621
622 for(G4int 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 1091 of file G4GDMLWriteSolids.cc.

1093{
1094 xercesc::DOMElement* propElement;
1096 auto pvec = ptable->GetProperties();
1097 auto cvec = ptable->GetConstProperties();
1098
1099 for(size_t i = 0; i < pvec.size(); ++i)
1100 {
1101 if(pvec[i] != nullptr) {
1102 propElement = NewElement("property");
1103 propElement->setAttributeNode(
1104 NewAttribute("name", ptable->GetMaterialPropertyNames()[i]));
1105 propElement->setAttributeNode(NewAttribute(
1106 "ref", GenerateName(ptable->GetMaterialPropertyNames()[i],
1107 pvec[i])));
1109 pvec[i]);
1110 optElement->appendChild(propElement);
1111 }
1112 }
1113 for(size_t i = 0; i < cvec.size(); ++i)
1114 {
1115 if (cvec[i].second == true) {
1116 propElement = NewElement("property");
1117 propElement->setAttributeNode(NewAttribute(
1118 "name", ptable->GetMaterialConstPropertyNames()[i]));
1119 propElement->setAttributeNode(NewAttribute(
1120 "ref", ptable->GetMaterialConstPropertyNames()[i]));
1121 xercesc::DOMElement* constElement = NewElement("constant");
1122 constElement->setAttributeNode(NewAttribute(
1123 "name", ptable->GetMaterialConstPropertyNames()[i]));
1124 constElement->setAttributeNode(NewAttribute("value", cvec[i].first));
1125 defineElement->appendChild(constElement);
1126 optElement->appendChild(propElement);
1127 }
1128 }
1129}
xercesc::DOMElement * defineElement
void PropertyVectorWrite(const G4String &, const G4PhysicsFreeVector *const)
const std::vector< G4String > & GetMaterialConstPropertyNames() const
const std::vector< std::pair< G4double, G4bool > > & GetConstProperties() const
const std::vector< G4MaterialPropertyVector * > & GetProperties() const
const std::vector< G4String > & GetMaterialPropertyNames() const

Referenced by OpticalSurfaceWrite().

◆ RZPointWrite()

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

Definition at line 1055 of file G4GDMLWriteSolids.cc.

1057{
1058 xercesc::DOMElement* rzpointElement = NewElement("rzpoint");
1059 rzpointElement->setAttributeNode(NewAttribute("r", r / mm));
1060 rzpointElement->setAttributeNode(NewAttribute("z", z / mm));
1061 element->appendChild(rzpointElement);
1062}

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 1132 of file G4GDMLWriteSolids.cc.

1133{
1134#ifdef G4VERBOSE
1135 G4cout << "G4GDML: Writing solids..." << G4endl;
1136#endif
1137 solidsElement = NewElement("solids");
1138 gdmlElement->appendChild(solidsElement);
1139
1140 solidList.clear();
1141}
#define G4endl
Definition G4ios.hh:67
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((G4int)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((G4int)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 if ( ! vertexMap.insert(std::make_pair(vertex, ref)).second )
724 {
725 G4ExceptionDescription description;
726 description << "Failed to insert [vertex, ref] " << vertex << ", "
727 << ref << " in map.";
728 G4Exception("G4GDMLWriteSolids::TessellatedWrite", "WriteError",
729 JustWarning, description);
730 }
731 AddPosition(ref, vertex); // ... add it to define section!
732 ++NumVertex;
733 }
734
735 // Now create association of the vertex with its facet
736 //
737 facetElement->setAttributeNode(NewAttribute(fname, ref));
738 }
739 }
740}
@ JustWarning
std::ostringstream G4ExceptionDescription
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 743 of file G4GDMLWriteSolids.cc.

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

Referenced by AddSolid().

◆ TorusWrite()

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

Definition at line 767 of file G4GDMLWriteSolids.cc.

769{
770 const G4String& name = GenerateName(torus->GetName(), torus);
771
772 xercesc::DOMElement* torusElement = NewElement("torus");
773 torusElement->setAttributeNode(NewAttribute("name", name));
774 torusElement->setAttributeNode(NewAttribute("rmin", torus->GetRmin() / mm));
775 torusElement->setAttributeNode(NewAttribute("rmax", torus->GetRmax() / mm));
776 torusElement->setAttributeNode(NewAttribute("rtor", torus->GetRtor() / mm));
777 torusElement->setAttributeNode(
778 NewAttribute("startphi", torus->GetSPhi() / degree));
779 torusElement->setAttributeNode(
780 NewAttribute("deltaphi", torus->GetDPhi() / degree));
781 torusElement->setAttributeNode(NewAttribute("aunit", "deg"));
782 torusElement->setAttributeNode(NewAttribute("lunit", "mm"));
783 solElement->appendChild(torusElement);
784}
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 819 of file G4GDMLWriteSolids.cc.

821{
822 const G4String& name = GenerateName(trap->GetName(), trap);
823
824 const G4ThreeVector& simaxis = trap->GetSymAxis();
825 const G4double phi = simaxis.phi();
826 const G4double theta = simaxis.theta();
827 const G4double alpha1 = std::atan(trap->GetTanAlpha1());
828 const G4double alpha2 = std::atan(trap->GetTanAlpha2());
829
830 xercesc::DOMElement* trapElement = NewElement("trap");
831 trapElement->setAttributeNode(NewAttribute("name", name));
832 trapElement->setAttributeNode(
833 NewAttribute("z", 2.0 * trap->GetZHalfLength() / mm));
834 trapElement->setAttributeNode(NewAttribute("theta", theta / degree));
835 trapElement->setAttributeNode(NewAttribute("phi", phi / degree));
836 trapElement->setAttributeNode(
837 NewAttribute("y1", 2.0 * trap->GetYHalfLength1() / mm));
838 trapElement->setAttributeNode(
839 NewAttribute("x1", 2.0 * trap->GetXHalfLength1() / mm));
840 trapElement->setAttributeNode(
841 NewAttribute("x2", 2.0 * trap->GetXHalfLength2() / mm));
842 trapElement->setAttributeNode(NewAttribute("alpha1", alpha1 / degree));
843 trapElement->setAttributeNode(
844 NewAttribute("y2", 2.0 * trap->GetYHalfLength2() / mm));
845 trapElement->setAttributeNode(
846 NewAttribute("x3", 2.0 * trap->GetXHalfLength3() / mm));
847 trapElement->setAttributeNode(
848 NewAttribute("x4", 2.0 * trap->GetXHalfLength4() / mm));
849 trapElement->setAttributeNode(NewAttribute("alpha2", alpha2 / degree));
850 trapElement->setAttributeNode(NewAttribute("aunit", "deg"));
851 trapElement->setAttributeNode(NewAttribute("lunit", "mm"));
852 solElement->appendChild(trapElement);
853}
const G4double alpha2
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 856 of file G4GDMLWriteSolids.cc.

858{
859 const G4String& name = GenerateName(trd->GetName(), trd);
860
861 xercesc::DOMElement* trdElement = NewElement("trd");
862 trdElement->setAttributeNode(NewAttribute("name", name));
863 trdElement->setAttributeNode(
864 NewAttribute("x1", 2.0 * trd->GetXHalfLength1() / mm));
865 trdElement->setAttributeNode(
866 NewAttribute("x2", 2.0 * trd->GetXHalfLength2() / mm));
867 trdElement->setAttributeNode(
868 NewAttribute("y1", 2.0 * trd->GetYHalfLength1() / mm));
869 trdElement->setAttributeNode(
870 NewAttribute("y2", 2.0 * trd->GetYHalfLength2() / mm));
871 trdElement->setAttributeNode(
872 NewAttribute("z", 2.0 * trd->GetZHalfLength() / mm));
873 trdElement->setAttributeNode(NewAttribute("lunit", "mm"));
874 solElement->appendChild(trdElement);
875}
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 878 of file G4GDMLWriteSolids.cc.

880{
881 const G4String& name = GenerateName(tube->GetName(), tube);
882
883 xercesc::DOMElement* tubeElement = NewElement("tube");
884 tubeElement->setAttributeNode(NewAttribute("name", name));
885 tubeElement->setAttributeNode(
886 NewAttribute("rmin", tube->GetInnerRadius() / mm));
887 tubeElement->setAttributeNode(
888 NewAttribute("rmax", tube->GetOuterRadius() / mm));
889 tubeElement->setAttributeNode(
890 NewAttribute("z", 2.0 * tube->GetZHalfLength() / mm));
891 tubeElement->setAttributeNode(
892 NewAttribute("startphi", tube->GetStartPhiAngle() / degree));
893 tubeElement->setAttributeNode(
894 NewAttribute("deltaphi", tube->GetDeltaPhiAngle() / degree));
895 tubeElement->setAttributeNode(NewAttribute("aunit", "deg"));
896 tubeElement->setAttributeNode(NewAttribute("lunit", "mm"));
897 solElement->appendChild(tubeElement);
898}
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 936 of file G4GDMLWriteSolids.cc.

938{
939 const G4String& name = GenerateName(twistedbox->GetName(), twistedbox);
940
941 xercesc::DOMElement* twistedboxElement = NewElement("twistedbox");
942 twistedboxElement->setAttributeNode(NewAttribute("name", name));
943 twistedboxElement->setAttributeNode(
944 NewAttribute("x", 2.0 * twistedbox->GetXHalfLength() / mm));
945 twistedboxElement->setAttributeNode(
946 NewAttribute("y", 2.0 * twistedbox->GetYHalfLength() / mm));
947 twistedboxElement->setAttributeNode(
948 NewAttribute("z", 2.0 * twistedbox->GetZHalfLength() / mm));
949 twistedboxElement->setAttributeNode(
950 NewAttribute("PhiTwist", twistedbox->GetPhiTwist() / degree));
951 twistedboxElement->setAttributeNode(NewAttribute("aunit", "deg"));
952 twistedboxElement->setAttributeNode(NewAttribute("lunit", "mm"));
953 solElement->appendChild(twistedboxElement);
954}
G4double GetPhiTwist() const
G4double GetXHalfLength() const
G4double GetZHalfLength() const
G4double GetYHalfLength() const

Referenced by AddSolid().

◆ TwistedtrapWrite()

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

Definition at line 957 of file G4GDMLWriteSolids.cc.

959{
960 const G4String& name = GenerateName(twistedtrap->GetName(), twistedtrap);
961
962 xercesc::DOMElement* twistedtrapElement = NewElement("twistedtrap");
963 twistedtrapElement->setAttributeNode(NewAttribute("name", name));
964 twistedtrapElement->setAttributeNode(
965 NewAttribute("y1", 2.0 * twistedtrap->GetY1HalfLength() / mm));
966 twistedtrapElement->setAttributeNode(
967 NewAttribute("x1", 2.0 * twistedtrap->GetX1HalfLength() / mm));
968 twistedtrapElement->setAttributeNode(
969 NewAttribute("x2", 2.0 * twistedtrap->GetX2HalfLength() / mm));
970 twistedtrapElement->setAttributeNode(
971 NewAttribute("y2", 2.0 * twistedtrap->GetY2HalfLength() / mm));
972 twistedtrapElement->setAttributeNode(
973 NewAttribute("x3", 2.0 * twistedtrap->GetX3HalfLength() / mm));
974 twistedtrapElement->setAttributeNode(
975 NewAttribute("x4", 2.0 * twistedtrap->GetX4HalfLength() / mm));
976 twistedtrapElement->setAttributeNode(
977 NewAttribute("z", 2.0 * twistedtrap->GetZHalfLength() / mm));
978 twistedtrapElement->setAttributeNode(
979 NewAttribute("Alph", twistedtrap->GetTiltAngleAlpha() / degree));
980 twistedtrapElement->setAttributeNode(
981 NewAttribute("Theta", twistedtrap->GetPolarAngleTheta() / degree));
982 twistedtrapElement->setAttributeNode(
983 NewAttribute("Phi", twistedtrap->GetAzimuthalAnglePhi() / degree));
984 twistedtrapElement->setAttributeNode(
985 NewAttribute("PhiTwist", twistedtrap->GetPhiTwist() / degree));
986 twistedtrapElement->setAttributeNode(NewAttribute("aunit", "deg"));
987 twistedtrapElement->setAttributeNode(NewAttribute("lunit", "mm"));
988
989 solElement->appendChild(twistedtrapElement);
990}
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 993 of file G4GDMLWriteSolids.cc.

995{
996 const G4String& name = GenerateName(twistedtrd->GetName(), twistedtrd);
997
998 xercesc::DOMElement* twistedtrdElement = NewElement("twistedtrd");
999 twistedtrdElement->setAttributeNode(NewAttribute("name", name));
1000 twistedtrdElement->setAttributeNode(
1001 NewAttribute("x1", 2.0 * twistedtrd->GetX1HalfLength() / mm));
1002 twistedtrdElement->setAttributeNode(
1003 NewAttribute("x2", 2.0 * twistedtrd->GetX2HalfLength() / mm));
1004 twistedtrdElement->setAttributeNode(
1005 NewAttribute("y1", 2.0 * twistedtrd->GetY1HalfLength() / mm));
1006 twistedtrdElement->setAttributeNode(
1007 NewAttribute("y2", 2.0 * twistedtrd->GetY2HalfLength() / mm));
1008 twistedtrdElement->setAttributeNode(
1009 NewAttribute("z", 2.0 * twistedtrd->GetZHalfLength() / mm));
1010 twistedtrdElement->setAttributeNode(
1011 NewAttribute("PhiTwist", twistedtrd->GetPhiTwist() / degree));
1012 twistedtrdElement->setAttributeNode(NewAttribute("aunit", "deg"));
1013 twistedtrdElement->setAttributeNode(NewAttribute("lunit", "mm"));
1014 solElement->appendChild(twistedtrdElement);
1015}
G4double GetX2HalfLength() const
G4double GetY2HalfLength() const
G4double GetY1HalfLength() const
G4double GetZHalfLength() const
G4double GetPhiTwist() const
G4double GetX1HalfLength() const

Referenced by AddSolid().

◆ TwistedtubsWrite()

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

Definition at line 1018 of file G4GDMLWriteSolids.cc.

1020{
1021 const G4String& name = GenerateName(twistedtubs->GetName(), twistedtubs);
1022
1023 xercesc::DOMElement* twistedtubsElement = NewElement("twistedtubs");
1024 twistedtubsElement->setAttributeNode(NewAttribute("name", name));
1025 twistedtubsElement->setAttributeNode(
1026 NewAttribute("twistedangle", twistedtubs->GetPhiTwist() / degree));
1027 twistedtubsElement->setAttributeNode(
1028 NewAttribute("midinnerrad", twistedtubs->GetInnerRadius() / mm));
1029 twistedtubsElement->setAttributeNode(
1030 NewAttribute("midouterrad", twistedtubs->GetOuterRadius() / mm));
1031 twistedtubsElement->setAttributeNode(
1032 NewAttribute("negativeEndz", twistedtubs->GetEndZ(0) / mm));
1033 twistedtubsElement->setAttributeNode(
1034 NewAttribute("positiveEndz", twistedtubs->GetEndZ(1) / mm));
1035 twistedtubsElement->setAttributeNode(
1036 NewAttribute("phi", twistedtubs->GetDPhi() / degree));
1037 twistedtubsElement->setAttributeNode(NewAttribute("aunit", "deg"));
1038 twistedtubsElement->setAttributeNode(NewAttribute("lunit", "mm"));
1039 solElement->appendChild(twistedtubsElement);
1040}
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 1043 of file G4GDMLWriteSolids.cc.

1046{
1047 xercesc::DOMElement* zplaneElement = NewElement("zplane");
1048 zplaneElement->setAttributeNode(NewAttribute("z", z / mm));
1049 zplaneElement->setAttributeNode(NewAttribute("rmin", rmin / mm));
1050 zplaneElement->setAttributeNode(NewAttribute("rmax", rmax / mm));
1051 element->appendChild(zplaneElement);
1052}

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

Member Data Documentation

◆ maxTransforms

const G4int G4GDMLWriteSolids::maxTransforms = 8
staticprotected

Definition at line 146 of file G4GDMLWriteSolids.hh.

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

◆ solidList

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

Definition at line 144 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: