VTK  9.0.1
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 =========================================================================*/
30 #ifndef vtkConvexPointSet_h
31 #define vtkConvexPointSet_h
32 
33 #include "vtkCell3D.h"
34 #include "vtkCommonDataModelModule.h" // For export macro
35 
37 class vtkCellArray;
38 class vtkTriangle;
39 class vtkTetra;
40 class vtkDoubleArray;
41 
42 class VTKCOMMONDATAMODEL_EXPORT vtkConvexPointSet : public vtkCell3D
43 {
44 public:
46  vtkTypeMacro(vtkConvexPointSet, vtkCell3D);
47  void PrintSelf(ostream& os, vtkIndent indent) override;
48 
52  virtual int HasFixedTopology() { return 0; }
53 
55 
59  void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
60  {
61  vtkWarningMacro(<< "vtkConvexPointSet::GetEdgePoints Not Implemented");
62  }
63  // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
64  VTK_LEGACY(void GetEdgePoints(int vtkNotUsed(edgeId), int*& vtkNotUsed(pts)) override {
65  vtkWarningMacro(<< "vtkConvexPointSet::GetEdgePoints Not Implemented. "
66  << "Also note that this signature is deprecated. "
67  << "Please use GetEdgePoints(vtkIdType, const vtkIdType*& instead");
68  });
69  vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(pts)) override
70  {
71  vtkWarningMacro(<< "vtkConvexPointSet::GetFacePoints Not Implemented");
72  return 0;
73  }
74  // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
75  VTK_LEGACY(void GetFacePoints(int vtkNotUsed(faceId), int*& vtkNotUsed(pts)) override {
76  vtkWarningMacro(<< "vtkConvexPointSet::GetFacePoints Not Implemented. "
77  << "Also note that this signature is deprecated. "
78  << "Please use GetFacePoints(vtkIdType, const vtkIdType*& instead");
79  });
81  vtkIdType vtkNotUsed(edgeId), const vtkIdType*& vtkNotUsed(pts)) override
82  {
83  vtkWarningMacro(<< "vtkConvexPointSet::GetEdgeToAdjacentFaces Not Implemented");
84  }
86  vtkIdType vtkNotUsed(faceId), const vtkIdType*& vtkNotUsed(faceIds)) override
87  {
88  vtkWarningMacro(<< "vtkConvexPointSet::GetFaceToAdjacentFaces Not Implemented");
89  return 0;
90  }
92  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(edgeIds)) override
93  {
94  vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentEdges Not Implemented");
95  return 0;
96  }
98  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(faceIds)) override
99  {
100  vtkWarningMacro(<< "vtkConvexPointSet::GetPointToIncidentFaces Not Implemented");
101  return 0;
102  }
104  vtkIdType vtkNotUsed(pointId), const vtkIdType*& vtkNotUsed(pts)) override
105  {
106  vtkWarningMacro(<< "vtkConvexPointSet::GetPointToOneRingPoints Not Implemented");
107  return 0;
108  }
109  bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
110  {
111  vtkWarningMacro(<< "vtkConvexPointSet::GetCentroid Not Implemented");
112  return 0;
113  }
115 
119  double* GetParametricCoords() override;
120 
124  int GetCellType() override { return VTK_CONVEX_POINT_SET; }
125 
129  int RequiresInitialization() override { return 1; }
130  void Initialize() override;
131 
133 
143  int GetNumberOfEdges() override { return 0; }
144  vtkCell* GetEdge(int) override { return nullptr; }
145  int GetNumberOfFaces() override;
146  vtkCell* GetFace(int faceId) override;
148 
153  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
154  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
155  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
156 
162  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
163  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
164  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
165 
171  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
172  double& dist2, double weights[]) override;
173 
177  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
178 
183  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
184  double pcoords[3], int& subId) override;
185 
189  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
190 
196  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
197 
203  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
204 
208  int GetParametricCenter(double pcoords[3]) override;
209 
214  int IsPrimaryCell() override { return 0; }
215 
217 
221  void InterpolateFunctions(const double pcoords[3], double* sf) override;
222  void InterpolateDerivs(const double pcoords[3], double* derivs) override;
224 
225 protected:
227  ~vtkConvexPointSet() override;
228 
233 
237 
238 private:
239  vtkConvexPointSet(const vtkConvexPointSet&) = delete;
240  void operator=(const vtkConvexPointSet&) = delete;
241 };
242 
243 //----------------------------------------------------------------------------
244 inline int vtkConvexPointSet::GetParametricCenter(double pcoords[3])
245 {
246  pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
247  return 0;
248 }
249 
250 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:34
vtkConvexPointSet::GetCentroid
bool GetCentroid(double vtkNotUsed(centroid)[3]) const override
Definition: vtkConvexPointSet.h:109
vtkConvexPointSet::Derivatives
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...
vtkConvexPointSet::RequiresInitialization
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
Definition: vtkConvexPointSet.h:129
vtkConvexPointSet::New
static vtkConvexPointSet * New()
vtkConvexPointSet::BoundaryTris
vtkCellArray * BoundaryTris
Definition: vtkConvexPointSet.h:234
vtkConvexPointSet::CellBoundary
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...
vtkConvexPointSet::Tetra
vtkTetra * Tetra
Definition: vtkConvexPointSet.h:229
vtkConvexPointSet
a 3D cell defined by a set of convex points
Definition: vtkConvexPointSet.h:43
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:32
vtkX3D::value
@ value
Definition: vtkX3D.h:226
vtkConvexPointSet::vtkConvexPointSet
vtkConvexPointSet()
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkConvexPointSet::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkConvexPointSet::GetFaceToAdjacentFaces
vtkIdType GetFaceToAdjacentFaces(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(faceIds)) override
Definition: vtkConvexPointSet.h:85
vtkConvexPointSet::GetParametricCenter
int GetParametricCenter(double pcoords[3]) override
Return the center of the cell in parametric coordinates.
Definition: vtkConvexPointSet.h:244
vtkCell3D::GetFacePoints
virtual vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts)=0
Get the list of vertices that define a face.
vtkConvexPointSet::ParametricCoords
vtkDoubleArray * ParametricCoords
Definition: vtkConvexPointSet.h:236
vtkConvexPointSet::Contour
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.
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
vtkCell3D.h
vtkConvexPointSet::GetPointToIncidentFaces
vtkIdType GetPointToIncidentFaces(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(faceIds)) override
Definition: vtkConvexPointSet.h:97
vtkConvexPointSet::TetraIds
vtkIdList * TetraIds
Definition: vtkConvexPointSet.h:230
vtkConvexPointSet::TetraScalars
vtkDoubleArray * TetraScalars
Definition: vtkConvexPointSet.h:232
vtkConvexPointSet::Clip
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.
vtkConvexPointSet::TetraPoints
vtkPoints * TetraPoints
Definition: vtkConvexPointSet.h:231
vtkConvexPointSet::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
vtkCell3D::GetEdgePoints
virtual void GetEdgePoints(vtkIdType edgeId, const vtkIdType *&pts)=0
Get the pair of vertices that define an edge.
vtkCell3D
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
VTK_CONVEX_POINT_SET
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:85
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:57
vtkConvexPointSet::InterpolateFunctions
void InterpolateFunctions(const double pcoords[3], double *sf) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
vtkConvexPointSet::Triangle
vtkTriangle * Triangle
Definition: vtkConvexPointSet.h:235
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:34
vtkConvexPointSet::EvaluateLocation
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
The inverse of EvaluatePosition.
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:180
vtkConvexPointSet::~vtkConvexPointSet
~vtkConvexPointSet() override
vtkConvexPointSet::GetEdgePoints
void GetEdgePoints(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
Definition: vtkConvexPointSet.h:59
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:31
vtkConvexPointSet::IsPrimaryCell
int IsPrimaryCell() override
A convex point set is triangulated prior to any operations on it so it is not a primary cell,...
Definition: vtkConvexPointSet.h:214
vtkConvexPointSet::HasFixedTopology
virtual int HasFixedTopology()
See vtkCell3D API for description of this method.
Definition: vtkConvexPointSet.h:52
vtkTriangle
a cell that represents a triangle
Definition: vtkTriangle.h:36
vtkConvexPointSet::InterpolateDerivs
void InterpolateDerivs(const double pcoords[3], double *derivs) override
vtkConvexPointSet::GetPointToIncidentEdges
vtkIdType GetPointToIncidentEdges(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(edgeIds)) override
Definition: vtkConvexPointSet.h:91
vtkConvexPointSet::GetEdge
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
Definition: vtkConvexPointSet.h:144
vtkConvexPointSet::GetCellType
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkConvexPointSet.h:124
vtkConvexPointSet::GetParametricCoords
double * GetParametricCoords() override
See vtkCell3D API for description of this method.
vtkConvexPointSet::GetFace
vtkCell * GetFace(int faceId) override
Return the face cell from the faceId of the cell.
vtkConvexPointSet::Initialize
void Initialize() override
vtkCell::GetParametricCenter
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:36
vtkConvexPointSet::GetPointToOneRingPoints
vtkIdType GetPointToOneRingPoints(vtkIdType vtkNotUsed(pointId), const vtkIdType *&vtkNotUsed(pts)) override
Definition: vtkConvexPointSet.h:103
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:93
vtkConvexPointSet::EvaluatePosition
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Satisfy the vtkCell API.
vtkTetra
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:42
vtkConvexPointSet::IntersectWithLine
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.
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkConvexPointSet::GetEdgeToAdjacentFaces
void GetEdgeToAdjacentFaces(vtkIdType vtkNotUsed(edgeId), const vtkIdType *&vtkNotUsed(pts)) override
Definition: vtkConvexPointSet.h:80
vtkConvexPointSet::GetFacePoints
vtkIdType GetFacePoints(vtkIdType vtkNotUsed(faceId), const vtkIdType *&vtkNotUsed(pts)) override
Definition: vtkConvexPointSet.h:69
vtkConvexPointSet::GetNumberOfEdges
int GetNumberOfEdges() override
A convex point set has no explicit cell edge or faces; however implicitly (after triangulation) it do...
Definition: vtkConvexPointSet.h:143
vtkConvexPointSet::Triangulate
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Triangulate using methods of vtkOrderedTriangulator.