Geant4 Cross Reference |
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 // Class Description: 27 // HepPolyhedron is an intermediate class between description of a shape 28 // and visualization systems. It is intended to provide some service like: 29 // - polygonization of shapes with triangulization (quadrilaterization) 30 // of complex polygons; 31 // - calculation of normals for faces and vertices; 32 // - finding result of boolean operation on polyhedra; 33 // 34 // Public constructors: 35 // 36 // HepPolyhedronBox (dx,dy,dz) 37 // - create polyhedron for Box; 38 // HepPolyhedronTrd1 (dx1,dx2,dy,dz) 39 // - create polyhedron for Trd1; 40 // HepPolyhedronTrd2 (dx1,dx2,dy1,dy2,dz) 41 // - create polyhedron for Trd2; 42 // HepPolyhedronTrap (dz,theta,phi, h1,bl1,tl1,alp1, h2,bl2,tl2,alp2) 43 // - create polyhedron for Trap; 44 // HepPolyhedronPara (dx,dy,dz,alpha,theta,phi) 45 // - create polyhedron for Para; 46 // HepPolyhedronTube (rmin,rmax,dz) 47 // - create polyhedron for Tube; 48 // HepPolyhedronTubs (rmin,rmax,dz,phi1,dphi) 49 // - create polyhedron for Tubs; 50 // HepPolyhedronCone (rmin1,rmax1,rmin2,rmax2,dz) 51 // - create polyhedron for Cone; 52 // HepPolyhedronCons (rmin1,rmax1,rmin2,rmax2,dz,phi1,dphi) 53 // - create polyhedron for Cons; 54 // HepPolyhedronPgon (phi,dphi,npdv,nz, z(*),rmin(*),rmax(*)) 55 // - create polyhedron for Pgon; 56 // HepPolyhedronPcon (phi,dphi,nz, z(*),rmin(*),rmax(*)) 57 // - create polyhedron for Pcon; 58 // HepPolyhedronSphere (rmin,rmax,phi,dphi,the,dthe) 59 // - create polyhedron for Sphere; 60 // HepPolyhedronTorus (rmin,rmax,rtor,phi,dphi) 61 // - create polyhedron for Torus; 62 // HepPolyhedronTet (p0[3],p1[3],p2[3],p3[3]) 63 // - create polyhedron for Tet; 64 // HepPolyhedronEllipsoid (dx,dy,dz,zcut1,zcut2) 65 // - create polyhedron for Ellipsoid; 66 // HepPolyhedronEllipticalCone(dx,dy,z,zcut1) 67 // - create polyhedron for Elliptical cone; 68 // HepPolyhedronParaboloid (r1,r2,dz,phi,dphi) 69 // - create polyhedron for Paraboloid; 70 // HepPolyhedronHype (r1,r2,tan1,tan2,halfz) 71 // - create polyhedron for Hype; 72 // HepPolyhedronHyperbolicMirror (a,h,r) 73 // - create polyhedron for Hyperbolic mirror; 74 // HepPolyhedronTetMesh (vector<p>) 75 // - create polyhedron for tetrahedron mesh; 76 // HepPolyhedronBoxMesh (sx,sy,sz,vector<p>) 77 // - create polyhedron for box mesh; 78 // Public functions: 79 // 80 // GetNoVertices () - returns number of vertices; 81 // GetNoFacets () - returns number of faces; 82 // GetNextVertexIndex (index,edgeFlag) - get vertex indices of the 83 // quadrilaterals in order; 84 // returns false when finished each face; 85 // GetVertex (index) - returns vertex by index; 86 // GetNextVertex (vertex,edgeFlag) - get vertices with edge visibility 87 // of the quadrilaterals in order; 88 // returns false when finished each face; 89 // GetNextVertex (vertex,edgeFlag,normal) - get vertices with edge 90 // visibility and normal of the quadrilaterals 91 // in order; returns false when finished each face; 92 // GetNextEdgeIndices (i1,i2,edgeFlag) - get indices of the next edge; 93 // returns false for the last edge; 94 // GetNextEdgeIndices (i1,i2,edgeFlag,iface1,iface2) - get indices of 95 // the next edge with indices of the faces 96 // to which the edge belongs; 97 // returns false for the last edge; 98 // GetNextEdge (p1,p2,edgeFlag) - get next edge; 99 // returns false for the last edge; 100 // GetNextEdge (p1,p2,edgeFlag,iface1,iface2) - get next edge with indices 101 // of the faces to which the edge belongs; 102 // returns false for the last edge; 103 // GetFacet (index,n,nodes,edgeFlags=0,normals=0) - get face by index; 104 // GetNextFacet (n,nodes,edgeFlags=0,normals=0) - get next face with normals 105 // at the nodes; returns false for the last face; 106 // GetNormal (index) - get normal of face given by index; 107 // GetUnitNormal (index) - get unit normal of face given by index; 108 // GetNextNormal (normal) - get normals of each face in order; 109 // returns false when finished all faces; 110 // GetNextUnitNormal (normal) - get normals of unit length of each face 111 // in order; returns false when finished all faces; 112 // GetSurfaceArea() - get surface area of the polyhedron; 113 // GetVolume() - get volume of the polyhedron; 114 // GetNumberOfRotationSteps() - get number of steps for whole circle; 115 // SetVertex(index, v) - set vertex; 116 // SetFacet(index,iv1,iv2,iv3,iv4) - set facet; 117 // SetReferences() - set references to neighbouring facets; 118 // JoinCoplanarFacets(tolerance) - join coplanar facets where it is possible 119 // InvertFacets() - invert the order on nodes in facets; 120 // SetNumberOfRotationSteps (n) - set number of steps for whole circle; 121 // ResetNumberOfRotationSteps() - reset number of steps for whole circle 122 // to default value; 123 // History: 124 // 125 // 20.06.96 Evgeni Chernyaev <Evgueni.Tcherniaev@cern.ch> - initial version 126 // 127 // 23.07.96 John Allison 128 // - added GetNoVertices, GetNoFacets, GetNextVertex, GetNextNormal 129 // 130 // 30.09.96 E.Chernyaev 131 // - added GetNextVertexIndex, GetVertex by Yasuhide Sawada 132 // - added GetNextUnitNormal, GetNextEdgeIndices, GetNextEdge 133 // - improvements: angles now expected in radians 134 // int -> G4int, double -> G4double 135 // - G4ThreeVector replaced by either G4Point3D or G4Normal3D 136 // 137 // 15.12.96 E.Chernyaev 138 // - private functions G4PolyhedronAlloc, G4PolyhedronPrism renamed 139 // to AllocateMemory and CreatePrism 140 // - added private functions GetNumberOfRotationSteps, RotateEdge, 141 // RotateAroundZ, SetReferences 142 // - rewritten G4PolyhedronCons; 143 // - added G4PolyhedronPara, ...Trap, ...Pgon, ...Pcon, ...Sphere, ...Torus, 144 // so full List of implemented shapes now looks like: 145 // BOX, TRD1, TRD2, TRAP, TUBE, TUBS, CONE, CONS, PARA, PGON, PCON, 146 // SPHERE, TORUS 147 // 148 // 01.06.97 E.Chernyaev 149 // - RotateAroundZ modified and SetSideFacets added to allow Rmin=Rmax 150 // in bodies of revolution 151 // 152 // 24.06.97 J.Allison 153 // - added static private member fNumberOfRotationSteps and static public 154 // functions void SetNumberOfRotationSteps (G4int n) and 155 // void ResetNumberOfRotationSteps (). Modified 156 // GetNumberOfRotationSteps() appropriately. Made all three functions 157 // inline (at end of this .hh file). 158 // Usage: 159 // G4Polyhedron::SetNumberOfRotationSteps 160 // (fpView -> GetViewParameters ().GetNoOfSides ()); 161 // pPolyhedron = solid.CreatePolyhedron (); 162 // G4Polyhedron::ResetNumberOfRotationSteps (); 163 // 164 // 19.03.00 E.Chernyaev 165 // - added boolean operations (add, subtract, intersect) on polyhedra; 166 // 167 // 25.05.01 E.Chernyaev 168 // - added GetSurfaceArea() and GetVolume(); 169 // 170 // 05.11.02 E.Chernyaev 171 // - added createTwistedTrap() and createPolyhedron(); 172 // 173 // 06.03.05 J.Allison 174 // - added IsErrorBooleanProcess 175 // 176 // 20.06.05 G.Cosmo 177 // - added HepPolyhedronEllipsoid 178 // 179 // 18.07.07 T.Nikitina 180 // - added HepPolyhedronParaboloid; 181 // 182 // 21.10.09 J.Allison 183 // - removed IsErrorBooleanProcess (now error is returned through argument) 184 // 185 // 22.02.20 E.Chernyaev 186 // - added HepPolyhedronTet, HepPolyhedronHyberbolicMirror 187 // 188 // 12.05.21 E.Chernyaev 189 // - added TriangulatePolygon(), RotateContourAroundZ() 190 // - added HepPolyhedronPgon, HepPolyhedronPcon given by rz-countour 191 // 192 // 26.03.22 E.Chernyaev 193 // - added HepPolyhedronTetMesh 194 // 195 // 04.04.22 E.Chernyaev 196 // - added JoinCoplanarFacets() 197 // 198 // 07.04.22 E.Chernyaev 199 // - added HepPolyhedronBoxMesh 200 201 #ifndef HEP_POLYHEDRON_HH 202 #define HEP_POLYHEDRON_HH 203 204 #include <vector> 205 #include "G4Types.hh" 206 #include "G4TwoVector.hh" 207 #include "G4ThreeVector.hh" 208 #include "G4Point3D.hh" 209 #include "G4Normal3D.hh" 210 #include "G4Transform3D.hh" 211 212 #ifndef DEFAULT_NUMBER_OF_STEPS 213 #define DEFAULT_NUMBER_OF_STEPS 24 214 #endif 215 216 class G4Facet { 217 friend class HepPolyhedron; 218 friend std::ostream& operator<<(std::ostream&, const G4Facet &facet); 219 220 private: 221 struct G4Edge { G4int v,f; }; 222 G4Edge edge[4]; 223 224 public: 225 G4Facet(G4int v1=0, G4int f1=0, G4int v2=0, G4int f2=0, 226 G4int v3=0, G4int f3=0, G4int v4=0, G4int f4=0) 227 { edge[0].v=v1; edge[0].f=f1; edge[1].v=v2; edge[1].f=f2; 228 edge[2].v=v3; edge[2].f=f3; edge[3].v=v4; edge[3].f=f4; } 229 }; 230 231 class HepPolyhedron { 232 friend std::ostream& operator<<(std::ostream&, const HepPolyhedron &ph); 233 234 protected: 235 static G4ThreadLocal G4int fNumberOfRotationSteps; 236 G4int nvert, nface; 237 G4Point3D *pV; 238 G4Facet *pF; 239 240 // Re-allocate memory for HepPolyhedron 241 void AllocateMemory(G4int Nvert, G4int Nface); 242 243 // Find neighbouring facet 244 G4int FindNeighbour(G4int iFace, G4int iNode, G4int iOrder) const; 245 246 // Find normal at node 247 G4Normal3D FindNodeNormal(G4int iFace, G4int iNode) const; 248 249 // Create HepPolyhedron for prism with quadrilateral base 250 void CreatePrism(); 251 252 // Generate facets by revolving an edge around Z-axis 253 void RotateEdge(G4int k1, G4int k2, G4double r1, G4double r2, 254 G4int v1, G4int v2, G4int vEdge, 255 G4bool ifWholeCircle, G4int ns, G4int &kface); 256 257 // Set side facets for the case of incomplete rotation 258 void SetSideFacets(G4int ii[4], G4int vv[4], 259 G4int *kk, G4double *r, 260 G4double dphi, G4int ns, G4int &kface); 261 262 // Create HepPolyhedron for body of revolution around Z-axis 263 void RotateAroundZ(G4int nstep, G4double phi, G4double dphi, 264 G4int np1, G4int np2, 265 const G4double *z, G4double *r, 266 G4int nodeVis, G4int edgeVis); 267 268 // Create HepPolyhedron for body of revolution around Z-axis 269 void RotateContourAroundZ(G4int nstep, G4double phi, G4double dphi, 270 const std::vector<G4TwoVector> &rz, 271 G4int nodeVis, G4int edgeVis); 272 273 // Triangulate closed polygon (contour) 274 G4bool TriangulatePolygon(const std::vector<G4TwoVector> &polygon, 275 std::vector<G4int> &result); 276 277 // Helper function for TriangulatePolygon() 278 G4bool CheckSnip(const std::vector<G4TwoVector> &contour, 279 G4int a, G4int b, G4int c, 280 G4int n, const G4int* V); 281 282 public: 283 // Default constructor 284 HepPolyhedron() : nvert(0), nface(0), pV(nullptr), pF(nullptr) {} 285 286 // Constructor with allocation of memory 287 HepPolyhedron(G4int Nvert, G4int Nface); 288 289 // Copy constructor 290 HepPolyhedron(const HepPolyhedron & from); 291 292 // Move constructor 293 HepPolyhedron(HepPolyhedron && from); 294 295 // Destructor 296 virtual ~HepPolyhedron() { delete [] pV; delete [] pF; } 297 298 // Assignment 299 HepPolyhedron & operator=(const HepPolyhedron & from); 300 301 // Move assignment 302 HepPolyhedron & operator=(HepPolyhedron && from); 303 304 // Get number of vertices 305 G4int GetNoVertices() const { return nvert; } 306 G4int GetNoVerteces() const { return nvert; } // Old spelling. 307 308 // Get number of facets 309 G4int GetNoFacets() const { return nface; } 310 311 // Transform the polyhedron 312 HepPolyhedron & Transform(const G4Transform3D & t); 313 314 // Get next vertex index of the quadrilateral 315 G4bool GetNextVertexIndex(G4int & index, G4int & edgeFlag) const; 316 317 // Get vertex by index 318 G4Point3D GetVertex(G4int index) const; 319 320 // Get next vertex + edge visibility of the quadrilateral 321 G4bool GetNextVertex(G4Point3D & vertex, G4int & edgeFlag) const; 322 323 // Get next vertex + edge visibility + normal of the quadrilateral 324 G4bool GetNextVertex(G4Point3D & vertex, G4int & edgeFlag, 325 G4Normal3D & normal) const; 326 327 // Get indices of the next edge with indices of the faces 328 G4bool GetNextEdgeIndices(G4int & i1, G4int & i2, G4int & edgeFlag, 329 G4int & iface1, G4int & iface2) const; 330 G4bool GetNextEdgeIndeces(G4int & i1, G4int & i2, G4int & edgeFlag, 331 G4int & iface1, G4int & iface2) const 332 {return GetNextEdgeIndices(i1,i2,edgeFlag,iface1,iface2);} // Old spelling 333 334 // Get indices of the next edge 335 G4bool GetNextEdgeIndices(G4int & i1, G4int & i2, G4int & edgeFlag) const; 336 G4bool GetNextEdgeIndeces(G4int & i1, G4int & i2, G4int & edgeFlag) const 337 {return GetNextEdgeIndices(i1,i2,edgeFlag);} // Old spelling. 338 339 // Get next edge 340 G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag) const; 341 342 // Get next edge 343 G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag, 344 G4int &iface1, G4int &iface2) const; 345 346 // Get face by index 347 void GetFacet(G4int iFace, G4int &n, G4int *iNodes, 348 G4int *edgeFlags = nullptr, G4int *iFaces = nullptr) const; 349 350 // Get face by index 351 void GetFacet(G4int iFace, G4int &n, G4Point3D *nodes, 352 G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const; 353 354 // Get next face with normals at the nodes 355 G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, 356 G4Normal3D *normals=nullptr) const; 357 358 // Get normal of the face given by index 359 G4Normal3D GetNormal(G4int iFace) const; 360 361 // Get unit normal of the face given by index 362 G4Normal3D GetUnitNormal(G4int iFace) const; 363 364 // Get normal of the next face 365 G4bool GetNextNormal(G4Normal3D &normal) const; 366 367 // Get normal of unit length of the next face 368 G4bool GetNextUnitNormal(G4Normal3D &normal) const; 369 370 // Boolean operations 371 HepPolyhedron add(const HepPolyhedron &p) const; 372 HepPolyhedron subtract(const HepPolyhedron &p) const; 373 HepPolyhedron intersect(const HepPolyhedron &p) const; 374 375 // Get area of the surface of the polyhedron 376 G4double GetSurfaceArea() const; 377 378 // Get volume of the polyhedron 379 G4double GetVolume() const; 380 381 // Get number of steps for whole circle 382 static G4int GetNumberOfRotationSteps(); 383 384 // Set vertex (1 <= index <= Nvert) 385 void SetVertex(G4int index, const G4Point3D& v); 386 387 // Set facet (1 <= index <= Nface) 388 void SetFacet(G4int index, G4int iv1, G4int iv2, G4int iv3, G4int iv4 = 0); 389 390 // For each edge set reference to neighbouring facet, 391 // call this after all vertices and facets have been set 392 void SetReferences(); 393 394 // Join couples of triangular facets to quadrangular facets 395 // where it is possible 396 void JoinCoplanarFacets(G4double tolerance); 397 398 // Invert the order on nodes in facets 399 void InvertFacets(); 400 401 // Set number of steps for whole circle 402 static void SetNumberOfRotationSteps(G4int n); 403 404 // Reset number of steps for whole circle to default value 405 static void ResetNumberOfRotationSteps(); 406 407 /** 408 * Creates polyhedron for twisted trapezoid. 409 * The trapezoid is given by two bases perpendicular to the z-axis. 410 * 411 * @param Dz half length in z 412 * @param xy1 1st base (at z = -Dz) 413 * @param xy2 2nd base (at z = +Dz) 414 * @return status of the operation - is non-zero in case of problem 415 */ 416 G4int createTwistedTrap(G4double Dz, 417 const G4double xy1[][2], const G4double xy2[][2]); 418 419 /** 420 * Creates user defined polyhedron. 421 * This function allows to the user to define arbitrary polyhedron. 422 * The faces of the polyhedron should be either triangles or planar 423 * quadrilateral. Nodes of a face are defined by indexes pointing to 424 * the elements in the xyz array. Numeration of the elements in the 425 * array starts from 1 (like in fortran). The indexes can be positive 426 * or negative. Negative sign means that the corresponding edge is 427 * invisible. The normal of the face should be directed to exterior 428 * of the polyhedron. 429 * 430 * @param Nnodes number of nodes 431 * @param Nfaces number of faces 432 * @param xyz nodes 433 * @param faces faces (quadrilaterals or triangles) 434 * @return status of the operation - is non-zero in case of problem 435 */ 436 G4int createPolyhedron(G4int Nnodes, G4int Nfaces, 437 const G4double xyz[][3], const G4int faces[][4]); 438 439 /** 440 * Calculate the unweighted mean of all the vertices in the polyhedron. Not to be 441 * confused with the polyhedron centre or centre of mass 442 * @return G4Point3D of the unweighted mean vertex position 443 */ 444 G4Point3D vertexUnweightedMean() const; 445 }; 446 447 class HepPolyhedronTrd2 : public HepPolyhedron 448 { 449 public: 450 HepPolyhedronTrd2(G4double Dx1, G4double Dx2, 451 G4double Dy1, G4double Dy2, G4double Dz); 452 ~HepPolyhedronTrd2() override; 453 }; 454 455 class HepPolyhedronTrd1 : public HepPolyhedronTrd2 456 { 457 public: 458 HepPolyhedronTrd1(G4double Dx1, G4double Dx2, 459 G4double Dy, G4double Dz); 460 ~HepPolyhedronTrd1() override; 461 }; 462 463 class HepPolyhedronBox : public HepPolyhedronTrd2 464 { 465 public: 466 HepPolyhedronBox(G4double Dx, G4double Dy, G4double Dz); 467 ~HepPolyhedronBox() override; 468 }; 469 470 class HepPolyhedronTrap : public HepPolyhedron 471 { 472 public: 473 HepPolyhedronTrap(G4double Dz, G4double Theta, G4double Phi, 474 G4double Dy1, 475 G4double Dx1, G4double Dx2, G4double Alp1, 476 G4double Dy2, 477 G4double Dx3, G4double Dx4, G4double Alp2); 478 ~HepPolyhedronTrap() override; 479 }; 480 481 class HepPolyhedronPara : public HepPolyhedronTrap 482 { 483 public: 484 HepPolyhedronPara(G4double Dx, G4double Dy, G4double Dz, 485 G4double Alpha, G4double Theta, G4double Phi); 486 ~HepPolyhedronPara() override; 487 }; 488 489 class HepPolyhedronParaboloid : public HepPolyhedron 490 { 491 public: 492 HepPolyhedronParaboloid(G4double r1, 493 G4double r2, 494 G4double dz, 495 G4double Phi1, 496 G4double Dphi); 497 ~HepPolyhedronParaboloid() override; 498 }; 499 500 class HepPolyhedronHype : public HepPolyhedron 501 { 502 public: 503 HepPolyhedronHype(G4double r1, 504 G4double r2, 505 G4double tan1, 506 G4double tan2, 507 G4double halfZ); 508 ~HepPolyhedronHype() override; 509 }; 510 511 class HepPolyhedronCons : public HepPolyhedron 512 { 513 public: 514 HepPolyhedronCons(G4double Rmn1, G4double Rmx1, 515 G4double Rmn2, G4double Rmx2, G4double Dz, 516 G4double Phi1, G4double Dphi); 517 ~HepPolyhedronCons() override; 518 }; 519 520 class HepPolyhedronCone : public HepPolyhedronCons 521 { 522 public: 523 HepPolyhedronCone(G4double Rmn1, G4double Rmx1, 524 G4double Rmn2, G4double Rmx2, G4double Dz); 525 ~HepPolyhedronCone() override; 526 }; 527 528 class HepPolyhedronTubs : public HepPolyhedronCons 529 { 530 public: 531 HepPolyhedronTubs(G4double Rmin, G4double Rmax, G4double Dz, 532 G4double Phi1, G4double Dphi); 533 ~HepPolyhedronTubs() override; 534 }; 535 536 class HepPolyhedronTube : public HepPolyhedronCons 537 { 538 public: 539 HepPolyhedronTube (G4double Rmin, G4double Rmax, G4double Dz); 540 ~HepPolyhedronTube() override; 541 }; 542 543 class HepPolyhedronPgon : public HepPolyhedron 544 { 545 public: 546 HepPolyhedronPgon(G4double phi, G4double dphi, G4int npdv, G4int nz, 547 const G4double *z, 548 const G4double *rmin, 549 const G4double *rmax); 550 HepPolyhedronPgon(G4double phi, G4double dphi, G4int npdv, 551 const std::vector<G4TwoVector> &rz); 552 ~HepPolyhedronPgon() override; 553 }; 554 555 class HepPolyhedronPcon : public HepPolyhedronPgon 556 { 557 public: 558 HepPolyhedronPcon(G4double phi, G4double dphi, G4int nz, 559 const G4double *z, 560 const G4double *rmin, 561 const G4double *rmax); 562 HepPolyhedronPcon(G4double phi, G4double dphi, 563 const std::vector<G4TwoVector> &rz); 564 ~HepPolyhedronPcon() override; 565 }; 566 567 class HepPolyhedronSphere : public HepPolyhedron 568 { 569 public: 570 HepPolyhedronSphere(G4double rmin, G4double rmax, 571 G4double phi, G4double dphi, 572 G4double the, G4double dthe); 573 ~HepPolyhedronSphere() override; 574 }; 575 576 class HepPolyhedronTorus : public HepPolyhedron 577 { 578 public: 579 HepPolyhedronTorus(G4double rmin, G4double rmax, G4double rtor, 580 G4double phi, G4double dphi); 581 ~HepPolyhedronTorus() override; 582 }; 583 584 class HepPolyhedronTet : public HepPolyhedron 585 { 586 public: 587 HepPolyhedronTet(const G4double p0[3], 588 const G4double p1[3], 589 const G4double p2[3], 590 const G4double p3[3]); 591 ~HepPolyhedronTet() override; 592 }; 593 594 class HepPolyhedronEllipsoid : public HepPolyhedron 595 { 596 public: 597 HepPolyhedronEllipsoid(G4double dx, G4double dy, G4double dz, 598 G4double zcut1, G4double zcut2); 599 ~HepPolyhedronEllipsoid() override; 600 }; 601 602 class HepPolyhedronEllipticalCone : public HepPolyhedron 603 { 604 public: 605 HepPolyhedronEllipticalCone(G4double dx, G4double dy, G4double z, 606 G4double zcut1); 607 ~HepPolyhedronEllipticalCone() override; 608 }; 609 610 class HepPolyhedronHyperbolicMirror : public HepPolyhedron 611 { 612 public: 613 HepPolyhedronHyperbolicMirror(G4double a, G4double h, G4double r); 614 ~HepPolyhedronHyperbolicMirror() override; 615 }; 616 617 class HepPolyhedronTetMesh : public HepPolyhedron 618 { 619 public: 620 HepPolyhedronTetMesh(const std::vector<G4ThreeVector>& tetrahedra); 621 ~HepPolyhedronTetMesh() override; 622 }; 623 624 class HepPolyhedronBoxMesh : public HepPolyhedron 625 { 626 public: 627 HepPolyhedronBoxMesh(G4double sizeX, G4double sizeY, G4double sizeZ, 628 const std::vector<G4ThreeVector>& positions); 629 ~HepPolyhedronBoxMesh() override; 630 }; 631 632 #endif /* HEP_POLYHEDRON_HH */ 633