Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
HepPolyhedron.h
이 파일의 문서화 페이지로 가기
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 //
27 // $Id: HepPolyhedron.h 106123 2017-09-13 12:53:03Z gcosmo $
28 //
29 //
30 // Class Description:
31 // HepPolyhedron is an intermediate class between description of a shape
32 // and visualization systems. It is intended to provide some service like:
33 // - polygonization of shapes with triangulization (quadrilaterization)
34 // of complex polygons;
35 // - calculation of normals for faces and vertices;
36 // - finding result of boolean operation on polyhedra;
37 //
38 // Public constructors:
39 //
40 // HepPolyhedronBox (dx,dy,dz)
41 // - create polyhedron for Box;
42 // HepPolyhedronTrd1 (dx1,dx2,dy,dz)
43 // - create polyhedron for Trd1;
44 // HepPolyhedronTrd2 (dx1,dx2,dy1,dy2,dz)
45 // - create polyhedron for Trd2;
46 // HepPolyhedronTrap (dz,theta,phi, h1,bl1,tl1,alp1, h2,bl2,tl2,alp2)
47 // - create polyhedron for Trap;
48 // HepPolyhedronPara (dx,dy,dz,alpha,theta,phi)
49 // - create polyhedron for Para;
50 // HepPolyhedronTube (rmin,rmax,dz)
51 // - create polyhedron for Tube;
52 // HepPolyhedronTubs (rmin,rmax,dz,phi1,dphi)
53 // - create polyhedron for Tubs;
54 // HepPolyhedronCone (rmin1,rmax1,rmin2,rmax2,dz)
55 // - create polyhedron for Cone;
56 // HepPolyhedronCons (rmin1,rmax1,rmin2,rmax2,dz,phi1,dphi)
57 // - create polyhedron for Cons;
58 // HepPolyhedronPgon (phi,dphi,npdv,nz, z(*),rmin(*),rmax(*))
59 // - create polyhedron for Pgon;
60 // HepPolyhedronPcon (phi,dphi,nz, z(*),rmin(*),rmax(*))
61 // - create polyhedron for Pcon;
62 // HepPolyhedronSphere (rmin,rmax,phi,dphi,the,dthe)
63 // - create polyhedron for Sphere;
64 // HepPolyhedronTorus (rmin,rmax,rtor,phi,dphi)
65 // - create polyhedron for Torus;
66 // HepPolyhedronEllipsoid (dx,dy,dz,zcut1,zcut2)
67 // - create polyhedron for Ellipsoid;
68 // Public functions:
69 //
70 // GetNoVertices () - returns number of vertices;
71 // GetNoFacets () - returns number of faces;
72 // GetNextVertexIndex (index,edgeFlag) - get vertex indices of the
73 // quadrilaterals in order;
74 // returns false when finished each face;
75 // GetVertex (index) - returns vertex by index;
76 // GetNextVertex (vertex,edgeFlag) - get vertices with edge visibility
77 // of the quadrilaterals in order;
78 // returns false when finished each face;
79 // GetNextVertex (vertex,edgeFlag,normal) - get vertices with edge
80 // visibility and normal of the quadrilaterals
81 // in order; returns false when finished each face;
82 // GetNextEdgeIndices (i1,i2,edgeFlag) - get indices of the next edge;
83 // returns false for the last edge;
84 // GetNextEdgeIndices (i1,i2,edgeFlag,iface1,iface2) - get indices of
85 // the next edge with indices of the faces
86 // to which the edge belongs;
87 // returns false for the last edge;
88 // GetNextEdge (p1,p2,edgeFlag) - get next edge;
89 // returns false for the last edge;
90 // GetNextEdge (p1,p2,edgeFlag,iface1,iface2) - get next edge with indices
91 // of the faces to which the edge belongs;
92 // returns false for the last edge;
93 // GetFacet (index,n,nodes,edgeFlags=0,normals=0) - get face by index;
94 // GetNextFacet (n,nodes,edgeFlags=0,normals=0) - get next face with normals
95 // at the nodes; returns false for the last face;
96 // GetNormal (index) - get normal of face given by index;
97 // GetUnitNormal (index) - get unit normal of face given by index;
98 // GetNextNormal (normal) - get normals of each face in order;
99 // returns false when finished all faces;
100 // GetNextUnitNormal (normal) - get normals of unit length of each face
101 // in order; returns false when finished all faces;
102 // GetSurfaceArea() - get surface area of the polyhedron;
103 // GetVolume() - get volume of the polyhedron;
104 // GetNumberOfRotationSteps() - get number of steps for whole circle;
105 // SetNumberOfRotationSteps (n) - set number of steps for whole circle;
106 // ResetNumberOfRotationSteps() - reset number of steps for whole circle
107 // to default value;
108 // History:
109 //
110 // 20.06.96 Evgeni Chernyaev <Evgueni.Tcherniaev@cern.ch> - initial version
111 //
112 // 23.07.96 John Allison
113 // - added GetNoVertices, GetNoFacets, GetNextVertex, GetNextNormal
114 //
115 // 30.09.96 E.Chernyaev
116 // - added GetNextVertexIndex, GetVertex by Yasuhide Sawada
117 // - added GetNextUnitNormal, GetNextEdgeIndices, GetNextEdge
118 // - improvements: angles now expected in radians
119 // int -> G4int, double -> G4double
120 // - G4ThreeVector replaced by either G4Point3D or G4Normal3D
121 //
122 // 15.12.96 E.Chernyaev
123 // - private functions G4PolyhedronAlloc, G4PolyhedronPrism renamed
124 // to AllocateMemory and CreatePrism
125 // - added private functions GetNumberOfRotationSteps, RotateEdge,
126 // RotateAroundZ, SetReferences
127 // - rewritten G4PolyhedronCons;
128 // - added G4PolyhedronPara, ...Trap, ...Pgon, ...Pcon, ...Sphere, ...Torus,
129 // so full List of implemented shapes now looks like:
130 // BOX, TRD1, TRD2, TRAP, TUBE, TUBS, CONE, CONS, PARA, PGON, PCON,
131 // SPHERE, TORUS
132 //
133 // 01.06.97 E.Chernyaev
134 // - RotateAroundZ modified and SetSideFacets added to allow Rmin=Rmax
135 // in bodies of revolution
136 //
137 // 24.06.97 J.Allison
138 // - added static private member fNumberOfRotationSteps and static public
139 // functions void SetNumberOfRotationSteps (G4int n) and
140 // void ResetNumberOfRotationSteps (). Modified
141 // GetNumberOfRotationSteps() appropriately. Made all three functions
142 // inline (at end of this .hh file).
143 // Usage:
144 // G4Polyhedron::SetNumberOfRotationSteps
145 // (fpView -> GetViewParameters ().GetNoOfSides ());
146 // pPolyhedron = solid.CreatePolyhedron ();
147 // G4Polyhedron::ResetNumberOfRotationSteps ();
148 //
149 // 19.03.00 E.Chernyaev
150 // - added boolean operations (add, subtract, intersect) on polyhedra;
151 //
152 // 25.05.01 E.Chernyaev
153 // - added GetSurfaceArea() and GetVolume();
154 //
155 // 05.11.02 E.Chernyaev
156 // - added createTwistedTrap() and createPolyhedron();
157 //
158 // 06.03.05 J.Allison
159 // - added IsErrorBooleanProcess
160 //
161 // 20.06.05 G.Cosmo
162 // - added HepPolyhedronEllipsoid
163 //
164 // 21.10.09 J.Allison
165 // - removed IsErrorBooleanProcess (now error is returned through argument)
166 //
167 
168 #ifndef HEP_POLYHEDRON_HH
169 #define HEP_POLYHEDRON_HH
170 
171 #include "G4Types.hh"
172 #include "G4Point3D.hh"
173 #include "G4Normal3D.hh"
174 #include "G4Transform3D.hh"
175 
176 #ifndef DEFAULT_NUMBER_OF_STEPS
177 #define DEFAULT_NUMBER_OF_STEPS 24
178 #endif
179 
180 class G4Facet {
181  friend class HepPolyhedron;
182  friend std::ostream& operator<<(std::ostream&, const G4Facet &facet);
183 
184  private:
185  struct G4Edge { G4int v,f; };
187 
188  public:
189  G4Facet(G4int v1=0, G4int f1=0, G4int v2=0, G4int f2=0,
190  G4int v3=0, G4int f3=0, G4int v4=0, G4int f4=0)
191  { edge[0].v=v1; edge[0].f=f1; edge[1].v=v2; edge[1].f=f2;
192  edge[2].v=v3; edge[2].f=f3; edge[3].v=v4; edge[3].f=f4; }
193 };
194 
196  friend std::ostream& operator<<(std::ostream&, const HepPolyhedron &ph);
197 
198  protected:
203 
204  // Re-allocate memory for HepPolyhedron
205  void AllocateMemory(G4int Nvert, G4int Nface);
206 
207  // Find neighbouring facet
208  G4int FindNeighbour(G4int iFace, G4int iNode, G4int iOrder) const;
209 
210  // Find normal at node
211  G4Normal3D FindNodeNormal(G4int iFace, G4int iNode) const;
212 
213  // Create HepPolyhedron for prism with quadrilateral base
214  void CreatePrism();
215 
216  // Generate facets by revolving an edge around Z-axis
217  void RotateEdge(G4int k1, G4int k2, G4double r1, G4double r2,
218  G4int v1, G4int v2, G4int vEdge,
219  G4bool ifWholeCircle, G4int ns, G4int &kface);
220 
221  // Set side facets for the case of incomplete rotation
222  void SetSideFacets(G4int ii[4], G4int vv[4],
223  G4int *kk, G4double *r,
224  G4double dphi, G4int ns, G4int &kface);
225 
226  // Create HepPolyhedron for body of revolution around Z-axis
227  void RotateAroundZ(G4int nstep, G4double phi, G4double dphi,
228  G4int np1, G4int np2,
229  const G4double *z, G4double *r,
230  G4int nodeVis, G4int edgeVis);
231 
232  // For each edge set reference to neighbouring facet
233  void SetReferences();
234 
235  // Invert the order on nodes in facets
236  void InvertFacets();
237 
238  public:
239  // Constructor
240  HepPolyhedron() : nvert(0), nface(0), pV(0), pF(0) {}
241 
242  // Copy constructor
243  HepPolyhedron(const HepPolyhedron & from);
244 
245  // Destructor
246  virtual ~HepPolyhedron() { delete [] pV; delete [] pF; }
247 
248  // Assignment
249  HepPolyhedron & operator=(const HepPolyhedron & from);
250 
251  // Get number of vertices
252  G4int GetNoVertices() const { return nvert; }
253  G4int GetNoVerteces() const { return nvert; } // Old spelling.
254 
255  // Get number of facets
256  G4int GetNoFacets() const { return nface; }
257 
258  // Transform the polyhedron
260 
261  // Get next vertex index of the quadrilateral
262  G4bool GetNextVertexIndex(G4int & index, G4int & edgeFlag) const;
263 
264  // Get vertex by index
265  G4Point3D GetVertex(G4int index) const;
266 
267  // Get next vertex + edge visibility of the quadrilateral
268  G4bool GetNextVertex(G4Point3D & vertex, G4int & edgeFlag) const;
269 
270  // Get next vertex + edge visibility + normal of the quadrilateral
271  G4bool GetNextVertex(G4Point3D & vertex, G4int & edgeFlag,
272  G4Normal3D & normal) const;
273 
274  // Get indices of the next edge with indices of the faces
275  G4bool GetNextEdgeIndices(G4int & i1, G4int & i2, G4int & edgeFlag,
276  G4int & iface1, G4int & iface2) const;
277  G4bool GetNextEdgeIndeces(G4int & i1, G4int & i2, G4int & edgeFlag,
278  G4int & iface1, G4int & iface2) const
279  {return GetNextEdgeIndices(i1,i2,edgeFlag,iface1,iface2);} // Old spelling
280 
281  // Get indices of the next edge
282  G4bool GetNextEdgeIndices(G4int & i1, G4int & i2, G4int & edgeFlag) const;
283  G4bool GetNextEdgeIndeces(G4int & i1, G4int & i2, G4int & edgeFlag) const
284  {return GetNextEdgeIndices(i1,i2,edgeFlag);} // Old spelling.
285 
286  // Get next edge
287  G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag) const;
288 
289  // Get next edge
290  G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag,
291  G4int &iface1, G4int &iface2) const;
292 
293  // Get face by index
294  void GetFacet(G4int iFace, G4int &n, G4int *iNodes,
295  G4int *edgeFlags = 0, G4int *iFaces = 0) const;
296 
297  // Get face by index
298  void GetFacet(G4int iFace, G4int &n, G4Point3D *nodes,
299  G4int *edgeFlags=0, G4Normal3D *normals=0) const;
300 
301  // Get next face with normals at the nodes
302  G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=0,
303  G4Normal3D *normals=0) const;
304 
305  // Get normal of the face given by index
306  G4Normal3D GetNormal(G4int iFace) const;
307 
308  // Get unit normal of the face given by index
309  G4Normal3D GetUnitNormal(G4int iFace) const;
310 
311  // Get normal of the next face
313 
314  // Get normal of unit length of the next face
316 
317  // Boolean operations
318  HepPolyhedron add(const HepPolyhedron &p) const;
319  HepPolyhedron subtract(const HepPolyhedron &p) const;
320  HepPolyhedron intersect(const HepPolyhedron &p) const;
321 
322  // Get area of the surface of the polyhedron
323  G4double GetSurfaceArea() const;
324 
325  // Get volume of the polyhedron
326  G4double GetVolume() const;
327 
328  // Get number of steps for whole circle
330 
331  // Set number of steps for whole circle
332  static void SetNumberOfRotationSteps(G4int n);
333 
334  // Reset number of steps for whole circle to default value
335  static void ResetNumberOfRotationSteps();
336 
347  const G4double xy1[][2], const G4double xy2[][2]);
348 
366  G4int createPolyhedron(G4int Nnodes, G4int Nfaces,
367  const G4double xyz[][3], const G4int faces[][4]);
368 };
369 
371 {
372  public:
374  G4double Dy1, G4double Dy2, G4double Dz);
375  virtual ~HepPolyhedronTrd2();
376 };
377 
379 {
380  public:
382  G4double Dy, G4double Dz);
383  virtual ~HepPolyhedronTrd1();
384 };
385 
387 {
388  public:
390  virtual ~HepPolyhedronBox();
391 };
392 
394 {
395  public:
397  G4double Dy1,
398  G4double Dx1, G4double Dx2, G4double Alp1,
399  G4double Dy2,
400  G4double Dx3, G4double Dx4, G4double Alp2);
401  virtual ~HepPolyhedronTrap();
402 };
403 
405 {
406  public:
408  G4double Alpha, G4double Theta, G4double Phi);
409  virtual ~HepPolyhedronPara();
410 };
411 
413 {
414  public:
416  G4double r2,
417  G4double dz,
418  G4double Phi1,
419  G4double Dphi);
420  virtual ~HepPolyhedronParaboloid();
421 };
422 
424 {
425  public:
427  G4double r2,
428  G4double tan1,
429  G4double tan2,
430  G4double halfZ);
431  virtual ~HepPolyhedronHype();
432 };
433 
435 {
436  public:
437  HepPolyhedronCons(G4double Rmn1, G4double Rmx1,
438  G4double Rmn2, G4double Rmx2, G4double Dz,
439  G4double Phi1, G4double Dphi);
440  virtual ~HepPolyhedronCons();
441 };
442 
444 {
445  public:
446  HepPolyhedronCone(G4double Rmn1, G4double Rmx1,
447  G4double Rmn2, G4double Rmx2, G4double Dz);
448  virtual ~HepPolyhedronCone();
449 };
450 
452 {
453  public:
455  G4double Phi1, G4double Dphi);
456  virtual ~HepPolyhedronTubs();
457 };
458 
460 {
461  public:
462  HepPolyhedronTube (G4double Rmin, G4double Rmax, G4double Dz);
463  virtual ~HepPolyhedronTube();
464 };
465 
467 {
468  public:
469  HepPolyhedronPgon(G4double phi, G4double dphi, G4int npdv, G4int nz,
470  const G4double *z,
471  const G4double *rmin,
472  const G4double *rmax);
473  virtual ~HepPolyhedronPgon();
474 };
475 
477 {
478  public:
480  const G4double *z,
481  const G4double *rmin,
482  const G4double *rmax);
483  virtual ~HepPolyhedronPcon();
484 };
485 
487 {
488  public:
490  G4double phi, G4double dphi,
491  G4double the, G4double dthe);
492  virtual ~HepPolyhedronSphere();
493 };
494 
496 {
497  public:
498  HepPolyhedronTorus(G4double rmin, G4double rmax, G4double rtor,
499  G4double phi, G4double dphi);
500  virtual ~HepPolyhedronTorus();
501 };
502 
504 {
505  public:
507  G4double zcut1, G4double zcut2);
508  virtual ~HepPolyhedronEllipsoid();
509 };
510 
512 {
513  public:
515  G4double zcut1);
517 };
518 
519 #endif /* HEP_POLYHEDRON_HH */
HepPolyhedronTubs(G4double Rmin, G4double Rmax, G4double Dz, G4double Phi1, G4double Dphi)
G4Point3D GetVertex(G4int index) const
Float_t f4
HepPolyhedronTorus(G4double rmin, G4double rmax, G4double rtor, G4double phi, G4double dphi)
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
HepPolyhedronPcon(G4double phi, G4double dphi, G4int nz, const G4double *z, const G4double *rmin, const G4double *rmax)
virtual ~HepPolyhedronBox()
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=0, G4Normal3D *normals=0) const
G4Normal3D FindNodeNormal(G4int iFace, G4int iNode) const
HepPolyhedron subtract(const HepPolyhedron &p) const
G4double GetSurfaceArea() const
virtual ~HepPolyhedronCone()
HepPolyhedronCone(G4double Rmn1, G4double Rmx1, G4double Rmn2, G4double Rmx2, G4double Dz)
static void ResetNumberOfRotationSteps()
HepPolyhedronTrap(G4double Dz, G4double Theta, G4double Phi, G4double Dy1, G4double Dx1, G4double Dx2, G4double Alp1, G4double Dy2, G4double Dx3, G4double Dx4, G4double Alp2)
HepPolyhedronEllipticalCone(G4double dx, G4double dy, G4double z, G4double zcut1)
void RotateAroundZ(G4int nstep, G4double phi, G4double dphi, G4int np1, G4int np2, const G4double *z, G4double *r, G4int nodeVis, G4int edgeVis)
const char * p
Definition: xmltok.h:285
Double_t z
virtual ~HepPolyhedronTrd1()
HepPolyhedronTube(G4double Rmin, G4double Rmax, G4double Dz)
HepPolyhedronBox(G4double Dx, G4double Dy, G4double Dz)
G4Edge edge[4]
static G4ThreadLocal G4int fNumberOfRotationSteps
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4int FindNeighbour(G4int iFace, G4int iNode, G4int iOrder) const
Float_t f2
friend std::ostream & operator<<(std::ostream &, const G4Facet &facet)
HepPolyhedronEllipsoid(G4double dx, G4double dy, G4double dz, G4double zcut1, G4double zcut2)
#define G4ThreadLocal
Definition: tls.hh:69
virtual ~HepPolyhedronTrd2()
static void SetNumberOfRotationSteps(G4int n)
virtual ~HepPolyhedronPgon()
G4int GetNoVerteces() const
virtual ~HepPolyhedronPcon()
G4Point3D * pV
G4bool GetNextVertexIndex(G4int &index, G4int &edgeFlag) const
virtual ~HepPolyhedron()
Float_t f3
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4Facet * pF
HepPolyhedron add(const HepPolyhedron &p) const
G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag) const
G4Facet(G4int v1=0, G4int f1=0, G4int v2=0, G4int f2=0, G4int v3=0, G4int f3=0, G4int v4=0, G4int f4=0)
void GetFacet(G4int iFace, G4int &n, G4int *iNodes, G4int *edgeFlags=0, G4int *iFaces=0) const
HepPolyhedronPgon(G4double phi, G4double dphi, G4int npdv, G4int nz, const G4double *z, const G4double *rmin, const G4double *rmax)
HepPolyhedron & Transform(const G4Transform3D &t)
void RotateEdge(G4int k1, G4int k2, G4double r1, G4double r2, G4int v1, G4int v2, G4int vEdge, G4bool ifWholeCircle, G4int ns, G4int &kface)
void SetReferences()
G4bool GetNextUnitNormal(G4Normal3D &normal) const
virtual ~HepPolyhedronTubs()
virtual ~HepPolyhedronTube()
G4bool GetNextVertex(G4Point3D &vertex, G4int &edgeFlag) const
void SetSideFacets(G4int ii[4], G4int vv[4], G4int *kk, G4double *r, G4double dphi, G4int ns, G4int &kface)
Float_t f1
virtual ~HepPolyhedronPara()
virtual ~HepPolyhedronTorus()
G4bool GetNextEdgeIndeces(G4int &i1, G4int &i2, G4int &edgeFlag, G4int &iface1, G4int &iface2) const
static G4int GetNumberOfRotationSteps()
G4int GetNoVertices() const
G4bool GetNextNormal(G4Normal3D &normal) const
G4Normal3D GetUnitNormal(G4int iFace) const
int G4int
Definition: G4Types.hh:78
G4int createTwistedTrap(G4double Dz, const G4double xy1[][2], const G4double xy2[][2])
friend std::ostream & operator<<(std::ostream &, const HepPolyhedron &ph)
G4Normal3D GetNormal(G4int iFace) const
virtual ~HepPolyhedronCons()
HepPolyhedronSphere(G4double rmin, G4double rmax, G4double phi, G4double dphi, G4double the, G4double dthe)
virtual ~HepPolyhedronTrap()
HepPolyhedronParaboloid(G4double r1, G4double r2, G4double dz, G4double Phi1, G4double Dphi)
Char_t n[5]
virtual ~HepPolyhedronSphere()
G4double GetVolume() const
void AllocateMemory(G4int Nvert, G4int Nface)
HepPolyhedronTrd1(G4double Dx1, G4double Dx2, G4double Dy, G4double Dz)
G4int GetNoFacets() const
G4bool GetNextEdgeIndices(G4int &i1, G4int &i2, G4int &edgeFlag, G4int &iface1, G4int &iface2) const
G4bool GetNextEdgeIndeces(G4int &i1, G4int &i2, G4int &edgeFlag) const
HepPolyhedron intersect(const HepPolyhedron &p) const
#define ns
Definition: xmlparse.cc:614
virtual ~HepPolyhedronHype()
HepPolyhedronPara(G4double Dx, G4double Dy, G4double Dz, G4double Alpha, G4double Theta, G4double Phi)
HepPolyhedron & operator=(const HepPolyhedron &from)
HepPolyhedronCons(G4double Rmn1, G4double Rmx1, G4double Rmn2, G4double Rmx2, G4double Dz, G4double Phi1, G4double Dphi)
HepPolyhedronTrd2(G4double Dx1, G4double Dx2, G4double Dy1, G4double Dy2, G4double Dz)
HepPolyhedronHype(G4double r1, G4double r2, G4double tan1, G4double tan2, G4double halfZ)