VTK  9.3.20240425
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
148VTK_ABI_NAMESPACE_BEGIN
149class vtkIdTypeArray;
150class vtkCellArray;
151class vtkTriangle;
152class vtkQuad;
153class vtkTetra;
154class vtkPolygon;
155class vtkLine;
156class vtkIdToIdVectorMapType;
157class vtkIdToIdMapType;
158class vtkEdgeTable;
159class vtkPolyData;
160class vtkCellLocator;
161class vtkGenericCell;
162class vtkPointLocator;
163
164class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
165{
166public:
167 typedef std::map<vtkIdType, vtkIdType> vtkPointIdMap;
168
170
174 vtkTypeMacro(vtkPolyhedron, vtkCell3D);
175 void PrintSelf(ostream& os, vtkIndent indent) override;
177
179
183 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
184 {
185 vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
186 }
187 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
188 {
189 vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
190 return 0;
191 }
193 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
194 {
195 vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
196 }
198 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
199 {
200 vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
201 return 0;
202 }
204 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
205 {
206 vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
207 return 0;
208 }
209 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
211 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
212 {
213 vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
214 return 0;
215 }
216 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
217 {
218 vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
219 return false;
220 }
222
226 double* GetParametricCoords() override;
227
231 int GetCellType() override { return VTK_POLYHEDRON; }
232
236 int RequiresInitialization() override { return 1; }
237
243 void Initialize() override;
244
246
250 int GetNumberOfEdges() override;
251 vtkCell* GetEdge(int) override;
252 int GetNumberOfFaces() override;
253 vtkCell* GetFace(int faceId) override;
255
261 void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
262 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
263 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
264
274 void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
275 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
276 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
277
285 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
286 double& dist2, double weights[]) override;
287
292 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
293
300 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
301 double pcoords[3], int& subId) override;
302
318 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
319
327
335
344 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
345
350 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
351
356 int GetParametricCenter(double pcoords[3]) override;
357
361 int IsPrimaryCell() override { return 1; }
362
364
369 void InterpolateFunctions(const double x[3], double* sf) override;
370 void InterpolateDerivs(const double x[3], double* derivs) override;
372
378 int RequiresExplicitFaceRepresentation() override { return 1; }
379
397 void SetFaces(vtkIdType* faces) override;
398
415 vtkIdType* GetFaces() override;
416
425
427
437
444 int IsInside(const double x[3], double tolerance);
445
452 bool IsConvex();
453
458
459protected:
461 ~vtkPolyhedron() override;
462
463 // Internal classes for supporting operations on this cell
469
470 // Filled with the SetFaces method.
471 // These faces are numbered in global id space
473
474 // Backward compatibility
476
477 // vtkCell has the data members Points (x,y,z coordinates) and PointIds (global cell ids).
478 // These data members are implicitly organized in canonical space, i.e., where the cell
479 // point ids are (0,1,...,npts-1).
480 // The PointIdMap is constructed during the call of the Initialize() method and maps global
481 // point ids to the canonical point ids.
483
484 // If edges are needed. Note that the edge numbering is in canonical space.
485 int EdgesGenerated; // true/false
486 vtkEdgeTable* EdgeTable; // keep track of all edges
487 vtkIdTypeArray* Edges; // edge pairs kept in this list, in canonical id space
488 vtkIdTypeArray* EdgeFaces; // face pairs that comprise each edge, with the
489 // same ordering as EdgeTable
490 int GenerateEdges(); // method populates the edge table and edge array
491
492 // Numerous methods needs faces to be numbered in the canonical space.
493 // This method uses PointIdMap to fill the Faces member (faces described
494 // with canonical IDs) from the GlobalFaces member (faces described with
495 // global IDs).
497 vtkCellArray* Faces; // These are numbered in canonical id space
498 int FacesGenerated; // True when Faces have been successfully constructed
499
500 // Bounds management
503 void ComputeParametricCoordinate(const double x[3], double pc[3]);
504 void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
505
507
508 // Members for supporting geometric operations
517
518 // Members used in GetPointToIncidentFaces
521
522private:
523 vtkPolyhedron(const vtkPolyhedron&) = delete;
524 void operator=(const vtkPolyhedron&) = delete;
525
527};
528
529//----------------------------------------------------------------------------
530inline int vtkPolyhedron::GetParametricCenter(double pcoords[3])
531{
532 pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
533 return 0;
534}
535
536VTK_ABI_NAMESPACE_END
537#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
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: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
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 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
vtkCellArray * GetCellFaces()
Get the faces of the polyhedron.
int SetCellFaces(vtkCellArray *faces)
Set the faces of the polyhedron.
vtkCellArray * GlobalFaces
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 ...
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])
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.
vtkCellArray * Faces
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
int TriangulateFaces(vtkCellArray *newFaces)
Triangulate each face of the polyhedron.
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
vtkIdTypeArray * LegacyGlobalFaces
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()
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.
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
@ VTK_POLYHEDRON
Definition vtkCellType.h:99
int vtkIdType
Definition vtkType.h:315