VTK  9.3.0
vtkPolyhedron.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
123#ifndef vtkPolyhedron_h
124#define vtkPolyhedron_h
125
126#include "vtkCell3D.h"
127#include "vtkCommonDataModelModule.h" // For export macro
128
129VTK_ABI_NAMESPACE_BEGIN
130class vtkIdTypeArray;
131class vtkCellArray;
132class vtkTriangle;
133class vtkQuad;
134class vtkTetra;
135class vtkPolygon;
136class vtkLine;
137class vtkIdToIdVectorMapType;
138class vtkIdToIdMapType;
139class vtkEdgeTable;
140class vtkPolyData;
141class vtkCellLocator;
142class vtkGenericCell;
143class vtkPointLocator;
144
145class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
146{
147public:
148 typedef std::map<vtkIdType, vtkIdType> vtkPointIdMap;
149
151
155 vtkTypeMacro(vtkPolyhedron, vtkCell3D);
156 void PrintSelf(ostream& os, vtkIndent indent) override;
158
160
164 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
165 {
166 vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
167 }
168 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
169 {
170 vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
171 return 0;
172 }
174 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
175 {
176 vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
177 }
179 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
180 {
181 vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
182 return 0;
183 }
185 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
186 {
187 vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
188 return 0;
189 }
190 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
192 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
193 {
194 vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
195 return 0;
196 }
197 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
198 {
199 vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
200 return false;
201 }
203
207 double* GetParametricCoords() override;
208
212 int GetCellType() override { return VTK_POLYHEDRON; }
213
217 int RequiresInitialization() override { return 1; }
218
224 void Initialize() override;
225
227
231 int GetNumberOfEdges() override;
232 vtkCell* GetEdge(int) override;
233 int GetNumberOfFaces() override;
234 vtkCell* GetFace(int faceId) override;
236
242 void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
243 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
244 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
245
255 void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
256 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
257 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
258
266 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
267 double& dist2, double weights[]) override;
268
273 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
274
281 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
282 double pcoords[3], int& subId) override;
283
299 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
300
308
317 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
318
323 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
324
329 int GetParametricCenter(double pcoords[3]) override;
330
334 int IsPrimaryCell() override { return 1; }
335
337
342 void InterpolateFunctions(const double x[3], double* sf) override;
343 void InterpolateDerivs(const double x[3], double* derivs) override;
345
351 int RequiresExplicitFaceRepresentation() override { return 1; }
352
369 void SetFaces(vtkIdType* faces) override;
370
387 vtkIdType* GetFaces() override;
388
395 int IsInside(const double x[3], double tolerance);
396
403 bool IsConvex();
404
409
410protected:
412 ~vtkPolyhedron() override;
413
414 // Internal classes for supporting operations on this cell
420
421 // Filled with the SetFaces method.
422 // These faces are numbered in global id space
423 // (in the legacy vtkCellArray form)
425
426 // Filled with the SetFaces method.
427 // Used to to point to each face in the GlobalFaces array.
429
430 // vtkCell has the data members Points (x,y,z coordinates) and PointIds (global cell ids).
431 // These data members are implicitly organized in canonical space, i.e., where the cell
432 // point ids are (0,1,...,npts-1).
433 // The PointIdMap is constructed during the call of the Initialize() method and maps global
434 // point ids to the canonical point ids.
436
437 // If edges are needed. Note that the edge numbering is in canonical space.
438 int EdgesGenerated; // true/false
439 vtkEdgeTable* EdgeTable; // keep track of all edges
440 vtkIdTypeArray* Edges; // edge pairs kept in this list, in canonical id space
441 vtkIdTypeArray* EdgeFaces; // face pairs that comprise each edge, with the
442 // same ordering as EdgeTable
443 int GenerateEdges(); // method populates the edge table and edge array
444
445 // Numerous methods needs faces to be numbered in the canonical space.
446 // This method uses PointIdMap to fill the Faces member (faces described
447 // with canonical IDs) from the GlobalFaces member (faces described with
448 // global IDs).
450 vtkIdTypeArray* Faces; // These are numbered in canonical id space
451 int FacesGenerated; // True when Faces have been successfully constructed
452
453 // Bounds management
456 void ComputeParametricCoordinate(const double x[3], double pc[3]);
457 void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
458
460
461 // Members for supporting geometric operations
471
472 // Members used in GetPointToIncidentFaces
475
476private:
477 vtkPolyhedron(const vtkPolyhedron&) = delete;
478 void operator=(const vtkPolyhedron&) = delete;
479
481};
482
483//----------------------------------------------------------------------------
484inline int vtkPolyhedron::GetParametricCenter(double pcoords[3])
485{
486 pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
487 return 0;
488}
489
490VTK_ABI_NAMESPACE_END
491#endif
abstract class to specify 3D cell interface
Definition vtkCell3D.h:28
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:31
octree-based spatial search object to quickly locate cells
abstract class to specify cell behavior
Definition vtkCell.h:50
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:23
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:29
cell represents a 1D line
Definition vtkLine.h:23
represent and manipulate point attribute data
quickly locate points in 3-space
represent and manipulate 3D points
Definition vtkPoints.h:29
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:80
a cell that represents an n-sided polygon
Definition vtkPolygon.h:32
vtkPolyhedron utilities
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
Satisfy the vtkCell API.
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
Set the faces of the polyhedron.
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
int TriangulateFaces(vtkIdList *newFaces)
Triangulate each face of the polyhedron.
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.
std::map< vtkIdType, vtkIdType > vtkPointIdMap
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
Get the faces of the polyhedron.
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
The Initialize method builds up internal structures of vtkPolyhedron.
a cell that represents a 2D quadrilateral
Definition vtkQuad.h:28
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:34
a cell that represents a triangle
Definition vtkTriangle.h:28
@ VTK_POLYHEDRON
Definition vtkCellType.h:80
int vtkIdType
Definition vtkType.h:315