Geant4_10
HepPolyhedron.h
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 //
27 // $Id: HepPolyhedron.h 66872 2013-01-15 01:25:57Z japost $
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 indeces 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 // GetNextEdgeIndeces (i1,i2,edgeFlag) - get indeces of the next edge;
83 // returns false for the last edge;
84 // GetNextEdgeIndeces (i1,i2,edgeFlag,iface1,iface2) - get indeces of
85 // the next edge with indeces 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 indeces
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, GetNextEdgeIndeces, 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; };
186  G4Edge edge[4];
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 
254  // Get number of facets
255  G4int GetNoFacets() const { return nface; }
256 
257  // Transform the polyhedron
259 
260  // Get next vertex index of the quadrilateral
261  G4bool GetNextVertexIndex(G4int & index, G4int & edgeFlag) const;
262 
263  // Get vertex by index
265 
266  // Get next vertex + edge visibility of the quadrilateral
267  G4bool GetNextVertex(G4Point3D & vertex, G4int & edgeFlag) const;
268 
269  // Get next vertex + edge visibility + normal of the quadrilateral
270  G4bool GetNextVertex(G4Point3D & vertex, G4int & edgeFlag,
271  G4Normal3D & normal) const;
272 
273  // Get indeces of the next edge with indeces of the faces
274  G4bool GetNextEdgeIndeces(G4int & i1, G4int & i2, G4int & edgeFlag,
275  G4int & iface1, G4int & iface2) const;
276 
277  // Get indeces of the next edge
278  G4bool GetNextEdgeIndeces(G4int & i1, G4int & i2, G4int & edgeFlag) const;
279 
280  // Get next edge
281  G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag) const;
282 
283  // Get next edge
284  G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag,
285  G4int &iface1, G4int &iface2) const;
286 
287  // Get face by index
288  void GetFacet(G4int iFace, G4int &n, G4int *iNodes,
289  G4int *edgeFlags = 0, G4int *iFaces = 0) const;
290 
291  // Get face by index
292  void GetFacet(G4int iFace, G4int &n, G4Point3D *nodes,
293  G4int *edgeFlags=0, G4Normal3D *normals=0) const;
294 
295  // Get next face with normals at the nodes
296  G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=0,
297  G4Normal3D *normals=0) const;
298 
299  // Get normal of the face given by index
300  G4Normal3D GetNormal(G4int iFace) const;
301 
302  // Get unit normal of the face given by index
303  G4Normal3D GetUnitNormal(G4int iFace) const;
304 
305  // Get normal of the next face
306  G4bool GetNextNormal(G4Normal3D &normal) const;
307 
308  // Get normal of unit length of the next face
309  G4bool GetNextUnitNormal(G4Normal3D &normal) const;
310 
311  // Boolean operations
312  HepPolyhedron add(const HepPolyhedron &p) const;
313  HepPolyhedron subtract(const HepPolyhedron &p) const;
314  HepPolyhedron intersect(const HepPolyhedron &p) const;
315 
316  // Get area of the surface of the polyhedron
317  G4double GetSurfaceArea() const;
318 
319  // Get volume of the polyhedron
320  G4double GetVolume() const;
321 
322  // Get number of steps for whole circle
324 
325  // Set number of steps for whole circle
326  static void SetNumberOfRotationSteps(G4int n);
327 
328  // Reset number of steps for whole circle to default value
329  static void ResetNumberOfRotationSteps();
330 
341  const G4double xy1[][2], const G4double xy2[][2]);
342 
360  G4int createPolyhedron(G4int Nnodes, G4int Nfaces,
361  const G4double xyz[][3], const G4int faces[][4]);
362 };
363 
365 {
366  public:
368  G4double Dy1, G4double Dy2, G4double Dz);
369  virtual ~HepPolyhedronTrd2();
370 };
371 
373 {
374  public:
376  G4double Dy, G4double Dz);
377  virtual ~HepPolyhedronTrd1();
378 };
379 
381 {
382  public:
384  virtual ~HepPolyhedronBox();
385 };
386 
388 {
389  public:
391  G4double Dy1,
392  G4double Dx1, G4double Dx2, G4double Alp1,
393  G4double Dy2,
394  G4double Dx3, G4double Dx4, G4double Alp2);
395  virtual ~HepPolyhedronTrap();
396 };
397 
399 {
400  public:
402  G4double Alpha, G4double Theta, G4double Phi);
403  virtual ~HepPolyhedronPara();
404 };
405 
407 {
408  public:
410  G4double r2,
411  G4double dz,
412  G4double Phi1,
413  G4double Dphi);
414  virtual ~HepPolyhedronParaboloid();
415 };
416 
418 {
419  public:
421  G4double r2,
422  G4double tan1,
423  G4double tan2,
424  G4double halfZ);
425  virtual ~HepPolyhedronHype();
426 };
427 
429 {
430  public:
431  HepPolyhedronCons(G4double Rmn1, G4double Rmx1,
432  G4double Rmn2, G4double Rmx2, G4double Dz,
433  G4double Phi1, G4double Dphi);
434  virtual ~HepPolyhedronCons();
435 };
436 
438 {
439  public:
440  HepPolyhedronCone(G4double Rmn1, G4double Rmx1,
441  G4double Rmn2, G4double Rmx2, G4double Dz);
442  virtual ~HepPolyhedronCone();
443 };
444 
446 {
447  public:
449  G4double Phi1, G4double Dphi);
450  virtual ~HepPolyhedronTubs();
451 };
452 
454 {
455  public:
456  HepPolyhedronTube (G4double Rmin, G4double Rmax, G4double Dz);
457  virtual ~HepPolyhedronTube();
458 };
459 
461 {
462  public:
463  HepPolyhedronPgon(G4double phi, G4double dphi, G4int npdv, G4int nz,
464  const G4double *z,
465  const G4double *rmin,
466  const G4double *rmax);
467  virtual ~HepPolyhedronPgon();
468 };
469 
471 {
472  public:
474  const G4double *z,
475  const G4double *rmin,
476  const G4double *rmax);
477  virtual ~HepPolyhedronPcon();
478 };
479 
481 {
482  public:
484  G4double phi, G4double dphi,
485  G4double the, G4double dthe);
486  virtual ~HepPolyhedronSphere();
487 };
488 
490 {
491  public:
492  HepPolyhedronTorus(G4double rmin, G4double rmax, G4double rtor,
493  G4double phi, G4double dphi);
494  virtual ~HepPolyhedronTorus();
495 };
496 
498 {
499  public:
501  G4double zcut1, G4double zcut2);
502  virtual ~HepPolyhedronEllipsoid();
503 };
504 
506 {
507  public:
509  G4double zcut1);
511 };
512 
513 #endif /* HEP_POLYHEDRON_HH */
virtual ~HepPolyhedronTorus()
void AllocateMemory(G4int Nvert, G4int Nface)
void GetFacet(G4int iFace, G4int &n, G4int *iNodes, G4int *edgeFlags=0, G4int *iFaces=0) const
HepPolyhedronEllipticalCone(G4double dx, G4double dy, G4double z, G4double zcut1)
HepPolyhedronPgon(G4double phi, G4double dphi, G4int npdv, G4int nz, const G4double *z, const G4double *rmin, const G4double *rmax)
G4Normal3D GetUnitNormal(G4int iFace) const
G4Facet * pF
virtual ~HepPolyhedronSphere()
G4int createTwistedTrap(G4double Dz, const G4double xy1[][2], const G4double xy2[][2])
HepPolyhedronParaboloid(G4double r1, G4double r2, G4double dz, G4double Phi1, G4double Dphi)
Int_t index
Definition: macro.C:9
G4bool GetNextUnitNormal(G4Normal3D &normal) const
void SetSideFacets(G4int ii[4], G4int vv[4], G4int *kk, G4double *r, G4double dphi, G4int ns, G4int &kface)
G4bool GetNextNormal(G4Normal3D &normal) const
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
G4bool GetNextEdgeIndeces(G4int &i1, G4int &i2, G4int &edgeFlag, G4int &iface1, G4int &iface2) const
G4bool GetNextEdge(G4Point3D &p1, G4Point3D &p2, G4int &edgeFlag) const
const char * p
Definition: xmltok.h:285
friend std::ostream & operator<<(std::ostream &, const HepPolyhedron &ph)
virtual ~HepPolyhedronPara()
virtual ~HepPolyhedronCone()
Float_t f4
virtual ~HepPolyhedronTrap()
HepPolyhedronTrd1(G4double Dx1, G4double Dx2, G4double Dy, G4double Dz)
HepPolyhedron subtract(const HepPolyhedron &p) const
HepPolyhedronSphere(G4double rmin, G4double rmax, G4double phi, G4double dphi, G4double the, G4double dthe)
HepPolyhedron intersect(const HepPolyhedron &p) const
G4Normal3D FindNodeNormal(G4int iFace, G4int iNode) const
#define G4ThreadLocal
Definition: tls.hh:52
HepPolyhedron & Transform(const G4Transform3D &t)
int G4int
Definition: G4Types.hh:78
HepPolyhedronTorus(G4double rmin, G4double rmax, G4double rtor, G4double phi, G4double dphi)
TFile f
Definition: plotHisto.C:6
virtual ~HepPolyhedronTubs()
HepPolyhedronEllipsoid(G4double dx, G4double dy, G4double dz, G4double zcut1, G4double zcut2)
G4int FindNeighbour(G4int iFace, G4int iNode, G4int iOrder) const
virtual ~HepPolyhedronHype()
Char_t n[5]
Float_t f1
HepPolyhedronCone(G4double Rmn1, G4double Rmx1, G4double Rmn2, G4double Rmx2, G4double Dz)
HepPolyhedronCons(G4double Rmn1, G4double Rmx1, G4double Rmn2, G4double Rmx2, G4double Dz, G4double Phi1, G4double Dphi)
static G4ThreadLocal G4int fNumberOfRotationSteps
Float_t f2
bool G4bool
Definition: G4Types.hh:79
HepPolyhedron add(const HepPolyhedron &p) const
G4Point3D * pV
G4Normal3D GetNormal(G4int iFace) const
virtual ~HepPolyhedronTube()
virtual ~HepPolyhedron()
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=0, G4Normal3D *normals=0) const
tuple v
Definition: test.py:18
HepPolyhedronPara(G4double Dx, G4double Dy, G4double Dz, G4double Alpha, G4double Theta, G4double Phi)
G4int GetNoVertices() const
Float_t f3
HepPolyhedronHype(G4double r1, G4double r2, G4double tan1, G4double tan2, G4double halfZ)
virtual ~HepPolyhedronTrd2()
HepPolyhedronBox(G4double Dx, G4double Dy, G4double Dz)
void RotateAroundZ(G4int nstep, G4double phi, G4double dphi, G4int np1, G4int np2, const G4double *z, G4double *r, G4int nodeVis, G4int edgeVis)
jump r
Definition: plot.C:36
void SetReferences()
static G4int GetNumberOfRotationSteps()
HepPolyhedronPcon(G4double phi, G4double dphi, G4int nz, const G4double *z, const G4double *rmin, const G4double *rmax)
virtual ~HepPolyhedronBox()
HepPolyhedronTubs(G4double Rmin, G4double Rmax, G4double Dz, G4double Phi1, G4double Dphi)
HepPolyhedronTrap(G4double Dz, G4double Theta, G4double Phi, G4double Dy1, G4double Dx1, G4double Dx2, G4double Alp1, G4double Dy2, G4double Dx3, G4double Dx4, G4double Alp2)
HepPolyhedron & operator=(const HepPolyhedron &from)
tuple z
Definition: test.py:28
G4bool GetNextVertexIndex(G4int &index, G4int &edgeFlag) const
HepPolyhedronTrd2(G4double Dx1, G4double Dx2, G4double Dy1, G4double Dy2, G4double Dz)
void RotateEdge(G4int k1, G4int k2, G4double r1, G4double r2, G4int v1, G4int v2, G4int vEdge, G4bool ifWholeCircle, G4int ns, G4int &kface)
virtual ~HepPolyhedronPgon()
static void SetNumberOfRotationSteps(G4int n)
friend std::ostream & operator<<(std::ostream &, const G4Facet &facet)
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)
virtual ~HepPolyhedronTrd1()
double G4double
Definition: G4Types.hh:76
G4double GetVolume() const
G4bool GetNextVertex(G4Point3D &vertex, G4int &edgeFlag) const
virtual ~HepPolyhedronCons()
G4double GetSurfaceArea() const
G4int GetNoFacets() const
#define ns
Definition: xmlparse.cc:597
HepPolyhedronTube(G4double Rmin, G4double Rmax, G4double Dz)
G4Point3D GetVertex(G4int index) const
virtual ~HepPolyhedronPcon()
static void ResetNumberOfRotationSteps()