VTK
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 =========================================================================*/
33 #ifndef vtkConvexPointSet_h
34 #define vtkConvexPointSet_h
35 
36 #include "vtkCommonDataModelModule.h" // For export macro
37 #include "vtkCell3D.h"
38 
40 class vtkCellArray;
41 class vtkTriangle;
42 class vtkTetra;
43 class vtkDoubleArray;
44 
45 class VTKCOMMONDATAMODEL_EXPORT vtkConvexPointSet : public vtkCell3D
46 {
47 public:
50  void PrintSelf(ostream& os, vtkIndent indent) override;
51 
55  virtual int HasFixedTopology() {return 0;}
56 
60  void GetEdgePoints(int vtkNotUsed(edgeId), int* &vtkNotUsed(pts)) override {}
61  void GetFacePoints(int vtkNotUsed(faceId), int* &vtkNotUsed(pts)) override {}
62  double *GetParametricCoords() override;
63 
67  int GetCellType() override {return VTK_CONVEX_POINT_SET;}
68 
72  int RequiresInitialization() override {return 1;}
73  void Initialize() override;
74 
76 
86  int GetNumberOfEdges() override {return 0;}
87  vtkCell *GetEdge(int) override {return nullptr;}
88  int GetNumberOfFaces() override;
89  vtkCell *GetFace(int faceId) override;
91 
96  void Contour(double value, vtkDataArray *cellScalars,
98  vtkCellArray *lines, vtkCellArray *polys,
99  vtkPointData *inPd, vtkPointData *outPd,
100  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override;
101 
107  void Clip(double value, vtkDataArray *cellScalars,
108  vtkIncrementalPointLocator *locator, vtkCellArray *connectivity,
109  vtkPointData *inPd, vtkPointData *outPd,
110  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
111  int insideOut) override;
112 
118  int EvaluatePosition(const double x[3], double closestPoint[3],
119  int& subId, double pcoords[3],
120  double& dist2, double weights[]) override;
121 
125  void EvaluateLocation(int& subId, const double pcoords[3], double x[3],
126  double *weights) override;
127 
132  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
133  double x[3], double pcoords[3], int& subId) override;
134 
138  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
139 
144  void Derivatives(int subId, const double pcoords[3], const double *values,
145  int dim, double *derivs) override;
146 
152  int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
153 
157  int GetParametricCenter(double pcoords[3]) override;
158 
163  int IsPrimaryCell() override {return 0;}
164 
166 
170  void InterpolateFunctions(const double pcoords[3], double *sf) override;
171  void InterpolateDerivs(const double pcoords[3], double *derivs) override;
173 
174 protected:
176  ~vtkConvexPointSet() override;
177 
182 
186 
187 private:
188  vtkConvexPointSet(const vtkConvexPointSet&) = delete;
189  void operator=(const vtkConvexPointSet&) = delete;
190 };
191 
192 //----------------------------------------------------------------------------
193 inline int vtkConvexPointSet::GetParametricCenter(double pcoords[3])
194 {
195  pcoords[0] = pcoords[1] = pcoords[2] = 0.5;
196  return 0;
197 }
198 
199 #endif
200 
201 
202 
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:40
vtkConvexPointSet::GetFacePoints
void GetFacePoints(int vtkNotUsed(faceId), int *&vtkNotUsed(pts)) override
Definition: vtkConvexPointSet.h:61
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:72
vtkConvexPointSet::New
static vtkConvexPointSet * New()
vtkConvexPointSet::BoundaryTris
vtkCellArray * BoundaryTris
Definition: vtkConvexPointSet.h:183
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:178
vtkConvexPointSet
a 3D cell defined by a set of convex points
Definition: vtkConvexPointSet.h:46
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:38
vtkX3D::value
@ value
Definition: vtkX3D.h:220
vtkConvexPointSet::vtkConvexPointSet
vtkConvexPointSet()
vtkIdType
int vtkIdType
Definition: vtkType.h:347
vtkConvexPointSet::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkConvexPointSet::GetParametricCenter
int GetParametricCenter(double pcoords[3]) override
Return the center of the cell in parametric coordinates.
Definition: vtkConvexPointSet.h:193
vtkConvexPointSet::GetEdgePoints
void GetEdgePoints(int vtkNotUsed(edgeId), int *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
Definition: vtkConvexPointSet.h:60
vtkConvexPointSet::ParametricCoords
vtkDoubleArray * ParametricCoords
Definition: vtkConvexPointSet.h:185
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:55
vtkCell3D.h
vtkConvexPointSet::TetraIds
vtkIdList * TetraIds
Definition: vtkConvexPointSet.h:179
vtkConvexPointSet::TetraScalars
vtkDoubleArray * TetraScalars
Definition: vtkConvexPointSet.h:181
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:180
vtkConvexPointSet::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
vtkCell3D
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
VTK_CONVEX_POINT_SET
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:84
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:60
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:184
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:39
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:40
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:51
vtkConvexPointSet::~vtkConvexPointSet
~vtkConvexPointSet() override
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:37
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:163
vtkConvexPointSet::HasFixedTopology
virtual int HasFixedTopology()
See vtkCell3D API for description of this method.
Definition: vtkConvexPointSet.h:55
vtkTriangle
a cell that represents a triangle
Definition: vtkTriangle.h:42
vtkConvexPointSet::InterpolateDerivs
void InterpolateDerivs(const double pcoords[3], double *derivs) override
vtkConvexPointSet::GetEdge
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
Definition: vtkConvexPointSet.h:87
vtkConvexPointSet::GetCellType
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkConvexPointSet.h:67
vtkConvexPointSet::GetParametricCoords
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
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:42
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:89
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:48
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:246
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:86
vtkConvexPointSet::Triangulate
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Triangulate using methods of vtkOrderedTriangulator.