VTK  9.4.20241221
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 "vtkCommonDataModelModule.h" // For export macro
147#include "vtkDeprecation.h" // For VTK_DEPRECATED
148#include "vtkNew.h" // For vtkNew
149
150VTK_ABI_NAMESPACE_BEGIN
151class vtkIdTypeArray;
152class vtkCellArray;
153class vtkTriangle;
154class vtkQuad;
155class vtkTetra;
156class vtkPolygon;
157class vtkLine;
158class vtkEdgeTable;
159class vtkPolyData;
160class vtkCellLocator;
161class vtkGenericCell;
162class vtkPointLocator;
164
165class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
166{
167public:
168 using vtkPointIdMap = std::map<vtkIdType, vtkIdType>;
169
171
175 vtkTypeMacro(vtkPolyhedron, vtkCell3D);
176 void PrintSelf(ostream& os, vtkIndent indent) override;
178
180
184 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
185 {
186 vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
187 }
188 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
189 {
190 vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
191 return 0;
192 }
194 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
195 {
196 vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
197 }
199 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
200 {
201 vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
202 return 0;
203 }
205 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
206 {
207 vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
208 return 0;
209 }
210 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
212 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
213 {
214 vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
215 return 0;
216 }
217 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
218 {
219 vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
220 return false;
221 }
223
227 double* GetParametricCoords() override;
228
232 int GetCellType() override { return VTK_POLYHEDRON; }
233
237 int RequiresInitialization() override { return 1; }
238
244 void Initialize() override;
245
247
251 int GetNumberOfEdges() override;
252 vtkCell* GetEdge(int) override;
253 int GetNumberOfFaces() override;
254 vtkCell* GetFace(int faceId) override;
256
262 void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
263 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
264 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
265
275 void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
276 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
277 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
278
286 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
287 double& dist2, double weights[]) override;
288
293 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
294
301 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
302 double pcoords[3], int& subId) override;
303
319 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
320
328
336
345 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
346
351 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
352
357 int GetParametricCenter(double pcoords[3]) override;
358
362 int IsPrimaryCell() override { return 1; }
363
365
370 void InterpolateFunctions(const double x[3], double* sf) override;
371 void InterpolateDerivs(const double x[3], double* derivs) override;
373
379 int RequiresExplicitFaceRepresentation() override { return 1; }
380
398 void SetFaces(vtkIdType* faces) override;
399
416 vtkIdType* GetFaces() override;
417
426
428
438
445 int IsInside(const double x[3], double tolerance);
446
453 bool IsConvex();
454
459
463 void ShallowCopy(vtkCell* c) override;
464
468 void DeepCopy(vtkCell* c) override;
469
470protected:
472 ~vtkPolyhedron() override;
473
474 // Internal classes for supporting operations on this cell
480
481 // Filled with the SetFaces method.
482 // These faces are numbered in global id space
484
485 // Backward compatibility
487
488 // If edges are needed. Note that the edge numbering is in canonical space.
489 int EdgesGenerated = 0; // true/false
490 vtkNew<vtkEdgeTable> EdgeTable; // keep track of all edges
491 vtkNew<vtkIdTypeArray> Edges; // edge pairs kept in this list, in canonical id space
492 vtkNew<vtkIdTypeArray> EdgeFaces; // face pairs that comprise each edge, with the
493 // same ordering as EdgeTable
494 int GenerateEdges(); // method populates the edge table and edge array
495
496 // Numerous methods needs faces to be numbered in the canonical space.
497 // This method uses PointIdMap to fill the Faces member (faces described
498 // with canonical IDs) from the GlobalFaces member (faces described with
499 // global IDs).
501 vtkNew<vtkCellArray> Faces; // These are numbered in canonical id space
502 int FacesGenerated = 0; // True when Faces have been successfully constructed
503
504 // Bounds management
505 int BoundsComputed = 0;
507 void ComputeParametricCoordinate(const double x[3], double pc[3]);
508 void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
509
510 VTK_DEPRECATED_IN_9_4_0("Use GeneratePointToIncidentFaces instead.")
511 void GeneratePointToIncidentFacesAndValenceAtPoint() { this->GeneratePointToIncidentFaces(); }
512
513 // Members for supporting geometric operations
514 int PolyDataConstructed = 0;
517 int LocatorConstructed = 0;
522
523private:
524 vtkPolyhedron(const vtkPolyhedron&) = delete;
525 void operator=(const vtkPolyhedron&) = delete;
526
528
529 // vtkCell has the data members Points (x,y,z coordinates) and PointIds (global cell ids).
530 // These data members are implicitly organized in canonical space, i.e., where the cell
531 // point ids are (0,1,...,npts-1).
532 // The PointIdMap is constructed during the call of the Initialize() method and maps global
533 // point ids to the canonical point ids.
534 vtkPointIdMap PointIdMap;
535
536 void GeneratePointToIncidentFaces();
537
538 // Members used in GetPointToIncidentFaces
539 std::vector<std::vector<vtkIdType>> PointToIncidentFaces;
540
542 std::atomic<bool> IsRandomSequenceSeedInitialized{ false };
543};
544
545VTK_ABI_NAMESPACE_END
546#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:130
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:131
vtkPolyhedron utilities
A 3D cell defined by a set of polygonal faces.
vtkNew< vtkLine > Line
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.
int RequiresExplicitFaceRepresentation() override
Satisfy the vtkCell API.
vtkCellArray * GetCellFaces()
Get the faces of the polyhedron.
vtkNew< vtkCellArray > GlobalFaces
int SetCellFaces(vtkCellArray *faces)
Set the faces of the polyhedron.
void ShallowCopy(vtkCell *c) override
Shallow copy of a polyhedron.
int GenerateEdges()
int GetNumberOfFaces() override
A polyhedron is represented internally by a set of polygonal faces.
vtkNew< vtkTetra > Tetra
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 ...
vtkNew< vtkTriangle > Triangle
std::map< vtkIdType, vtkIdType > vtkPointIdMap
void GetCellFaces(vtkCellArray *faces)
Get the faces of the polyhedron.
~vtkPolyhedron() override
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Use vtkOrderedTriangulator to tetrahedralize the polyhedron mesh.
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])
vtkNew< vtkPolyData > PolyData
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.
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).
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.
vtkNew< vtkIdTypeArray > LegacyGlobalFaces
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).
bool IsConvex()
Determine whether or not a polyhedron is convex.
vtkNew< vtkIdTypeArray > EdgeFaces
void ConstructPolyData()
int GetNumberOfEdges() override
A polyhedron is represented internally by a set of polygonal faces.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard new methods.
int IsPrimaryCell() override
A polyhedron is a full-fledged primary cell.
vtkPolyData * GetPolyData()
Construct polydata if no one exist, then return this->PolyData.
vtkNew< vtkCellArray > Faces
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 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.
bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
See vtkCell3D API for description of these methods.
vtkNew< vtkIdList > CellIds
vtkIdType * GetFaces() override
Get the faces of the polyhedron.
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
void ComputePositionFromParametricCoordinate(const double pc[3], double x[3])
static vtkPolyhedron * New()
Standard new methods.
vtkNew< vtkCellLocator > CellLocator
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
@ VTK_POLYHEDRON
Definition vtkCellType.h:80
#define VTK_DEPRECATED_IN_9_4_0(reason)
int vtkIdType
Definition vtkType.h:315