VTK  9.1.0
vtkVertex.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkVertex.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 =========================================================================*/
22 #ifndef vtkVertex_h
23 #define vtkVertex_h
24 
25 #include "vtkCell.h"
26 #include "vtkCommonDataModelModule.h" // For export macro
27 
29 
30 class VTKCOMMONDATAMODEL_EXPORT vtkVertex : public vtkCell
31 {
32 public:
33  static vtkVertex* New();
34  vtkTypeMacro(vtkVertex, vtkCell);
35  void PrintSelf(ostream& os, vtkIndent indent) override;
36 
42 
45  int GetCellType() override { return VTK_VERTEX; }
46  int GetCellDimension() override { return 0; }
47  int GetNumberOfEdges() override { return 0; }
48  int GetNumberOfFaces() override { return 0; }
49  vtkCell* GetEdge(int) override { return nullptr; }
50  vtkCell* GetFace(int) override { return nullptr; }
51  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
52  vtkCellArray* pts, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId,
53  vtkCellData* outCd, int insideOut) override;
54  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
55  double& dist2, double weights[]) override;
56  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
57  double* GetParametricCoords() override;
59 
65  int Inflate(double) override { return 0; }
66 
74  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
75 
82  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
83  vtkCellArray* verts1, vtkCellArray* lines, vtkCellArray* verts2, vtkPointData* inPd,
84  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
85 
89  int GetParametricCenter(double pcoords[3]) override;
90 
96  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
97  double pcoords[3], int& subId) override;
98 
103  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
104 
110  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
111 
112  static void InterpolationFunctions(const double pcoords[3], double weights[1]);
113  static void InterpolationDerivs(const double pcoords[3], double derivs[3]);
115 
119  void InterpolateFunctions(const double pcoords[3], double weights[1]) override
120  {
121  vtkVertex::InterpolationFunctions(pcoords, weights);
122  }
123  void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
124  {
125  vtkVertex::InterpolationDerivs(pcoords, derivs);
126  }
128 
129 protected:
131  ~vtkVertex() override = default;
132 
133 private:
134  vtkVertex(const vtkVertex&) = delete;
135  void operator=(const vtkVertex&) = delete;
136 };
137 
138 //----------------------------------------------------------------------------
139 inline int vtkVertex::GetParametricCenter(double pcoords[3])
140 {
141  pcoords[0] = pcoords[1] = pcoords[2] = 0.0;
142  return 0;
143 }
144 
145 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:181
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:58
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
list of point or cell ids
Definition: vtkIdList.h:31
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:34
represent and manipulate point attribute data
Definition: vtkPointData.h:33
represent and manipulate 3D points
Definition: vtkPoints.h:34
a cell that represents a 3D point
Definition: vtkVertex.h:31
void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts1, vtkCellArray *lines, vtkCellArray *verts2, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
Generate contouring primitives.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Make a new vtkVertex object with the same information as this object.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Intersect with a ray.
static void InterpolationFunctions(const double pcoords[3], double weights[1])
int GetNumberOfFaces() override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:48
int GetCellDimension() override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:46
vtkCell * GetEdge(int) override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:49
int GetParametricCenter(double pcoords[3]) override
Return the center of the triangle in parametric coordinates.
Definition: vtkVertex.h:139
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Triangulate the vertex.
double * GetParametricCoords() override
Make a new vtkVertex object with the same information as this object.
~vtkVertex() override=default
vtkCell * GetFace(int) override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:50
int Inflate(double) override
This method does nothing.
Definition: vtkVertex.h:65
void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkVertex.h:123
static vtkVertex * New()
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Make a new vtkVertex object with the same information as this object.
int GetNumberOfEdges() override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:47
int GetCellType() override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:45
static void InterpolationDerivs(const double pcoords[3], double derivs[3])
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Get the derivative of the vertex.
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *pts, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Make a new vtkVertex object with the same information as this object.
void InterpolateFunctions(const double pcoords[3], double weights[1]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkVertex.h:119
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_VERTEX
Definition: vtkCellType.h:47
int vtkIdType
Definition: vtkType.h:332