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
40class vtkCellArray;
41class vtkTriangle;
42class vtkTetra;
43class vtkDoubleArray;
44
45class VTKCOMMONDATAMODEL_EXPORT vtkConvexPointSet : public vtkCell3D
46{
47public:
50 void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
51
55 virtual int HasFixedTopology() {return 0;}
56
60 void GetEdgePoints(int vtkNotUsed(edgeId), int* &vtkNotUsed(pts)) VTK_OVERRIDE {}
61 void GetFacePoints(int vtkNotUsed(faceId), int* &vtkNotUsed(pts)) VTK_OVERRIDE {}
62 double *GetParametricCoords() VTK_OVERRIDE;
63
67 int GetCellType() VTK_OVERRIDE {return VTK_CONVEX_POINT_SET;}
68
72 int RequiresInitialization() VTK_OVERRIDE {return 1;}
73 void Initialize() VTK_OVERRIDE;
74
76
86 int GetNumberOfEdges() VTK_OVERRIDE {return 0;}
87 vtkCell *GetEdge(int) VTK_OVERRIDE {return NULL;}
88 int GetNumberOfFaces() VTK_OVERRIDE;
89 vtkCell *GetFace(int faceId) VTK_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) VTK_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) VTK_OVERRIDE;
112
118 int EvaluatePosition(double x[3], double* closestPoint,
119 int& subId, double pcoords[3],
120 double& dist2, double *weights) VTK_OVERRIDE;
121
125 void EvaluateLocation(int& subId, double pcoords[3], double x[3],
126 double *weights) VTK_OVERRIDE;
127
132 int IntersectWithLine(double p1[3], double p2[3], double tol, double& t,
133 double x[3], double pcoords[3], int& subId) VTK_OVERRIDE;
134
138 int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) VTK_OVERRIDE;
139
144 void Derivatives(int subId, double pcoords[3], double *values,
145 int dim, double *derivs) VTK_OVERRIDE;
146
152 int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) VTK_OVERRIDE;
153
157 int GetParametricCenter(double pcoords[3]) VTK_OVERRIDE;
158
163 int IsPrimaryCell() VTK_OVERRIDE {return 0;}
164
166
170 void InterpolateFunctions(double pcoords[3], double *sf) VTK_OVERRIDE;
171 void InterpolateDerivs(double pcoords[3], double *derivs) VTK_OVERRIDE;
173
174protected:
176 ~vtkConvexPointSet() VTK_OVERRIDE;
177
178 vtkTetra *Tetra;
179 vtkIdList *TetraIds;
180 vtkPoints *TetraPoints;
181 vtkDoubleArray *TetraScalars;
182
183 vtkCellArray *BoundaryTris;
184 vtkTriangle *Triangle;
185 vtkDoubleArray *ParametricCoords;
186
187private:
188 vtkConvexPointSet(const vtkConvexPointSet&) VTK_DELETE_FUNCTION;
189 void operator=(const vtkConvexPointSet&) VTK_DELETE_FUNCTION;
190};
191
192//----------------------------------------------------------------------------
193inline 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
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
object to represent cell connectivity
Definition: vtkCellArray.h:51
represent and manipulate cell attribute data
Definition: vtkCellData.h:39
abstract class to specify cell behavior
Definition: vtkCell.h:60
a 3D cell defined by a set of convex points
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int RequiresInitialization() override
This cell requires that it be initialized prior to access.
static vtkConvexPointSet * New()
void InterpolateFunctions(double pcoords[3], double *sf) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
~vtkConvexPointSet() override
void Initialize() override
int GetNumberOfFaces() override
Return the number of faces in the cell.
void InterpolateDerivs(double pcoords[3], double *derivs) override
void GetEdgePoints(int vtkNotUsed(edgeId), int *&vtkNotUsed(pts)) override
See vtkCell3D API for description of these methods.
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
virtual int HasFixedTopology()
See vtkCell3D API for description of this method.
void GetFacePoints(int vtkNotUsed(faceId), int *&vtkNotUsed(pts)) override
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:37
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:40
represent and manipulate point attribute data
Definition: vtkPointData.h:38
represent and manipulate 3D points
Definition: vtkPoints.h:40
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:48
a cell that represents a triangle
Definition: vtkTriangle.h:42
dataset represents arbitrary combinations of all possible cell types
int Contour(vtkDataSet *input, vtkPolyData *output, vtkDataArray *field, float isoValue, bool computeScalars)
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:84
int vtkIdType
Definition: vtkType.h:287