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=========================================================================*/
91#ifndef vtkVertex_h
92#define vtkVertex_h
93
94#include "vtkCell.h"
95#include "vtkCommonDataModelModule.h" // For export macro
96
98
99class VTKCOMMONDATAMODEL_EXPORT vtkVertex : public vtkCell
100{
101public:
102 static vtkVertex* New();
103 vtkTypeMacro(vtkVertex, vtkCell);
104 void PrintSelf(ostream& os, vtkIndent indent) override;
105
111
114 int GetCellType() override { return VTK_VERTEX; }
115 int GetCellDimension() override { return 0; }
116 int GetNumberOfEdges() override { return 0; }
117 int GetNumberOfFaces() override { return 0; }
118 vtkCell* GetEdge(int) override { return nullptr; }
119 vtkCell* GetFace(int) override { return nullptr; }
120 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
121 vtkCellArray* pts, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId,
122 vtkCellData* outCd, int insideOut) override;
123 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
124 double& dist2, double weights[]) override;
125 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
126 double* GetParametricCoords() override;
128
134 int Inflate(double) override { return 0; }
135
143 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
144
151 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
152 vtkCellArray* verts1, vtkCellArray* lines, vtkCellArray* verts2, vtkPointData* inPd,
153 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
154
158 int GetParametricCenter(double pcoords[3]) override;
159
165 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
166 double pcoords[3], int& subId) override;
167
172 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
173
179 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
180
181 static void InterpolationFunctions(const double pcoords[3], double weights[1]);
182 static void InterpolationDerivs(const double pcoords[3], double derivs[3]);
184
188 void InterpolateFunctions(const double pcoords[3], double weights[1]) override
189 {
190 vtkVertex::InterpolationFunctions(pcoords, weights);
191 }
192 void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
193 {
194 vtkVertex::InterpolationDerivs(pcoords, derivs);
195 }
197
198protected:
200 ~vtkVertex() override = default;
201
202private:
203 vtkVertex(const vtkVertex&) = delete;
204 void operator=(const vtkVertex&) = delete;
205};
206
207//----------------------------------------------------------------------------
208inline int vtkVertex::GetParametricCenter(double pcoords[3])
209{
210 pcoords[0] = pcoords[1] = pcoords[2] = 0.0;
211 return 0;
212}
213
214#endif
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
list of point or cell ids
Definition: vtkIdList.h:140
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:113
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
a cell that represents a 3D point
Definition: vtkVertex.h:100
double * GetParametricCoords() override
Make a new vtkVertex object with the same information as this object.
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.
vtkCell * GetEdge(int) override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:118
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:117
vtkCell * GetFace(int) override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:119
int GetCellDimension() override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:115
int GetParametricCenter(double pcoords[3]) override
Return the center of the triangle in parametric coordinates.
Definition: vtkVertex.h:208
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Triangulate the vertex.
~vtkVertex() override=default
int Inflate(double) override
This method does nothing.
Definition: vtkVertex.h:134
void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkVertex.h:192
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:116
static vtkVertex * New()
int GetCellType() override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:114
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:188
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_VERTEX
Definition: vtkCellType.h:86
int vtkIdType
Definition: vtkType.h:332