VTK  9.1.0
vtkConvexPointSet.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkConvexPointSet.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=========================================================================*/
39#ifndef vtkConvexPointSet_h
40#define vtkConvexPointSet_h
41
42#include "vtkCell3D.h"
43#include "vtkCommonDataModelModule.h" // For export macro
44
46class vtkCellArray;
47class vtkTriangle;
48class vtkTetra;
49class vtkDoubleArray;
50
51class VTKCOMMONDATAMODEL_EXPORT vtkConvexPointSet : public vtkCell3D
52{
53public:
56 void PrintSelf(ostream& os, vtkIndent indent) override;
57
61 virtual int HasFixedTopology() { return 0; }
62
64
68 void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
69 {
70 vtkWarningMacro(<< "vtkConvexPointSet::GetEdgePoints Not Implemented");
71 }
72 // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
74 "Replaced by vtkConvexPointSet::GetEdgePoints(vtkIdType, const vtkIdType*&)")
75 void GetEdgePoints(int vtkNotUsed(edgeId), int*& vtkNotUsed(pts)) override
76 {
77 vtkErrorMacro(<< "vtkConvexPointSet::GetEdgePoints Not Implemented. "
78 "Also note that this signature is deprecated. "
79 "Please use GetEdgePoints(vtkIdType, const vtkIdType*& instead");
80 }
81 vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
82 {
83 vtkWarningMacro(<< "vtkConvexPointSet::GetFacePoints Not Implemented");
84 return 0;
85 }
86 // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
88 "Replaced by vtkConvexPointSet::GetFacePoints(vtkIdType, const vtkIdType*&)")
89 void GetFacePoints(int vtkNotUsed(faceId), int*& vtkNotUsed(pts)) override
90 {
91 vtkErrorMacro(<< "vtkConvexPointSet::GetFacePoints Not Implemented. "
92 "Also note that this signature is deprecated. "
93 "Please use GetFacePoints(vtkIdType, const vtkIdType*& instead");
94 }
96 vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
97 {
98 vtkWarningMacro(<< "vtkConvexPointSet::GetEdgeToAdjacentFaces Not Implemented");
99 }
101 vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
102 {
103 vtkWarningMacro(<< "vtkConvexPointSet::GetFaceToAdjacentFaces Not Implemented");
104 return 0;
105 }
107 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
108 {
109 vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentEdges Not Implemented");
110 return 0;
111 }
113 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(faceIds)) override
114 {
115 vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentFaces Not Implemented");
116 return 0;
117 }
119 vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
120 {
121 vtkWarningMacro(<< "vtkConvexPointSet::GetPointToOneRingPoints Not Implemented");
122 return 0;
123 }
124 bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
125 {
126 vtkWarningMacro(<< "vtkConvexPointSet::GetCentroid Not Implemented");
127 return false;
128 }
130
134 double* GetParametricCoords() override;
135
139 int GetCellType() override { return VTK_CONVEX_POINT_SET; }
140
144 int RequiresInitialization() override { return 1; }
145 void Initialize() override;
146
148
158 int GetNumberOfEdges() override { return 0; }
159 vtkCell* GetEdge(int) override { return nullptr; }
160 int GetNumberOfFaces() override;
161 vtkCell* GetFace(int faceId) override;
163
168 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
169 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
170 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
171
177 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
178 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
179 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
180
186 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
187 double& dist2, double weights[]) override;
188
192 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
193
198 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
199 double pcoords[3], int& subId) override;
200
204 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
205
211 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
212
218 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
219
223 int GetParametricCenter(double pcoords[3]) override;
224
229 int IsPrimaryCell() override { return 0; }
230
232
236 void InterpolateFunctions(const double pcoords[3], double* sf) override;
237 void InterpolateDerivs(const double pcoords[3], double* derivs) override;
239
240protected:
243
248
252
253private:
254 vtkConvexPointSet(const vtkConvexPointSet&) = delete;
255 void operator=(const vtkConvexPointSet&) = delete;
256};
257
258//----------------------------------------------------------------------------
259inline int vtkConvexPointSet::GetParametricCenter(double pcoords[3])
260{
261 pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
262 return 0;
263}
264
265#endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:40
object to represent cell connectivity
Definition: vtkCellArray.h:190
represent and manipulate cell attribute data
Definition: vtkCellData.h:42
abstract class to specify cell behavior
Definition: vtkCell.h:67
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
a 3D cell defined by a set of convex points
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Triangulate using methods of vtkOrderedTriangulator.
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 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.
vtkTriangle * Triangle
void InterpolateDerivs(const double pcoords[3], double *derivs) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
vtkCellArray * BoundaryTris
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 IsPrimaryCell() override
A convex point set is triangulated prior to any operations on it so it is not a primary cell,...
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.
virtual int HasFixedTopology()
See vtkCell3D API for description of this method.
vtkIdType GetPointToOneRingPoints(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:59
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:40
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:43
represent and manipulate point attribute data
Definition: vtkPointData.h:42
represent and manipulate 3D points
Definition: vtkPoints.h:43
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:51
a cell that represents a triangle
Definition: vtkTriangle.h:45
dataset represents arbitrary combinations of all possible cell types
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:86
#define VTK_DEPRECATED_IN_9_0_0(reason)
int vtkIdType
Definition: vtkType.h:332