VTK  9.3.0
vtkQuadraticEdge.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
31#ifndef vtkQuadraticEdge_h
32#define vtkQuadraticEdge_h
33
34#include "vtkCommonDataModelModule.h" // For export macro
35#include "vtkNonLinearCell.h"
36
37VTK_ABI_NAMESPACE_BEGIN
38class vtkLine;
39class vtkDoubleArray;
40
41class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticEdge : public vtkNonLinearCell
42{
43public:
46 void PrintSelf(ostream& os, vtkIndent indent) override;
47
52 int GetCellType() override { return VTK_QUADRATIC_EDGE; }
53 int GetCellDimension() override { return 1; }
54 int GetNumberOfEdges() override { return 0; }
55 int GetNumberOfFaces() override { return 0; }
56 vtkCell* GetEdge(int) override { return nullptr; }
57 vtkCell* GetFace(int) override { return nullptr; }
58
59 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
60 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
61 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
62 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
63 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
64 double& dist2, double weights[]) override;
65 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
66 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
68 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
69 double* GetParametricCoords() override;
70
75 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
76 vtkCellArray* lines, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
77 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
78
83 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
84 double pcoords[3], int& subId) override;
85
89 int GetParametricCenter(double pcoords[3]) override;
90
91 static void InterpolationFunctions(const double pcoords[3], double weights[3]);
92 static void InterpolationDerivs(const double pcoords[3], double derivs[3]);
94
98 void InterpolateFunctions(const double pcoords[3], double weights[3]) override
99 {
101 }
102 void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
103 {
105 }
107
108protected:
111
113 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
114
115private:
116 vtkQuadraticEdge(const vtkQuadraticEdge&) = delete;
117 void operator=(const vtkQuadraticEdge&) = delete;
118};
119//----------------------------------------------------------------------------
120inline int vtkQuadraticEdge::GetParametricCenter(double pcoords[3])
121{
122 pcoords[0] = 0.5;
123 pcoords[1] = pcoords[2] = 0.;
124 return 0;
125}
126
127VTK_ABI_NAMESPACE_END
128#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:40
abstract class to specify cell behavior
Definition vtkCell.h:59
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:38
cell represents a 1D line
Definition vtkLine.h:32
abstract superclass for non-linear cells
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:38
cell represents a parabolic, isoparametric edge
int GetCellType() override
Implement the vtkCell API.
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
static void InterpolationFunctions(const double pcoords[3], double weights[3])
void InterpolateFunctions(const double pcoords[3], double weights[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
static vtkQuadraticEdge * New()
vtkDoubleArray * Scalars
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
void InterpolateDerivs(const double pcoords[3], double derivs[3]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
static void InterpolationDerivs(const double pcoords[3], double derivs[3])
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
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...
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic tetra in parametric coordinates.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Line-edge intersection.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
~vtkQuadraticEdge() override
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
int GetNumberOfEdges() override
Return the number of edges in the cell.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
int GetNumberOfFaces() override
Return the number of faces in the cell.
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *lines, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this edge using scalar value provided.
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
Generate contouring primitives.
@ VTK_QUADRATIC_EDGE
Definition vtkCellType.h:56
int vtkIdType
Definition vtkType.h:315