VTK  9.5.20250617
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
142#ifndef vtkPolyhedron_h
143#define vtkPolyhedron_h
144
145#include "vtkCell3D.h"
146#include "vtkCellStatus.h" // For enum.
147#include "vtkCommonDataModelModule.h" // For export macro
148#include "vtkDeprecation.h" // VTK_DEPRECATED_IN_9_6_0()
149#include "vtkNew.h" // For vtkNew
150
151VTK_ABI_NAMESPACE_BEGIN
152class vtkIdTypeArray;
153class vtkCellArray;
154class vtkTriangle;
155class vtkQuad;
156class vtkTetra;
157class vtkPolygon;
158class vtkLine;
159class vtkEdgeTable;
160class vtkPolyData;
161class vtkCellLocator;
162class vtkGenericCell;
163class vtkPointLocator;
165
166class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
167{
168public:
169 using vtkPointIdMap = std::map<vtkIdType, vtkIdType>;
170
173
175
179 vtkTypeMacro(vtkPolyhedron, vtkCell3D);
180 void PrintSelf(ostream& os, vtkIndent indent) override;
182
184
188 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
189 {
190 vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
191 }
192 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
193 {
194 vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
195 return 0;
196 }
198 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
199 {
200 vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
201 }
203 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
204 {
205 vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
206 return 0;
207 }
209 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
210 {
211 vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
212 return 0;
213 }
214 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
216 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
217 {
218 vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
219 return 0;
220 }
221 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
222 {
223 vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
224 return false;
225 }
227
231 double* GetParametricCoords() override;
232
236 int GetCellType() override { return VTK_POLYHEDRON; }
237
241 int RequiresInitialization() override { return 1; }
242
248 void Initialize() override;
249
251
255 int GetNumberOfEdges() override;
256 vtkCell* GetEdge(int) override;
257 int GetNumberOfFaces() override;
258 vtkCell* GetFace(int faceId) override;
260
266 void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
267 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
268 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
269
279 void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
280 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
281 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
282
290 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
291 double& dist2, double weights[]) override;
292
297 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
298
305 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
306 double pcoords[3], int& subId) override;
307
323 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
324
332
340
349 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
350
355 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
356
361 int GetParametricCenter(double pcoords[3]) override;
362
366 int IsPrimaryCell() VTK_FUTURE_CONST override { return 1; }
367
369
374 void InterpolateFunctions(const double x[3], double* sf) override;
375 void InterpolateDerivs(const double x[3], double* derivs) override;
377
383 int RequiresExplicitFaceRepresentation() VTK_FUTURE_CONST override { return 1; }
384
402 VTK_DEPRECATED_IN_9_6_0("This function is deprecated, use SetCellFaces")
403 void SetFaces(vtkIdType* faces);
404
421 VTK_DEPRECATED_IN_9_6_0("This function is deprecated, use GetCellFaces")
422 vtkIdType* GetFaces();
423
431 int SetCellFaces(vtkCellArray* faces);
432
434
441 vtkCellArray* GetCellFaces();
442 void GetCellFaces(vtkCellArray* faces);
444
451 int IsInside(const double x[3], double tolerance);
452
469 bool IsConvex();
470 Status IsConvex(double planarThreshold);
471
475 vtkPolyData* GetPolyData();
476
480 void ShallowCopy(vtkCell* c) override;
481
485 void DeepCopy(vtkCell* c) override;
486
487protected:
489 ~vtkPolyhedron() override;
490
491 // Internal classes for supporting operations on this cell
497
498 // Filled with the SetFaces method.
499 // These faces are numbered in global id space
500 vtkNew<vtkCellArray> GlobalFaces;
501
502 // VTK_DEPRECATED_IN_9_6_0()
503 // Backward compatibility
504 vtkNew<vtkIdTypeArray> LegacyGlobalFaces;
505
506 // If edges are needed. Note that the edge numbering is in canonical space.
507 int EdgesGenerated = 0; // true/false
508 vtkNew<vtkEdgeTable> EdgeTable; // keep track of all edges
509 vtkNew<vtkIdTypeArray> Edges; // edge pairs kept in this list, in canonical id space
510 vtkNew<vtkIdTypeArray> EdgeFaces; // face pairs that comprise each edge, with the
511 // same ordering as EdgeTable
512 int GenerateEdges(); // method populates the edge table and edge array
513
514 // Numerous methods needs faces to be numbered in the canonical space.
515 // This method uses PointIdMap to fill the Faces member (faces described
516 // with canonical IDs) from the GlobalFaces member (faces described with
517 // global IDs).
518 void GenerateFaces();
519 vtkNew<vtkCellArray> Faces; // These are numbered in canonical id space
520 int FacesGenerated = 0; // True when Faces have been successfully constructed
521
522 // Bounds management
523 int BoundsComputed = 0;
524 void ComputeBounds();
525 void ComputeParametricCoordinate(const double x[3], double pc[3]);
526 void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
527
528 // Members for supporting geometric operations
529 int PolyDataConstructed = 0;
531 void ConstructPolyData();
532 int LocatorConstructed = 0;
533 vtkNew<vtkCellLocator> CellLocator;
534 void ConstructLocator();
537
538private:
539 vtkPolyhedron(const vtkPolyhedron&) = delete;
540 void operator=(const vtkPolyhedron&) = delete;
541
543
544 // vtkCell has the data members Points (x,y,z coordinates) and PointIds (global cell ids).
545 // These data members are implicitly organized in canonical space, i.e., where the cell
546 // point ids are (0,1,...,npts-1).
547 // The PointIdMap is constructed during the call of the Initialize() method and maps global
548 // point ids to the canonical point ids.
549 vtkPointIdMap PointIdMap;
550
551 void GeneratePointToIncidentFaces();
552
553 // Members used in GetPointToIncidentFaces
554 std::vector<std::vector<vtkIdType>> PointToIncidentFaces;
555
557 std::atomic<bool> IsRandomSequenceSeedInitialized{ false };
558};
559
560VTK_ABI_NAMESPACE_END
561#endif
abstract class to specify 3D cell interface
Definition vtkCell3D.h:28
object to represent cell connectivity
represent and manipulate cell attribute data
octree-based spatial search object to quickly locate cells
abstract class to specify cell behavior
Definition vtkCell.h:129
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:133
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
Park and Miller Sequence of pseudo random numbers.
Allocate and hold a VTK object.
Definition vtkNew.h:167
represent and manipulate point attribute data
quickly locate points in 3-space
concrete dataset represents vertices, lines, polygons, and triangle strips
a cell that represents an n-sided polygon
Definition vtkPolygon.h:132
vtkPolyhedron utilities
A 3D cell defined by a set of polygonal faces.
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.
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
int GetNumberOfFaces() override
A polyhedron is represented internally by a set of polygonal 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 ...
std::map< vtkIdType, vtkIdType > vtkPointIdMap
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh.
int TriangulateFaces(vtkIdList *newFaces)
Triangulate each face of the polyhedron.
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
See vtkCell3D API for description of these methods.
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).
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.
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.
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
void InterpolateDerivs(const double x[3], double *derivs) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
int IsPrimaryCell() VTK_FUTURE_CONST override
A polyhedron is a full-fledged primary cell.
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.
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.
vtkCell * GetFace(int faceId) override
A polyhedron is represented internally by a set of polygonal faces.
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.
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
static vtkPolyhedron * New()
Standard new methods.
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...
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
vtkCellStatus
Diagnostic values indicating how well-specified a cell is.
@ VTK_POLYHEDRON
Definition vtkCellType.h:80
#define VTK_DEPRECATED_IN_9_6_0(reason)
int vtkIdType
Definition vtkType.h:332