VTK  9.6.20260517
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
146
147#ifndef vtkPolyhedron_h
148#define vtkPolyhedron_h
149
150#include "vtkCell3D.h"
151#include "vtkCellStatus.h" // For enum.
152#include "vtkCommonDataModelModule.h" // For export macro
153#include "vtkDeprecation.h" // VTK_DEPRECATED_IN_9_6_0()
154#include "vtkNew.h" // For vtkNew
155
156VTK_ABI_NAMESPACE_BEGIN
157class vtkIdTypeArray;
158class vtkCellArray;
159class vtkTriangle;
160class vtkQuad;
161class vtkTetra;
162class vtkPolygon;
163class vtkLine;
164class vtkEdgeTable;
165class vtkPolyData;
166class vtkCellLocator;
167class vtkGenericCell;
170
171class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
172{
173public:
174 using vtkPointIdMap = std::map<vtkIdType, vtkIdType>;
175
178
180
184 vtkTypeMacro(vtkPolyhedron, vtkCell3D);
185 void PrintSelf(ostream& os, vtkIndent indent) override;
187
189
193 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
194 {
195 vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
196 }
197 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
199 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
200 {
201 vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
202 }
204 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
205 {
206 vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
207 return 0;
208 }
210 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
211 {
212 vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
213 return 0;
214 }
215 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
217 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
218 {
219 vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
220 return 0;
221 }
222 bool GetCentroid(double centroid[3]) const override;
223
231
235 int GetCellType() override { return VTK_POLYHEDRON; }
236
240 int RequiresInitialization() override { return 1; }
241
247 void Initialize() override;
248
250
254 int GetNumberOfEdges() override;
255 vtkCell* GetEdge(int) override;
256 int GetNumberOfFaces() override;
257 vtkCell* GetFace(int faceId) override;
259
265 void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
266 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
267 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
268
278 void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
279 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
280 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
281
289 void ClipWithContext(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
290 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
291 vtkIdType cellId, vtkCellData* outCd, int insideOut, vtkCellArray* outFaces,
292 vtkCellArray* outFaceLocs);
293
298 void ClipWithContext(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
299 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
300 vtkIdType cellId, vtkCellData* outCd, int insideOut, vtkUnstructuredGrid* outUG);
301
309 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
310 double& dist2, double weights[]) override;
311
316 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
317
324 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
325 double pcoords[3], int& subId) override;
326
342 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
343
351
359
368 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
369
374 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
375
380 int GetParametricCenter(double pcoords[3]) override;
381
385 int IsPrimaryCell() VTK_FUTURE_CONST override { return 1; }
386
388
393 void InterpolateFunctions(const double x[3], double* sf) override;
394 void InterpolateDerivs(const double x[3], double* derivs) override;
396
402 int RequiresExplicitFaceRepresentation() VTK_FUTURE_CONST override { return 1; }
403
421 VTK_DEPRECATED_IN_9_6_0("This function is deprecated, use SetCellFaces")
422 void SetFaces(vtkIdType* faces);
423
440 VTK_DEPRECATED_IN_9_6_0("This function is deprecated, use GetCellFaces")
442
451
453
462
469 int IsInside(const double x[3], double tolerance);
470
487 bool IsConvex();
488 Status IsConvex(double planarThreshold);
489
494
498 void ShallowCopy(vtkCell* c) override;
499
503 void DeepCopy(vtkCell* c) override;
504
505protected:
507 ~vtkPolyhedron() override;
508
509 // Internal classes for supporting operations on this cell
515
516 // Filled with the SetFaces method.
517 // These faces are numbered in global id space
519
520 // VTK_DEPRECATED_IN_9_6_0()
521 // Backward compatibility
523
524 // If edges are needed. Note that the edge numbering is in canonical space.
525 int EdgesGenerated = 0; // true/false
526 vtkNew<vtkEdgeTable> EdgeTable; // keep track of all edges
527 vtkNew<vtkIdTypeArray> Edges; // edge pairs kept in this list, in canonical id space
528 vtkNew<vtkIdTypeArray> EdgeFaces; // face pairs that comprise each edge, with the
529 // same ordering as EdgeTable
530 int GenerateEdges(); // method populates the edge table and edge array
531
532 // Numerous methods needs faces to be numbered in the canonical space.
533 // This method uses PointIdMap to fill the Faces member (faces described
534 // with canonical IDs) from the GlobalFaces member (faces described with
535 // global IDs).
537 vtkNew<vtkCellArray> Faces; // These are numbered in canonical id space
538 int FacesGenerated = 0; // True when Faces have been successfully constructed
539
540 // Bounds management
543 void ComputeParametricCoordinate(const double x[3], double pc[3]);
544 void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
545
546 // Members for supporting geometric operations
555
556private:
557 vtkPolyhedron(const vtkPolyhedron&) = delete;
558 void operator=(const vtkPolyhedron&) = delete;
559
562
563 // vtkCell has the data members Points (x,y,z coordinates) and PointIds (global cell ids).
564 // These data members are implicitly organized in canonical space, i.e., where the cell
565 // point ids are (0,1,...,npts-1).
566 // The PointIdMap is constructed during the call of the Initialize() method and maps global
567 // point ids to the canonical point ids.
568 vtkPointIdMap PointIdMap;
569
570 void GeneratePointToIncidentFaces();
571
572 // This variant of GenerateEdges() always regenerates edges but also populates
573 // \a unevenCoedges with a map from edge ID to the signed number of mismatched coedges.
574 int GenerateEdges(std::map<vtkIdType, int>& unevenCoedges);
575
576 // Members used in GetPointToIncidentFaces
577 std::vector<std::vector<vtkIdType>> PointToIncidentFaces;
578};
579
580VTK_ABI_NAMESPACE_END
581#endif
object to represent cell connectivity
represent and manipulate cell attribute data
octree-based spatial search object to quickly locate cells
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:135
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:108
cell represents a 1D line
Definition vtkLine.h:132
Allocate and hold a VTK object.
Definition vtkNew.h:168
represent and manipulate point attribute data
concrete dataset represents vertices, lines, polygons, and triangle strips
a cell that represents an n-sided polygon
Definition vtkPolygon.h:134
Isosurface extraction using López polygon tracing algorithm.
vtkNew< vtkLine > Line
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
void GetEdgePoints(vtkIdType edgeId, const vtkIdType *&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.
vtkCellArray * GetCellFaces()
Get the faces of the polyhedron.
vtkNew< vtkCellArray > GlobalFaces
int SetCellFaces(vtkCellArray *faces)
Set the faces of the polyhedron.
void SetFaces(vtkIdType *faces)
Set the faces of the polyhedron.
void ShallowCopy(vtkCell *c) override
Shallow copy of a polyhedron.
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
int GenerateEdges()
int GetNumberOfFaces() override
A polyhedron is represented internally by a set of polygonal faces.
vtkNew< vtkTetra > Tetra
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 ...
vtkCellStatus Status
Adopt vtkCellStatus to describe degenerate polyhedral cells.
vtkNew< vtkTriangle > Triangle
std::map< vtkIdType, vtkIdType > vtkPointIdMap
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh.
void GenerateFaces()
friend class vtkPolyhedronContour
void ComputeParametricCoordinate(const double x[3], double pc[3])
vtkNew< vtkPolyData > PolyData
int TriangulateFaces(vtkIdList *newFaces)
Triangulate each face of the polyhedron.
void ComputeBounds()
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
See vtkCell3D API for description of these methods.
vtkNew< vtkQuad > Quad
vtkNew< vtkGenericCell > Cell
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
The inverse of EvaluatePosition.
vtkNew< vtkPolygon > Polygon
void InterpolateFunctions(const double x[3], double *sf) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
int RequiresExplicitFaceRepresentation() VTK_FUTURE_CONST override
Satisfy the vtkCell API.
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.
vtkIdType * GetFaces()
Get the faces of the polyhedron.
int IsInside(const double x[3], double tolerance)
A method particular to vtkPolyhedron.
vtkNew< vtkIdTypeArray > LegacyGlobalFaces
bool GetCentroid(double centroid[3]) const override
See vtkCell3D API for description of these methods.
int GetParametricCenter(double pcoords[3]) override
Return the center of the cell in parametric coordinates.
int TriangulateFaces(vtkCellArray *newFaces)
Triangulate each face of the polyhedron.
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.
int IsPrimaryCell() VTK_FUTURE_CONST override
A polyhedron is a full-fledged primary cell.
vtkNew< vtkIdTypeArray > EdgeFaces
void ConstructPolyData()
void ClipWithContext(double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut, vtkCellArray *outFaces, vtkCellArray *outFaceLocs)
Clip this polyhedron and write faces directly into outFaces and outFaceLocs, bypassing the embedded f...
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 GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
vtkPolyData * GetPolyData()
Construct polydata if no one exist, then return this->PolyData.
vtkNew< vtkCellArray > Faces
vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
vtkCell * GetEdge(int) override
A polyhedron is represented internally by a set of polygonal faces.
void ClipWithContext(double value, vtkDataArray *scalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut, vtkUnstructuredGrid *outUG)
Convenience overload of ClipWithContext that extracts outFaces and outFaceLocs from the given output ...
vtkCell * GetFace(int faceId) override
A polyhedron is represented internally by a set of polygonal faces.
void DeepCopy(vtkCell *c) override
Deep copy of a polyhedron.
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.
vtkNew< vtkIdList > CellIds
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
friend class vtkPolyhedronUtilities
void ComputePositionFromParametricCoordinate(const double pc[3], double x[3])
static vtkPolyhedron * New()
Standard new methods.
vtkNew< vtkCellLocator > CellLocator
double ComputeVolume()
Compute the volume of the polyhedron using the divergence theorem.
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...
vtkNew< vtkEdgeTable > EdgeTable
vtkNew< vtkIdTypeArray > Edges
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:87
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:113
a cell that represents a triangle
dataset represents arbitrary combinations of all possible cell types
vtkCellStatus
Diagnostic values indicating how well-specified a cell is.
@ VTK_POLYHEDRON
Definition vtkCellType.h:80
#define vtkDataArray
#define VTK_DEPRECATED_IN_9_6_0(reason)
int vtkIdType
Definition vtkType.h:363