Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLWriteSolids.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4GDMLWriteSolids implementation
27//
28// Author: Zoltan Torzsok, November 2007
29// --------------------------------------------------------------------
30
31#include "G4GDMLWriteSolids.hh"
32
33#include "G4SystemOfUnits.hh"
34#include "G4BooleanSolid.hh"
35#include "G4ScaledSolid.hh"
36#include "G4Box.hh"
37#include "G4Cons.hh"
38#include "G4Ellipsoid.hh"
39#include "G4EllipticalCone.hh"
40#include "G4EllipticalTube.hh"
41#include "G4ExtrudedSolid.hh"
42#include "G4Hype.hh"
43#include "G4Orb.hh"
44#include "G4Para.hh"
45#include "G4Paraboloid.hh"
47#include "G4Polycone.hh"
48#include "G4GenericPolycone.hh"
49#include "G4Polyhedra.hh"
50#include "G4ReflectedSolid.hh"
51#include "G4Sphere.hh"
52#include "G4SubtractionSolid.hh"
53#include "G4GenericTrap.hh"
54#include "G4TessellatedSolid.hh"
55#include "G4Tet.hh"
56#include "G4Torus.hh"
57#include "G4Trap.hh"
58#include "G4Trd.hh"
59#include "G4Tubs.hh"
60#include "G4CutTubs.hh"
61#include "G4TwistedBox.hh"
62#include "G4TwistedTrap.hh"
63#include "G4TwistedTrd.hh"
64#include "G4TwistedTubs.hh"
65#include "G4UnionSolid.hh"
66#include "G4OpticalSurface.hh"
67#include "G4SurfaceProperty.hh"
69
70// --------------------------------------------------------------------
73{
74}
75
76// --------------------------------------------------------------------
78{
79}
80
81// --------------------------------------------------------------------
82void G4GDMLWriteSolids::MultiUnionWrite(xercesc::DOMElement* solElement,
83 const G4MultiUnion* const munionSolid)
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
105 G4ThreeVector pos = transl.getTranslation();
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}
139
140// --------------------------------------------------------------------
141void G4GDMLWriteSolids::BooleanWrite(xercesc::DOMElement* solElement,
142 const G4BooleanSolid* const boolean)
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}
257
258// --------------------------------------------------------------------
259void G4GDMLWriteSolids::ScaledWrite(xercesc::DOMElement* solElement,
260 const G4ScaledSolid* const scaled)
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}
290
291// --------------------------------------------------------------------
292void G4GDMLWriteSolids::BoxWrite(xercesc::DOMElement* solElement,
293 const G4Box* const box)
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}
308
309// --------------------------------------------------------------------
310void G4GDMLWriteSolids::ConeWrite(xercesc::DOMElement* solElement,
311 const G4Cons* const cone)
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}
335
336// --------------------------------------------------------------------
337void G4GDMLWriteSolids::ElconeWrite(xercesc::DOMElement* solElement,
338 const G4EllipticalCone* const elcone)
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}
354
355// --------------------------------------------------------------------
356void G4GDMLWriteSolids::EllipsoidWrite(xercesc::DOMElement* solElement,
357 const G4Ellipsoid* const ellipsoid)
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}
376
377// --------------------------------------------------------------------
378void G4GDMLWriteSolids::EltubeWrite(xercesc::DOMElement* solElement,
379 const G4EllipticalTube* const eltube)
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}
391
392// --------------------------------------------------------------------
393void G4GDMLWriteSolids::XtruWrite(xercesc::DOMElement* solElement,
394 const G4ExtrudedSolid* const xtru)
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}
436
437// --------------------------------------------------------------------
438void G4GDMLWriteSolids::HypeWrite(xercesc::DOMElement* solElement,
439 const G4Hype* const hype)
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}
459
460// --------------------------------------------------------------------
461void G4GDMLWriteSolids::OrbWrite(xercesc::DOMElement* solElement,
462 const G4Orb* const orb)
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}
472
473// --------------------------------------------------------------------
474void G4GDMLWriteSolids::ParaWrite(xercesc::DOMElement* solElement,
475 const G4Para* const para)
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}
499
500// --------------------------------------------------------------------
501void G4GDMLWriteSolids::ParaboloidWrite(xercesc::DOMElement* solElement,
502 const G4Paraboloid* const paraboloid)
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}
517
518// --------------------------------------------------------------------
519void G4GDMLWriteSolids::PolyconeWrite(xercesc::DOMElement* solElement,
520 const G4Polycone* const polycone)
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}
545
546// --------------------------------------------------------------------
548 xercesc::DOMElement* solElement, const G4GenericPolycone* const polycone)
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}
570
571// --------------------------------------------------------------------
572void G4GDMLWriteSolids::PolyhedraWrite(xercesc::DOMElement* solElement,
573 const G4Polyhedra* const polyhedra)
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}
630
631// --------------------------------------------------------------------
632void G4GDMLWriteSolids::SphereWrite(xercesc::DOMElement* solElement,
633 const G4Sphere* const sphere)
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}
655
656// --------------------------------------------------------------------
658 xercesc::DOMElement* solElement, const G4TessellatedSolid* const tessellated)
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 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}
734
735// --------------------------------------------------------------------
736void G4GDMLWriteSolids::TetWrite(xercesc::DOMElement* solElement,
737 const G4Tet* const tet)
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}
758
759// --------------------------------------------------------------------
760void G4GDMLWriteSolids::TorusWrite(xercesc::DOMElement* solElement,
761 const G4Torus* const torus)
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}
778
779// --------------------------------------------------------------------
780void G4GDMLWriteSolids::GenTrapWrite(xercesc::DOMElement* solElement,
781 const G4GenericTrap* const gtrap)
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}
810
811// --------------------------------------------------------------------
812void G4GDMLWriteSolids::TrapWrite(xercesc::DOMElement* solElement,
813 const G4Trap* const trap)
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}
847
848// --------------------------------------------------------------------
849void G4GDMLWriteSolids::TrdWrite(xercesc::DOMElement* solElement,
850 const G4Trd* const trd)
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}
869
870// --------------------------------------------------------------------
871void G4GDMLWriteSolids::TubeWrite(xercesc::DOMElement* solElement,
872 const G4Tubs* const tube)
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}
892
893// --------------------------------------------------------------------
894void G4GDMLWriteSolids::CutTubeWrite(xercesc::DOMElement* solElement,
895 const G4CutTubs* const cuttube)
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}
927
928// --------------------------------------------------------------------
929void G4GDMLWriteSolids::TwistedboxWrite(xercesc::DOMElement* solElement,
930 const G4TwistedBox* const twistedbox)
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}
948
949// --------------------------------------------------------------------
950void G4GDMLWriteSolids::TwistedtrapWrite(xercesc::DOMElement* solElement,
951 const G4TwistedTrap* const twistedtrap)
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}
984
985// --------------------------------------------------------------------
986void G4GDMLWriteSolids::TwistedtrdWrite(xercesc::DOMElement* solElement,
987 const G4TwistedTrd* const twistedtrd)
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}
1009
1010// --------------------------------------------------------------------
1011void G4GDMLWriteSolids::TwistedtubsWrite(xercesc::DOMElement* solElement,
1012 const G4TwistedTubs* const twistedtubs)
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}
1034
1035// --------------------------------------------------------------------
1036void G4GDMLWriteSolids::ZplaneWrite(xercesc::DOMElement* element,
1037 const G4double& z, const G4double& rmin,
1038 const G4double& rmax)
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}
1046
1047// --------------------------------------------------------------------
1048void G4GDMLWriteSolids::RZPointWrite(xercesc::DOMElement* element,
1049 const G4double& r, const G4double& z)
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}
1056
1057// --------------------------------------------------------------------
1058void G4GDMLWriteSolids::OpticalSurfaceWrite(xercesc::DOMElement* solElement,
1059 const G4OpticalSurface* const surf)
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}
1082
1083// --------------------------------------------------------------------
1084void G4GDMLWriteSolids::PropertyWrite(xercesc::DOMElement* optElement,
1085 const G4OpticalSurface* const surf)
1086{
1087 xercesc::DOMElement* propElement;
1089 auto pvec = ptable->GetProperties();
1090 auto cvec = ptable->GetConstProperties();
1091
1092 for(size_t i = 0; i < pvec.size(); ++i)
1093 {
1094 if(pvec[i] != nullptr) {
1095 propElement = NewElement("property");
1096 propElement->setAttributeNode(
1097 NewAttribute("name", ptable->GetMaterialPropertyNames()[i]));
1098 propElement->setAttributeNode(NewAttribute(
1099 "ref", GenerateName(ptable->GetMaterialPropertyNames()[i],
1100 pvec[i])));
1102 pvec[i]);
1103 optElement->appendChild(propElement);
1104 }
1105 }
1106 for(size_t i = 0; i < cvec.size(); ++i)
1107 {
1108 if (cvec[i].second == true) {
1109 propElement = NewElement("property");
1110 propElement->setAttributeNode(NewAttribute(
1111 "name", ptable->GetMaterialConstPropertyNames()[i]));
1112 propElement->setAttributeNode(NewAttribute(
1113 "ref", ptable->GetMaterialConstPropertyNames()[i]));
1114 xercesc::DOMElement* constElement = NewElement("constant");
1115 constElement->setAttributeNode(NewAttribute(
1116 "name", ptable->GetMaterialConstPropertyNames()[i]));
1117 constElement->setAttributeNode(NewAttribute("value", cvec[i].first));
1118 defineElement->appendChild(constElement);
1119 optElement->appendChild(propElement);
1120 }
1121 }
1122}
1123
1124// --------------------------------------------------------------------
1125void G4GDMLWriteSolids::SolidsWrite(xercesc::DOMElement* gdmlElement)
1126{
1127#ifdef G4VERBOSE
1128 G4cout << "G4GDML: Writing solids..." << G4endl;
1129#endif
1130 solidsElement = NewElement("solids");
1131 gdmlElement->appendChild(solidsElement);
1132
1133 solidList.clear();
1134}
1135
1136// --------------------------------------------------------------------
1137void G4GDMLWriteSolids::AddSolid(const G4VSolid* const solidPtr)
1138{
1139 for(std::size_t i = 0; i < solidList.size(); ++i) // Check if solid is
1140 { // already in the list!
1141 if(solidList[i] == solidPtr)
1142 {
1143 return;
1144 }
1145 }
1146
1147 solidList.push_back(solidPtr);
1148
1149 if(const G4BooleanSolid* const booleanPtr =
1150 dynamic_cast<const G4BooleanSolid*>(solidPtr))
1151 {
1152 BooleanWrite(solidsElement, booleanPtr);
1153 }
1154 else if(const G4ScaledSolid* const scaledPtr =
1155 dynamic_cast<const G4ScaledSolid*>(solidPtr))
1156 {
1157 ScaledWrite(solidsElement, scaledPtr);
1158 }
1159 else if(solidPtr->GetEntityType() == "G4MultiUnion")
1160 {
1161 const G4MultiUnion* const munionPtr =
1162 static_cast<const G4MultiUnion*>(solidPtr);
1163 MultiUnionWrite(solidsElement, munionPtr);
1164 }
1165 else if(solidPtr->GetEntityType() == "G4Box")
1166 {
1167 const G4Box* const boxPtr = static_cast<const G4Box*>(solidPtr);
1168 BoxWrite(solidsElement, boxPtr);
1169 }
1170 else if(solidPtr->GetEntityType() == "G4Cons")
1171 {
1172 const G4Cons* const conePtr = static_cast<const G4Cons*>(solidPtr);
1173 ConeWrite(solidsElement, conePtr);
1174 }
1175 else if(solidPtr->GetEntityType() == "G4EllipticalCone")
1176 {
1177 const G4EllipticalCone* const elconePtr =
1178 static_cast<const G4EllipticalCone*>(solidPtr);
1179 ElconeWrite(solidsElement, elconePtr);
1180 }
1181 else if(solidPtr->GetEntityType() == "G4Ellipsoid")
1182 {
1183 const G4Ellipsoid* const ellipsoidPtr =
1184 static_cast<const G4Ellipsoid*>(solidPtr);
1185 EllipsoidWrite(solidsElement, ellipsoidPtr);
1186 }
1187 else if(solidPtr->GetEntityType() == "G4EllipticalTube")
1188 {
1189 const G4EllipticalTube* const eltubePtr =
1190 static_cast<const G4EllipticalTube*>(solidPtr);
1191 EltubeWrite(solidsElement, eltubePtr);
1192 }
1193 else if(solidPtr->GetEntityType() == "G4ExtrudedSolid")
1194 {
1195 const G4ExtrudedSolid* const xtruPtr =
1196 static_cast<const G4ExtrudedSolid*>(solidPtr);
1197 XtruWrite(solidsElement, xtruPtr);
1198 }
1199 else if(solidPtr->GetEntityType() == "G4Hype")
1200 {
1201 const G4Hype* const hypePtr = static_cast<const G4Hype*>(solidPtr);
1202 HypeWrite(solidsElement, hypePtr);
1203 }
1204 else if(solidPtr->GetEntityType() == "G4Orb")
1205 {
1206 const G4Orb* const orbPtr = static_cast<const G4Orb*>(solidPtr);
1207 OrbWrite(solidsElement, orbPtr);
1208 }
1209 else if(solidPtr->GetEntityType() == "G4Para")
1210 {
1211 const G4Para* const paraPtr = static_cast<const G4Para*>(solidPtr);
1212 ParaWrite(solidsElement, paraPtr);
1213 }
1214 else if(solidPtr->GetEntityType() == "G4Paraboloid")
1215 {
1216 const G4Paraboloid* const paraboloidPtr =
1217 static_cast<const G4Paraboloid*>(solidPtr);
1218 ParaboloidWrite(solidsElement, paraboloidPtr);
1219 }
1220 else if(solidPtr->GetEntityType() == "G4Polycone")
1221 {
1222 const G4Polycone* const polyconePtr =
1223 static_cast<const G4Polycone*>(solidPtr);
1224 PolyconeWrite(solidsElement, polyconePtr);
1225 }
1226 else if(solidPtr->GetEntityType() == "G4GenericPolycone")
1227 {
1228 const G4GenericPolycone* const genpolyconePtr =
1229 static_cast<const G4GenericPolycone*>(solidPtr);
1230 GenericPolyconeWrite(solidsElement, genpolyconePtr);
1231 }
1232 else if(solidPtr->GetEntityType() == "G4Polyhedra")
1233 {
1234 const G4Polyhedra* const polyhedraPtr =
1235 static_cast<const G4Polyhedra*>(solidPtr);
1236 PolyhedraWrite(solidsElement, polyhedraPtr);
1237 }
1238 else if(solidPtr->GetEntityType() == "G4Sphere")
1239 {
1240 const G4Sphere* const spherePtr = static_cast<const G4Sphere*>(solidPtr);
1241 SphereWrite(solidsElement, spherePtr);
1242 }
1243 else if(solidPtr->GetEntityType() == "G4TessellatedSolid")
1244 {
1245 const G4TessellatedSolid* const tessellatedPtr =
1246 static_cast<const G4TessellatedSolid*>(solidPtr);
1247 TessellatedWrite(solidsElement, tessellatedPtr);
1248 }
1249 else if(solidPtr->GetEntityType() == "G4Tet")
1250 {
1251 const G4Tet* const tetPtr = static_cast<const G4Tet*>(solidPtr);
1252 TetWrite(solidsElement, tetPtr);
1253 }
1254 else if(solidPtr->GetEntityType() == "G4Torus")
1255 {
1256 const G4Torus* const torusPtr = static_cast<const G4Torus*>(solidPtr);
1257 TorusWrite(solidsElement, torusPtr);
1258 }
1259 else if(solidPtr->GetEntityType() == "G4GenericTrap")
1260 {
1261 const G4GenericTrap* const gtrapPtr =
1262 static_cast<const G4GenericTrap*>(solidPtr);
1263 GenTrapWrite(solidsElement, gtrapPtr);
1264 }
1265 else if(solidPtr->GetEntityType() == "G4Trap")
1266 {
1267 const G4Trap* const trapPtr = static_cast<const G4Trap*>(solidPtr);
1268 TrapWrite(solidsElement, trapPtr);
1269 }
1270 else if(solidPtr->GetEntityType() == "G4Trd")
1271 {
1272 const G4Trd* const trdPtr = static_cast<const G4Trd*>(solidPtr);
1273 TrdWrite(solidsElement, trdPtr);
1274 }
1275 else if(solidPtr->GetEntityType() == "G4Tubs")
1276 {
1277 const G4Tubs* const tubePtr = static_cast<const G4Tubs*>(solidPtr);
1278 TubeWrite(solidsElement, tubePtr);
1279 }
1280 else if(solidPtr->GetEntityType() == "G4CutTubs")
1281 {
1282 const G4CutTubs* const cuttubePtr = static_cast<const G4CutTubs*>(solidPtr);
1283 CutTubeWrite(solidsElement, cuttubePtr);
1284 }
1285 else if(solidPtr->GetEntityType() == "G4TwistedBox")
1286 {
1287 const G4TwistedBox* const twistedboxPtr =
1288 static_cast<const G4TwistedBox*>(solidPtr);
1289 TwistedboxWrite(solidsElement, twistedboxPtr);
1290 }
1291 else if(solidPtr->GetEntityType() == "G4TwistedTrap")
1292 {
1293 const G4TwistedTrap* const twistedtrapPtr =
1294 static_cast<const G4TwistedTrap*>(solidPtr);
1295 TwistedtrapWrite(solidsElement, twistedtrapPtr);
1296 }
1297 else if(solidPtr->GetEntityType() == "G4TwistedTrd")
1298 {
1299 const G4TwistedTrd* const twistedtrdPtr =
1300 static_cast<const G4TwistedTrd*>(solidPtr);
1301 TwistedtrdWrite(solidsElement, twistedtrdPtr);
1302 }
1303 else if(solidPtr->GetEntityType() == "G4TwistedTubs")
1304 {
1305 const G4TwistedTubs* const twistedtubsPtr =
1306 static_cast<const G4TwistedTubs*>(solidPtr);
1307 TwistedtubsWrite(solidsElement, twistedtubsPtr);
1308 }
1309 else
1310 {
1311 G4String error_msg = "Unknown solid: " + solidPtr->GetName() +
1312 "; Type: " + solidPtr->GetEntityType();
1313 G4Exception("G4GDMLWriteSolids::AddSolid()", "WriteError", FatalException,
1314 error_msg);
1315 }
1316}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4OpticalSurfaceModel
@ glisur
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4double alpha2
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
double z() const
double phi() const
double getZ() const
double theta() const
double x() const
double y() const
double getX() const
double getY() const
Definition: G4Box.hh:56
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
Definition: G4Cons.hh:78
G4double GetOuterRadiusPlusZ() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4double GetInnerRadiusMinusZ() const
G4double GetInnerRadiusPlusZ() const
G4double GetOuterRadiusMinusZ() const
G4double GetZHalfLength() const
G4ThreeVector GetHighNorm() const
G4double GetStartPhiAngle() const
G4double GetZHalfLength() const
G4double GetInnerRadius() const
G4ThreeVector GetLowNorm() const
G4double GetDeltaPhiAngle() const
G4double GetOuterRadius() const
G4double GetSemiAxisMax(G4int i) const
G4double GetZTopCut() const
G4double GetZBottomCut() const
G4double GetSemiAxisX() const
G4double GetSemiAxisY() const
G4double GetZMax() const
G4double GetZTopCut() const
G4double GetDy() const
G4double GetDx() const
G4double GetDz() const
ZSection GetZSection(G4int index) const
G4int GetNofZSections() const
G4int GetNofVertices() const
G4TwoVector GetVertex(G4int index) const
xercesc::DOMElement * defineElement
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 AddPosition(const G4String &name, const G4ThreeVector &pos)
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 &)
void ScaleWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
void PropertyVectorWrite(const G4String &, const G4PhysicsFreeVector *const)
void TwistedtrdWrite(xercesc::DOMElement *, const G4TwistedTrd *const)
void TorusWrite(xercesc::DOMElement *, const G4Torus *const)
void RZPointWrite(xercesc::DOMElement *, const G4double &, const G4double &)
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)
virtual void AddSolid(const G4VSolid *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)
virtual void SolidsWrite(xercesc::DOMElement *)
xercesc::DOMElement * solidsElement
virtual ~G4GDMLWriteSolids()
static const G4int maxTransforms
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 PropertyWrite(xercesc::DOMElement *, const G4OpticalSurface *const)
void TrdWrite(xercesc::DOMElement *, const G4Trd *const)
void XtruWrite(xercesc::DOMElement *, const G4ExtrudedSolid *const)
void OpticalSurfaceWrite(xercesc::DOMElement *, const G4OpticalSurface *const)
void ZplaneWrite(xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:192
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:132
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:155
G4double GetStartPhi() const
G4double GetEndPhi() const
G4int GetNumRZCorner() const
G4PolyconeSideRZ GetCorner(G4int index) const
G4double GetZHalfLength() const
const std::vector< G4TwoVector > & GetVertices() const
Definition: G4Hype.hh:69
G4double GetInnerStereo() const
G4double GetZHalfLength() const
G4double GetOuterStereo() const
G4double GetOuterRadius() const
G4double GetInnerRadius() 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
const G4Transform3D & GetTransformation(G4int index) const
G4int GetNumberOfSolids() const
G4VSolid * GetSolid(G4int index) const
G4OpticalSurfaceModel GetModel() const
G4double GetSigmaAlpha() const
G4OpticalSurfaceFinish GetFinish() const
G4double GetPolish() const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Orb.hh:56
G4double GetRadius() const
Definition: G4Para.hh:79
G4double GetTanAlpha() const
G4ThreeVector GetSymAxis() const
G4double GetYHalfLength() const
G4double GetZHalfLength() const
G4double GetXHalfLength() const
G4double GetRadiusPlusZ() const
G4double GetRadiusMinusZ() const
G4double GetZHalfLength() const
G4PolyconeHistorical * GetOriginalParameters() const
G4int GetNumRZCorner() const
G4PolyhedraHistorical * GetOriginalParameters() const
G4PolyhedraSideRZ GetCorner(const G4int index) const
G4bool IsGeneric() const
G4VSolid * GetUnscaledSolid() const
G4Scale3D GetScaleTransform() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetDeltaThetaAngle() const
G4double GetStartThetaAngle() const
const G4String & GetName() const
const G4SurfaceType & GetType() const
G4int GetNumberOfFacets() const
G4VFacet * GetFacet(G4int i) const
Definition: G4Tet.hh:56
void GetVertices(G4ThreeVector &anchor, G4ThreeVector &p1, G4ThreeVector &p2, G4ThreeVector &p3) const
Definition: G4Tet.cc:286
G4double GetDPhi() const
G4double GetRmin() const
G4double GetRtor() const
G4double GetRmax() const
G4double GetSPhi() const
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
Definition: G4Trd.hh:63
G4double GetXHalfLength2() const
G4double GetYHalfLength2() const
G4double GetXHalfLength1() const
G4double GetYHalfLength1() const
G4double GetZHalfLength() const
Definition: G4Tubs.hh:75
G4double GetZHalfLength() const
G4double GetInnerRadius() const
G4double GetOuterRadius() const
G4double GetStartPhiAngle() const
G4double GetDeltaPhiAngle() const
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
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
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
G4double GetOuterRadius() const
G4double GetPhiTwist() const
G4double GetEndZ(G4int i) const
G4double GetInnerRadius() const
G4double GetDPhi() const
virtual G4ThreeVector GetVertex(G4int i) const =0
virtual G4int GetNumberOfVertices() const =0
G4String GetName() const
virtual G4GeometryType GetEntityType() const =0
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
Definition: xmlparse.c:284