VTK  9.1.0
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=========================================================================*/
74#ifndef vtkPolyhedron_h
75#define vtkPolyhedron_h
76
77#include "vtkCell3D.h"
78#include "vtkCommonDataModelModule.h" // For export macro
79
80class vtkIdTypeArray;
81class vtkCellArray;
82class vtkTriangle;
83class vtkQuad;
84class vtkTetra;
85class vtkPolygon;
86class vtkLine;
87class vtkPointIdMap;
88class vtkIdToIdVectorMapType;
89class vtkIdToIdMapType;
90class vtkEdgeTable;
91class vtkPolyData;
92class vtkCellLocator;
93class vtkGenericCell;
94class vtkPointLocator;
95
96class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
97{
98public:
100
104 vtkTypeMacro(vtkPolyhedron, vtkCell3D);
105 void PrintSelf(ostream& os, vtkIndent indent) override;
107
109
113 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
114 {
115 vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
116 }
117 // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
118 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkPolyhedron::GetEdgePoints(vtkIdType, const vtkIdType*&)")
119 void GetEdgePoints(int vtkNotUsed(edgeId), int*& vtkNotUsed(pts)) override
120 {
121 vtkErrorMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented. "
122 "Also note that this signature is deprecated. "
123 "Please use GetEdgePoints(vtkIdType, const vtkIdType*& instead");
124 };
125 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
126 {
127 vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
128 return 0;
129 }
130 // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
131 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkPolyhedron::GetFacePoints(vtkIdType, const vtkIdType*&)")
132 void GetFacePoints(int vtkNotUsed(faceId), int*& vtkNotUsed(pts)) override
133 {
134 vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented. "
135 "Also note that this signature is deprecated. "
136 "Please use GetFacePoints(vtkIdType, const vtkIdType*& instead");
137 };
139 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
140 {
141 vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
142 }
144 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
145 {
146 vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
147 return 0;
148 }
150 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
151 {
152 vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
153 return 0;
154 }
155 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
157 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
158 {
159 vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
160 return 0;
161 }
162 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
163 {
164 vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
165 return false;
166 }
168
172 double* GetParametricCoords() override;
173
177 int GetCellType() override { return VTK_POLYHEDRON; }
178
182 int RequiresInitialization() override { return 1; }
183 void Initialize() override;
184
186
190 int GetNumberOfEdges() override;
191 vtkCell* GetEdge(int) override;
192 int GetNumberOfFaces() override;
193 vtkCell* GetFace(int faceId) override;
195
201 void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
202 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
203 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
204
214 void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
215 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
216 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
217
225 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
226 double& dist2, double weights[]) override;
227
232 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
233
240 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
241 double pcoords[3], int& subId) override;
242
258 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
259
268 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
269
274 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
275
280 int GetParametricCenter(double pcoords[3]) override;
281
285 int IsPrimaryCell() override { return 1; }
286
288
293 void InterpolateFunctions(const double x[3], double* sf) override;
294 void InterpolateDerivs(const double x[3], double* derivs) override;
296
298
306 int RequiresExplicitFaceRepresentation() override { return 1; }
307 void SetFaces(vtkIdType* faces) override;
308 vtkIdType* GetFaces() override;
310
317 int IsInside(const double x[3], double tolerance);
318
325 bool IsConvex();
326
331
332protected:
334 ~vtkPolyhedron() override;
335
336 // Internal classes for supporting operations on this cell
342 vtkIdTypeArray* GlobalFaces; // these are numbered in global id space
344
345 // vtkCell has the data members Points (x,y,z coordinates) and PointIds
346 // (global cell ids corresponding to cell canonical numbering (0,1,2,....)).
347 // These data members are implicitly organized in canonical space, i.e., where
348 // the cell point ids are (0,1,...,npts-1). The PointIdMap maps global point id
349 // back to these canonoical point ids.
350 vtkPointIdMap* PointIdMap;
351
352 // If edges are needed. Note that the edge numbering is in
353 // canonical space.
354 int EdgesGenerated; // true/false
355 vtkEdgeTable* EdgeTable; // keep track of all edges
356 vtkIdTypeArray* Edges; // edge pairs kept in this list, in canonical id space
357 vtkIdTypeArray* EdgeFaces; // face pairs that comprise each edge, with the
358 // same ordering as EdgeTable
359 int GenerateEdges(); // method populates the edge table and edge array
360
361 // If faces need renumbering into canonical numbering space these members
362 // are used. When initiallly loaded, the face numbering uses global dataset
363 // ids. Once renumbered, they are converted to canonical space.
364 vtkIdTypeArray* Faces; // these are numbered in canonical id space
367
368 // Bounds management
371 void ComputeParametricCoordinate(const double x[3], double pc[3]);
372 void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
373
375
376 // Members for supporting geometric operations
386
387 // Members used in GetPointToIncidentFaces
390
391private:
392 vtkPolyhedron(const vtkPolyhedron&) = delete;
393 void operator=(const vtkPolyhedron&) = delete;
394};
395
396//----------------------------------------------------------------------------
397inline int vtkPolyhedron::GetParametricCenter(double pcoords[3])
398{
399 pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
400 return 0;
401}
402
403#endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:40
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
octree-based spatial search object to quickly locate cells
abstract class to specify cell behavior
Definition: vtkCell.h:147
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
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:140
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:113
cell represents a 1D line
Definition: vtkLine.h:140
represent and manipulate point attribute data
Definition: vtkPointData.h:142
quickly locate points in 3-space
represent and manipulate 3D points
Definition: vtkPoints.h:143
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:149
a 3D cell defined by a set of polygonal faces
Definition: vtkPolyhedron.h:97
vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
vtkPointIdMap * PointIdMap
int RequiresExplicitFaceRepresentation() override
Methods supporting the definition of faces.
vtkIdList * CellIds
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh.
int GenerateEdges()
int GetNumberOfFaces() override
A polyhedron is represented internally by a set of polygonal faces.
void SetFaces(vtkIdType *faces) override
Methods supporting the definition of faces.
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
See vtkCell3D API for description of these methods.
vtkQuad * Quad
void ComputeParametricCoordinate(const double x[3], double pc[3])
vtkIdType * ValenceAtPoint
vtkEdgeTable * EdgeTable
void ComputeBounds()
void GetEdgeToAdjacentFaces(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
vtkIdType GetPointToOneRingPoints(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
vtkIdType GetFaceToAdjacentFaces(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(faceIds)) override
See vtkCell3D API for description of these methods.
void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
vtkPolygon * Polygon
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
The inverse of EvaluatePosition.
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
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
vtkIdTypeArray * Edges
void InterpolateDerivs(const double x[3], double *derivs) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
bool IsConvex()
Determine whether or not a polyhedron is convex.
void ConstructPolyData()
vtkIdType ** PointToIncidentFaces
int GetNumberOfEdges() override
A polyhedron is represented internally by a set of polygonal faces.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard new methods.
void GeneratePointToIncidentFacesAndValenceAtPoint()
vtkIdTypeArray * Faces
vtkPolyData * PolyData
int IsPrimaryCell() override
A polyhedron is a full-fledged primary cell.
vtkPolyData * GetPolyData()
Construct polydata if no one exist, then return this->PolyData.
vtkLine * Line
vtkCell * GetEdge(int) override
A polyhedron is represented internally by a set of polygonal faces.
vtkTetra * Tetra
vtkCell * GetFace(int faceId) override
A polyhedron is represented internally by a set of polygonal faces.
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
See vtkCell3D API for description of these methods.
vtkIdType * GetFaces() override
Methods supporting the definition of faces.
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
vtkGenericCell * Cell
void ComputePositionFromParametricCoordinate(const double pc[3], double x[3])
static vtkPolyhedron * New()
Standard new methods.
vtkCellArray * Polys
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:95
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:111
a cell that represents a triangle
Definition: vtkTriangle.h:145
@ VTK_POLYHEDRON
Definition: vtkCellType.h:128
#define VTK_DEPRECATED_IN_9_0_0(reason)
int vtkIdType
Definition: vtkType.h:332