VTK  9.5.2
vtkConvexPointSet.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
20
21#ifndef vtkConvexPointSet_h
22#define vtkConvexPointSet_h
23
24#include "vtkCell3D.h"
25#include "vtkCommonDataModelModule.h" // For export macro
26#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_5_0
27
28VTK_ABI_NAMESPACE_BEGIN
30class vtkCellArray;
31class vtkTriangle;
32class vtkTetra;
33class vtkDoubleArray;
34
35class VTKCOMMONDATAMODEL_EXPORT vtkConvexPointSet : public vtkCell3D
36{
37public:
40 void PrintSelf(ostream& os, vtkIndent indent) override;
41
45#ifndef VTK_LEGACY_REMOVE
46 VTK_DEPRECATED_IN_9_5_0("HasFixedTopology() is always 0 and will be removed")
47 virtual vtkTypeBool HasFixedTopology() { return 0; }
48#endif
49
51
55 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
56 {
57 vtkWarningMacro(<< "vtkConvexPointSet::GetEdgePoints Not Implemented");
58 }
59 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
60 {
61 vtkWarningMacro(<< "vtkConvexPointSet::GetFacePoints Not Implemented");
62 return 0;
63 }
65 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
66 {
67 vtkWarningMacro(<< "vtkConvexPointSet::GetEdgeToAdjacentFaces Not Implemented");
68 }
70 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
71 {
72 vtkWarningMacro(<< "vtkConvexPointSet::GetFaceToAdjacentFaces Not Implemented");
73 return 0;
74 }
76 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
77 {
78 vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentEdges Not Implemented");
79 return 0;
80 }
82 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(faceIds)) override
83 {
84 vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentFaces Not Implemented");
85 return 0;
86 }
88 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
89 {
90 vtkWarningMacro(<< "vtkConvexPointSet::GetPointToOneRingPoints Not Implemented");
91 return 0;
92 }
93 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
94 {
95 vtkWarningMacro(<< "vtkConvexPointSet::GetCentroid Not Implemented");
96 return false;
97 }
98
99
103 double* GetParametricCoords() override;
104
108 int GetCellType() override { return VTK_CONVEX_POINT_SET; }
109
113 int RequiresInitialization() override { return 1; }
114 void Initialize() override;
115
117
127 int GetNumberOfEdges() override { return 0; }
128 vtkCell* GetEdge(int) override { return nullptr; }
129 int GetNumberOfFaces() override;
130 vtkCell* GetFace(int faceId) override;
132
137 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
138 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
139 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
140
146 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
147 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
148 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
149
155 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
156 double& dist2, double weights[]) override;
157
161 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
162
167 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
168 double pcoords[3], int& subId) override;
169
173 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
174
180 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
181
187 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
188
192 int GetParametricCenter(double pcoords[3]) override;
193
198 int IsPrimaryCell() VTK_FUTURE_CONST override { return 0; }
199
201
205 void InterpolateFunctions(const double pcoords[3], double* sf) override;
206 void InterpolateDerivs(const double pcoords[3], double* derivs) override;
208
209protected:
212
217
221
222private:
223 vtkConvexPointSet(const vtkConvexPointSet&) = delete;
224 void operator=(const vtkConvexPointSet&) = delete;
225};
226
227//----------------------------------------------------------------------------
228inline int vtkConvexPointSet::GetParametricCenter(double pcoords[3])
229{
230 pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
231 return 0;
232}
233
234VTK_ABI_NAMESPACE_END
235#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:32
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
vtkDoubleArray * ParametricCoords
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkCell * GetFace(int faceId) override
A convex point set has no explicit cell edge or faces; however implicitly (after triangulation) it do...
bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
See vtkCell3D API for description of these methods.
int GetNumberOfEdges() override
A convex point set has no explicit cell edge or faces; however implicitly (after triangulation) it do...
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
vtkDoubleArray * TetraScalars
static vtkConvexPointSet * New()
void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
Satisfy the vtkCell API.
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Returns the set of points forming a face of the triangulation of these points that are on the boundar...
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Satisfy the vtkCell API.
vtkIdType GetPointToIncidentEdges(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(edgeIds)) override
See vtkCell3D API for description of these methods.
vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
~vtkConvexPointSet() override
void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
void Initialize() override
int IsPrimaryCell() VTK_FUTURE_CONST override
A convex point set is triangulated prior to any operations on it so it is not a primary cell,...
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Triangulates the cells and then intersects them to determine the intersection point.
void InterpolateDerivs(const double pcoords[3], double *derivs) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
vtkCellArray * BoundaryTris
virtual vtkTypeBool HasFixedTopology()
See vtkCell3D API for description of this method.
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Triangulate using methods of vtkOrderedTriangulator.
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Satisfy the vtkCell API.
int GetParametricCenter(double pcoords[3]) override
Return the center of the cell in parametric coordinates.
vtkIdType GetFaceToAdjacentFaces(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(faceIds)) override
See vtkCell3D API for description of these methods.
void GetEdgeToAdjacentFaces(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
int GetCellType() override
See the vtkCell API for descriptions of these methods.
int GetNumberOfFaces() override
A convex point set has no explicit cell edge or faces; however implicitly (after triangulation) it do...
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Computes derivatives by triangulating and from subId and pcoords, evaluating derivatives on the resul...
void InterpolateFunctions(const double pcoords[3], double *sf) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
vtkCell * GetEdge(int) override
A convex point set has no explicit cell edge or faces; however implicitly (after triangulation) it do...
vtkIdType GetPointToIncidentFaces(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(faceIds)) 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.
vtkIdType GetPointToOneRingPoints(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:24
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:29
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:30
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:34
a cell that represents a triangle
Definition vtkTriangle.h:28
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:64
@ VTK_CONVEX_POINT_SET
Definition vtkCellType.h:77
#define vtkDataArray
#define VTK_DEPRECATED_IN_9_5_0(reason)
int vtkIdType
Definition vtkType.h:332