VTK  9.2.6
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
41class vtkIdTypeArray;
42class vtkCellArray;
43class vtkTriangle;
44class vtkQuad;
45class vtkTetra;
46class vtkPolygon;
47class vtkLine;
48class vtkPointIdMap;
49class vtkIdToIdVectorMapType;
50class vtkIdToIdMapType;
51class vtkEdgeTable;
52class vtkPolyData;
53class vtkCellLocator;
54class vtkGenericCell;
55class vtkPointLocator;
56
57class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
58{
59public:
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 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
79 {
80 vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
81 return 0;
82 }
84 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
85 {
86 vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
87 }
89 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
90 {
91 vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
92 return 0;
93 }
95 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
96 {
97 vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
98 return 0;
99 }
100 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
102 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
103 {
104 vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
105 return 0;
106 }
107 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
108 {
109 vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
110 return false;
111 }
113
117 double* GetParametricCoords() override;
118
122 int GetCellType() override { return VTK_POLYHEDRON; }
123
127 int RequiresInitialization() override { return 1; }
128 void Initialize() override;
129
131
135 int GetNumberOfEdges() override;
136 vtkCell* GetEdge(int) override;
137 int GetNumberOfFaces() override;
138 vtkCell* GetFace(int faceId) override;
140
146 void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
147 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
148 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
149
159 void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
160 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
161 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
162
170 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
171 double& dist2, double weights[]) override;
172
177 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
178
185 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
186 double pcoords[3], int& subId) override;
187
203 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
204
213 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
214
219 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
220
225 int GetParametricCenter(double pcoords[3]) override;
226
230 int IsPrimaryCell() override { return 1; }
231
233
238 void InterpolateFunctions(const double x[3], double* sf) override;
239 void InterpolateDerivs(const double x[3], double* derivs) override;
241
243
251 int RequiresExplicitFaceRepresentation() override { return 1; }
252 void SetFaces(vtkIdType* faces) override;
253 vtkIdType* GetFaces() override;
255
262 int IsInside(const double x[3], double tolerance);
263
270 bool IsConvex();
271
276
277protected:
279 ~vtkPolyhedron() override;
280
281 // Internal classes for supporting operations on this cell
287 vtkIdTypeArray* GlobalFaces; // these are numbered in global id space
289
290 // vtkCell has the data members Points (x,y,z coordinates) and PointIds
291 // (global cell ids corresponding to cell canonical numbering (0,1,2,....)).
292 // These data members are implicitly organized in canonical space, i.e., where
293 // the cell point ids are (0,1,...,npts-1). The PointIdMap maps global point id
294 // back to these canonoical point ids.
295 vtkPointIdMap* PointIdMap;
296
297 // If edges are needed. Note that the edge numbering is in
298 // canonical space.
299 int EdgesGenerated; // true/false
300 vtkEdgeTable* EdgeTable; // keep track of all edges
301 vtkIdTypeArray* Edges; // edge pairs kept in this list, in canonical id space
302 vtkIdTypeArray* EdgeFaces; // face pairs that comprise each edge, with the
303 // same ordering as EdgeTable
304 int GenerateEdges(); // method populates the edge table and edge array
305
306 // If faces need renumbering into canonical numbering space these members
307 // are used. When initiallly loaded, the face numbering uses global dataset
308 // ids. Once renumbered, they are converted to canonical space.
309 vtkIdTypeArray* Faces; // these are numbered in canonical id space
312
313 // Bounds management
316 void ComputeParametricCoordinate(const double x[3], double pc[3]);
317 void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
318
320
321 // Members for supporting geometric operations
331
332 // Members used in GetPointToIncidentFaces
335
336private:
337 vtkPolyhedron(const vtkPolyhedron&) = delete;
338 void operator=(const vtkPolyhedron&) = delete;
339};
340
341//----------------------------------------------------------------------------
342inline int vtkPolyhedron::GetParametricCenter(double pcoords[3])
343{
344 pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
345 return 0;
346}
347
348#endif
abstract class to specify 3D cell interface
Definition vtkCell3D.h:39
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:36
octree-based spatial search object to quickly locate cells
abstract class to specify cell behavior
Definition vtkCell.h:58
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
keep track of edges (edge is pair of integer id's)
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:31
represent and manipulate point attribute data
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
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.
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.
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:36
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:42
a cell that represents a triangle
Definition vtkTriangle.h:36
@ VTK_POLYHEDRON
Definition vtkCellType.h:89
int vtkIdType
Definition vtkType.h:332