VTK  9.0.3
vtkPolyhedron.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPolyhedron.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
35 #ifndef vtkPolyhedron_h
36 #define vtkPolyhedron_h
37 
38 #include "vtkCell3D.h"
39 #include "vtkCommonDataModelModule.h" // For export macro
40 
41 class vtkIdTypeArray;
42 class vtkCellArray;
43 class vtkTriangle;
44 class vtkQuad;
45 class vtkTetra;
46 class vtkPolygon;
47 class vtkLine;
48 class vtkPointIdMap;
49 class vtkIdToIdVectorMapType;
50 class vtkIdToIdMapType;
51 class vtkEdgeTable;
52 class vtkPolyData;
53 class vtkCellLocator;
54 class vtkGenericCell;
55 class vtkPointLocator;
56 
57 class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
58 {
59 public:
61 
64  static vtkPolyhedron* New();
65  vtkTypeMacro(vtkPolyhedron, vtkCell3D);
66  void PrintSelf(ostream& os, vtkIndent indent) override;
68 
70 
74  void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
75  {
76  vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
77  }
78  // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
79  VTK_LEGACY(void GetEdgePoints(int vtkNotUsed(edgeId), int*& vtkNotUsed(pts)) override {
80  vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented. "
81  << "Also note that this signature is deprecated. "
82  << "Please use GetEdgePoints(vtkIdType, const vtkIdType*& instead");
83  });
84  vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
85  {
86  vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
87  return 0;
88  }
89  // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
90  VTK_LEGACY(void GetFacePoints(int vtkNotUsed(faceId), int*& vtkNotUsed(pts)) override {
91  vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented. "
92  << "Also note that this signature is deprecated. "
93  << "Please use GetFacePoints(vtkIdType, const vtkIdType*& instead");
94  });
96  vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
97  {
98  vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
99  }
101  vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
102  {
103  vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
104  return 0;
105  }
107  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
108  {
109  vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
110  return 0;
111  }
113  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(faceIds)) override
114  {
115  vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentFaces Not Implemented");
116  return 0;
117  }
119  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
120  {
121  vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
122  return 0;
123  }
124  bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
125  {
126  vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
127  return 0;
128  }
130 
134  double* GetParametricCoords() override;
135 
139  int GetCellType() override { return VTK_POLYHEDRON; }
140 
144  int RequiresInitialization() override { return 1; }
145  void Initialize() override;
146 
148 
152  int GetNumberOfEdges() override;
153  vtkCell* GetEdge(int) override;
154  int GetNumberOfFaces() override;
155  vtkCell* GetFace(int faceId) override;
157 
163  void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
164  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
165  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
166 
176  void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
177  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
178  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
179 
187  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
188  double& dist2, double weights[]) override;
189 
194  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
195 
202  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
203  double pcoords[3], int& subId) override;
204 
220  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
221 
230  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
231 
236  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
237 
242  int GetParametricCenter(double pcoords[3]) override;
243 
247  int IsPrimaryCell() override { return 1; }
248 
250 
255  void InterpolateFunctions(const double x[3], double* sf) override;
256  void InterpolateDerivs(const double x[3], double* derivs) override;
258 
260 
268  int RequiresExplicitFaceRepresentation() override { return 1; }
269  void SetFaces(vtkIdType* faces) override;
270  vtkIdType* GetFaces() override;
272 
279  int IsInside(const double x[3], double tolerance);
280 
287  bool IsConvex();
288 
293 
294 protected:
296  ~vtkPolyhedron() override;
297 
298  // Internal classes for supporting operations on this cell
304  vtkIdTypeArray* GlobalFaces; // these are numbered in global id space
306 
307  // vtkCell has the data members Points (x,y,z coordinates) and PointIds
308  // (global cell ids corresponding to cell canonical numbering (0,1,2,....)).
309  // These data members are implicitly organized in canonical space, i.e., where
310  // the cell point ids are (0,1,...,npts-1). The PointIdMap maps global point id
311  // back to these canonoical point ids.
312  vtkPointIdMap* PointIdMap;
313 
314  // If edges are needed. Note that the edge numbering is in
315  // canonical space.
316  int EdgesGenerated; // true/false
317  vtkEdgeTable* EdgeTable; // keep track of all edges
318  vtkIdTypeArray* Edges; // edge pairs kept in this list, in canonical id space
319  vtkIdTypeArray* EdgeFaces; // face pairs that comprise each edge, with the
320  // same ordering as EdgeTable
321  int GenerateEdges(); // method populates the edge table and edge array
322 
323  // If faces need renumbering into canonical numbering space these members
324  // are used. When initiallly loaded, the face numbering uses global dataset
325  // ids. Once renumbered, they are converted to canonical space.
326  vtkIdTypeArray* Faces; // these are numbered in canonical id space
329 
330  // Bounds management
333  void ComputeParametricCoordinate(const double x[3], double pc[3]);
334  void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
335 
336  // Members for supporting geometric operations
346 
347 private:
348  vtkPolyhedron(const vtkPolyhedron&) = delete;
349  void operator=(const vtkPolyhedron&) = delete;
350 };
351 
352 //----------------------------------------------------------------------------
353 inline int vtkPolyhedron::GetParametricCenter(double pcoords[3])
354 {
355  pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
356  return 0;
357 }
358 
359 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
virtual vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts)=0
Get the list of vertices that define a face.
virtual void GetEdgePoints(vtkIdType edgeId, const vtkIdType *&pts)=0
Get the pair of vertices that define an edge.
object to represent cell connectivity
Definition: vtkCellArray.h:180
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
octree-based spatial search object to quickly locate cells
abstract class to specify cell behavior
Definition: vtkCell.h:57
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
keep track of edges (edge is pair of integer id's)
Definition: vtkEdgeTable.h:41
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:31
dynamic, self-adjusting array of vtkIdType
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:34
cell represents a 1D line
Definition: vtkLine.h:30
represent and manipulate point attribute data
Definition: vtkPointData.h:32
quickly locate points in 3-space
represent and manipulate 3D points
Definition: vtkPoints.h:34
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:85
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:40
a 3D cell defined by a set of polygonal faces
Definition: vtkPolyhedron.h:58
vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(pts)) override
Definition: vtkPolyhedron.h:84
vtkPointIdMap * PointIdMap
int RequiresExplicitFaceRepresentation() override
Methods supporting the definition of faces.
vtkIdList * CellIds
vtkPolyData * GetPolyData()
Construct polydata if no one exist, then return this->PolyData.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh.
int GenerateEdges()
int GetNumberOfFaces() override
Return the number of faces in the cell.
void SetFaces(vtkIdType *faces) override
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Satisfy the vtkCell API.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Intersect the line (p1,p2) with a given tolerance tol to determine a point of intersection x[3] with ...
~vtkPolyhedron() override
void GenerateFaces()
vtkIdType GetPointToIncidentEdges(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(edgeIds)) override
vtkQuad * Quad
void ComputeParametricCoordinate(const double x[3], double pc[3])
vtkEdgeTable * EdgeTable
void ComputeBounds()
void GetEdgeToAdjacentFaces(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
Definition: vtkPolyhedron.h:95
vtkIdType GetPointToOneRingPoints(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(pts)) override
vtkIdType GetFaceToAdjacentFaces(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(faceIds)) override
vtkIdType GetPointToIncidentFaces(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(faceIds)) override
void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
Definition: vtkPolyhedron.h:74
vtkPolygon * Polygon
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
The inverse of EvaluatePosition.
vtkCell * GetFace(int faceId) override
Return the face cell from the faceId of the cell.
static vtkPolyhedron * New()
Standard new methods.
void InterpolateFunctions(const double x[3], double *sf) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
vtkCellLocator * CellLocator
void Contour(double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
Satisfy the vtkCell API.
int IsInside(const double x[3], double tolerance)
A method particular to vtkPolyhedron.
vtkIdTypeArray * EdgeFaces
vtkTriangle * Triangle
int GetParametricCenter(double pcoords[3]) override
Return the center of the cell in parametric coordinates.
vtkIdTypeArray * GlobalFaces
vtkIdTypeArray * Edges
void InterpolateDerivs(const double x[3], double *derivs) override
bool IsConvex()
Determine whether or not a polyhedron is convex.
void ConstructPolyData()
int GetNumberOfEdges() override
A polyhedron is represented internally by a set of polygonal faces.
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
vtkIdTypeArray * Faces
vtkPolyData * PolyData
int IsPrimaryCell() override
A polyhedron is a full-fledged primary cell.
vtkLine * Line
vtkTetra * Tetra
void ConstructLocator()
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Computes derivatives at the point specified by the parameter coordinate.
bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
vtkGenericCell * Cell
void ComputePositionFromParametricCoordinate(const double pc[3], double x[3])
vtkCellArray * Polys
vtkIdType * GetFaces() override
int GetCellType() override
See the vtkCell API for descriptions of these methods.
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Find the boundary face closest to the point defined by the pcoords[3] and subId of the cell (subId ca...
vtkIdTypeArray * FaceLocations
void Clip(double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Satisfy the vtkCell API.
void Initialize() override
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:36
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:42
a cell that represents a triangle
Definition: vtkTriangle.h:36
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_POLYHEDRON
Definition: vtkCellType.h:88
int vtkIdType
Definition: vtkType.h:338