VTK  9.2.5
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=========================================================================*/
25#ifndef vtkVertex_h
26#define vtkVertex_h
27
28#include "vtkCell.h"
29#include "vtkCommonDataModelModule.h" // For export macro
30
32
33class VTKCOMMONDATAMODEL_EXPORT vtkVertex : public vtkCell
34{
35public:
36 static vtkVertex* New();
37 vtkTypeMacro(vtkVertex, vtkCell);
38 void PrintSelf(ostream& os, vtkIndent indent) override;
39
45
48 int GetCellType() override { return VTK_VERTEX; }
49 int GetCellDimension() override { return 0; }
50 int GetNumberOfEdges() override { return 0; }
51 int GetNumberOfFaces() override { return 0; }
52 vtkCell* GetEdge(int) override { return nullptr; }
53 vtkCell* GetFace(int) override { return nullptr; }
54 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
55 vtkCellArray* pts, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId,
56 vtkCellData* outCd, int insideOut) override;
57 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
58 double& dist2, double weights[]) override;
59 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
60 double* GetParametricCoords() override;
62
68 int Inflate(double) override { return 0; }
69
77 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
78
85 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
86 vtkCellArray* verts1, vtkCellArray* lines, vtkCellArray* verts2, vtkPointData* inPd,
87 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
88
92 int GetParametricCenter(double pcoords[3]) override;
93
99 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
100 double pcoords[3], int& subId) override;
101
106 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
107
113 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
114
115 static void InterpolationFunctions(const double pcoords[3], double weights[1]);
116 static void InterpolationDerivs(const double pcoords[3], double derivs[3]);
118
122 void InterpolateFunctions(const double pcoords[3], double weights[1]) override
123 {
124 vtkVertex::InterpolationFunctions(pcoords, weights);
125 }
126 void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
127 {
128 vtkVertex::InterpolationDerivs(pcoords, derivs);
129 }
131
132protected:
134 ~vtkVertex() override = default;
135
136private:
137 vtkVertex(const vtkVertex&) = delete;
138 void operator=(const vtkVertex&) = delete;
139};
140
141//----------------------------------------------------------------------------
142inline int vtkVertex::GetParametricCenter(double pcoords[3])
143{
144 pcoords[0] = pcoords[1] = pcoords[2] = 0.0;
145 return 0;
146}
147
148#endif
object to represent cell connectivity
Definition: vtkCellArray.h:187
represent and manipulate cell attribute data
Definition: vtkCellData.h:42
abstract class to specify cell behavior
Definition: vtkCell.h:61
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:56
list of point or cell ids
Definition: vtkIdList.h:34
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:42
represent and manipulate 3D points
Definition: vtkPoints.h:40
a cell that represents a 3D point
Definition: vtkVertex.h:34
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:52
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:51
vtkCell * GetFace(int) override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:53
int GetCellDimension() 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:142
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:68
void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkVertex.h:126
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:50
static vtkVertex * New()
int GetCellType() override
Make a new vtkVertex object with the same information as this object.
Definition: vtkVertex.h:48
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:122
@ VTK_VERTEX
Definition: vtkCellType.h:47
int vtkIdType
Definition: vtkType.h:332