VTK  9.3.20240328
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
40 #ifndef vtkConvexPointSet_h
41 #define vtkConvexPointSet_h
42 
43 #include "vtkCell3D.h"
44 #include "vtkCommonDataModelModule.h" // For export macro
45 
46 VTK_ABI_NAMESPACE_BEGIN
48 class vtkCellArray;
49 class vtkTriangle;
50 class vtkTetra;
51 class vtkDoubleArray;
52 
53 class VTKCOMMONDATAMODEL_EXPORT vtkConvexPointSet : public vtkCell3D
54 {
55 public:
57  vtkTypeMacro(vtkConvexPointSet, vtkCell3D);
58  void PrintSelf(ostream& os, vtkIndent indent) override;
59 
63 #ifndef VTK_LEGACY_REMOVE
64  virtual vtkTypeBool HasFixedTopology() { return 0; }
65 #endif
66 
68 
72  void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
73  {
74  vtkWarningMacro(<< "vtkConvexPointSet::GetEdgePoints Not Implemented");
75  }
76  vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
77  {
78  vtkWarningMacro(<< "vtkConvexPointSet::GetFacePoints Not Implemented");
79  return 0;
80  }
82  vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
83  {
84  vtkWarningMacro(<< "vtkConvexPointSet::GetEdgeToAdjacentFaces Not Implemented");
85  }
87  vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
88  {
89  vtkWarningMacro(<< "vtkConvexPointSet::GetFaceToAdjacentFaces Not Implemented");
90  return 0;
91  }
93  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
94  {
95  vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentEdges Not Implemented");
96  return 0;
97  }
99  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(faceIds)) override
100  {
101  vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentFaces Not Implemented");
102  return 0;
103  }
105  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
106  {
107  vtkWarningMacro(<< "vtkConvexPointSet::GetPointToOneRingPoints Not Implemented");
108  return 0;
109  }
110  bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
111  {
112  vtkWarningMacro(<< "vtkConvexPointSet::GetCentroid Not Implemented");
113  return false;
114  }
116 
120  double* GetParametricCoords() override;
121 
125  int GetCellType() override { return VTK_CONVEX_POINT_SET; }
126 
130  int RequiresInitialization() override { return 1; }
131  void Initialize() override;
132 
134 
144  int GetNumberOfEdges() override { return 0; }
145  vtkCell* GetEdge(int) override { return nullptr; }
146  int GetNumberOfFaces() override;
147  vtkCell* GetFace(int faceId) override;
149 
154  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
155  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
156  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
157 
163  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
164  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
165  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
166 
172  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
173  double& dist2, double weights[]) override;
174 
178  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
179 
184  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
185  double pcoords[3], int& subId) override;
186 
190  int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
191 
197  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
198 
204  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
205 
209  int GetParametricCenter(double pcoords[3]) override;
210 
215  int IsPrimaryCell() override { return 0; }
216 
218 
222  void InterpolateFunctions(const double pcoords[3], double* sf) override;
223  void InterpolateDerivs(const double pcoords[3], double* derivs) override;
225 
226 protected:
228  ~vtkConvexPointSet() override;
229 
234 
238 
239 private:
240  vtkConvexPointSet(const vtkConvexPointSet&) = delete;
241  void operator=(const vtkConvexPointSet&) = delete;
242 };
243 
244 //----------------------------------------------------------------------------
245 inline int vtkConvexPointSet::GetParametricCenter(double pcoords[3])
246 {
247  pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
248  return 0;
249 }
250 
251 VTK_ABI_NAMESPACE_END
252 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:28
object to represent cell connectivity
Definition: vtkCellArray.h:285
represent and manipulate cell attribute data
Definition: vtkCellData.h:140
abstract class to specify cell behavior
Definition: vtkCell.h:130
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
a 3D cell defined by a set of convex points
vtkDoubleArray * ParametricCoords
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
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
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
vtkCell * GetFace(int faceId) override
A convex point set has no explicit cell edge or faces; however implicitly (after triangulation) it do...
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
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.
static vtkConvexPointSet * New()
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)
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.
vtkCell * GetEdge(int) override
A convex point set has no explicit cell edge or faces; however implicitly (after triangulation) it do...
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:154
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:132
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:108
represent and manipulate point attribute data
Definition: vtkPointData.h:139
represent and manipulate 3D points
Definition: vtkPoints.h:138
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:113
a cell that represents a triangle
Definition: vtkTriangle.h:137
dataset represents arbitrary combinations of all possible cell types
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
int vtkTypeBool
Definition: vtkABI.h:64
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:96
int vtkIdType
Definition: vtkType.h:315