VTK  9.2.5
vtkPolyhedron.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPolyhedron.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
77 #ifndef vtkPolyhedron_h
78 #define vtkPolyhedron_h
79 
80 #include "vtkCell3D.h"
81 #include "vtkCommonDataModelModule.h" // For export macro
82 
83 class vtkIdTypeArray;
84 class vtkCellArray;
85 class vtkTriangle;
86 class vtkQuad;
87 class vtkTetra;
88 class vtkPolygon;
89 class vtkLine;
90 class vtkPointIdMap;
91 class vtkIdToIdVectorMapType;
92 class vtkIdToIdMapType;
93 class vtkEdgeTable;
94 class vtkPolyData;
95 class vtkCellLocator;
96 class vtkGenericCell;
97 class vtkPointLocator;
98 
99 class VTKCOMMONDATAMODEL_EXPORT vtkPolyhedron : public vtkCell3D
100 {
101 public:
103 
106  static vtkPolyhedron* New();
107  vtkTypeMacro(vtkPolyhedron, vtkCell3D);
108  void PrintSelf(ostream& os, vtkIndent indent) override;
110 
112 
116  void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
117  {
118  vtkWarningMacro(<< "vtkPolyhedron::GetEdgePoints Not Implemented");
119  }
120  vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
121  {
122  vtkWarningMacro(<< "vtkPolyhedron::GetFacePoints Not Implemented");
123  return 0;
124  }
126  vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
127  {
128  vtkWarningMacro(<< "vtkPolyhedron::GetEdgeToAdjacentFaces Not Implemented");
129  }
131  vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
132  {
133  vtkWarningMacro(<< "vtkPolyhedron::GetFaceToAdjacentFaces Not Implemented");
134  return 0;
135  }
137  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
138  {
139  vtkWarningMacro(<< "vtkPolyhedron::GetPointToIncidentEdges Not Implemented");
140  return 0;
141  }
142  vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
144  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
145  {
146  vtkWarningMacro(<< "vtkPolyhedron::GetPointToOneRingPoints Not Implemented");
147  return 0;
148  }
149  bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
150  {
151  vtkWarningMacro(<< "vtkPolyhedron::GetCentroid Not Implemented");
152  return false;
153  }
155 
159  double* GetParametricCoords() override;
160 
164  int GetCellType() override { return VTK_POLYHEDRON; }
165 
169  int RequiresInitialization() override { return 1; }
170  void Initialize() override;
171 
173 
177  int GetNumberOfEdges() override;
178  vtkCell* GetEdge(int) override;
179  int GetNumberOfFaces() override;
180  vtkCell* GetFace(int faceId) override;
182 
188  void Contour(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
189  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
190  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
191 
201  void Clip(double value, vtkDataArray* scalars, vtkIncrementalPointLocator* locator,
202  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
203  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
204 
212  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
213  double& dist2, double weights[]) override;
214 
219  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
220 
227  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
228  double pcoords[3], int& subId) override;
229 
245  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
246 
255  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
256 
261  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
262 
267  int GetParametricCenter(double pcoords[3]) override;
268 
272  int IsPrimaryCell() override { return 1; }
273 
275 
280  void InterpolateFunctions(const double x[3], double* sf) override;
281  void InterpolateDerivs(const double x[3], double* derivs) override;
283 
285 
293  int RequiresExplicitFaceRepresentation() override { return 1; }
294  void SetFaces(vtkIdType* faces) override;
295  vtkIdType* GetFaces() override;
297 
304  int IsInside(const double x[3], double tolerance);
305 
312  bool IsConvex();
313 
318 
319 protected:
321  ~vtkPolyhedron() override;
322 
323  // Internal classes for supporting operations on this cell
329  vtkIdTypeArray* GlobalFaces; // these are numbered in global id space
331 
332  // vtkCell has the data members Points (x,y,z coordinates) and PointIds
333  // (global cell ids corresponding to cell canonical numbering (0,1,2,....)).
334  // These data members are implicitly organized in canonical space, i.e., where
335  // the cell point ids are (0,1,...,npts-1). The PointIdMap maps global point id
336  // back to these canonoical point ids.
337  vtkPointIdMap* PointIdMap;
338 
339  // If edges are needed. Note that the edge numbering is in
340  // canonical space.
341  int EdgesGenerated; // true/false
342  vtkEdgeTable* EdgeTable; // keep track of all edges
343  vtkIdTypeArray* Edges; // edge pairs kept in this list, in canonical id space
344  vtkIdTypeArray* EdgeFaces; // face pairs that comprise each edge, with the
345  // same ordering as EdgeTable
346  int GenerateEdges(); // method populates the edge table and edge array
347 
348  // If faces need renumbering into canonical numbering space these members
349  // are used. When initiallly loaded, the face numbering uses global dataset
350  // ids. Once renumbered, they are converted to canonical space.
351  vtkIdTypeArray* Faces; // these are numbered in canonical id space
354 
355  // Bounds management
358  void ComputeParametricCoordinate(const double x[3], double pc[3]);
359  void ComputePositionFromParametricCoordinate(const double pc[3], double x[3]);
360 
362 
363  // Members for supporting geometric operations
373 
374  // Members used in GetPointToIncidentFaces
377 
378 private:
379  vtkPolyhedron(const vtkPolyhedron&) = delete;
380  void operator=(const vtkPolyhedron&) = delete;
381 };
382 
383 //----------------------------------------------------------------------------
384 inline int vtkPolyhedron::GetParametricCenter(double pcoords[3])
385 {
386  pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
387  return 0;
388 }
389 
390 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
object to represent cell connectivity
Definition: vtkCellArray.h:296
represent and manipulate cell attribute data
Definition: vtkCellData.h:151
octree-based spatial search object to quickly locate cells
abstract class to specify cell behavior
Definition: vtkCell.h:150
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
keep track of edges (edge is pair of integer id's)
Definition: vtkEdgeTable.h:41
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:143
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:119
cell represents a 1D line
Definition: vtkLine.h:143
represent and manipulate point attribute data
Definition: vtkPointData.h:151
quickly locate points in 3-space
represent and manipulate 3D points
Definition: vtkPoints.h:149
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:200
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:152
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
Methods supporting the definition of faces.
vtkIdList * CellIds
vtkPolyData * GetPolyData()
Construct polydata if no one exist, then return this->PolyData.
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
Methods supporting the definition of 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 ...
~vtkPolyhedron() override
void GenerateFaces()
vtkIdType GetPointToIncidentEdges(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(edgeIds)) override
See vtkCell3D API for description of these methods.
vtkQuad * Quad
void ComputeParametricCoordinate(const double x[3], double pc[3])
vtkIdType * ValenceAtPoint
vtkEdgeTable * EdgeTable
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.
vtkCell * GetFace(int faceId) override
A polyhedron is represented internally by a set of polygonal faces.
static vtkPolyhedron * New()
Standard new methods.
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.
vtkIdTypeArray * GlobalFaces
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.
vtkCell * GetEdge(int) override
A polyhedron is represented internally by a set of polygonal faces.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard new methods.
void GeneratePointToIncidentFacesAndValenceAtPoint()
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
vtkIdTypeArray * Faces
vtkPolyData * PolyData
int IsPrimaryCell() override
A polyhedron is a full-fledged primary cell.
vtkLine * Line
vtkTetra * Tetra
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.
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
vtkGenericCell * Cell
void ComputePositionFromParametricCoordinate(const double pc[3], double x[3])
vtkCellArray * Polys
vtkIdType * GetFaces() override
Methods supporting the definition of faces.
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
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:98
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:114
a cell that represents a triangle
Definition: vtkTriangle.h:148
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_POLYHEDRON
Definition: vtkCellType.h:128
int vtkIdType
Definition: vtkType.h:332