VTK  9.4.20241221
vtkQuadraticQuad.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
43#ifndef vtkQuadraticQuad_h
44#define vtkQuadraticQuad_h
45
46#include "vtkCommonDataModelModule.h" // For export macro
47#include "vtkNonLinearCell.h"
48
49VTK_ABI_NAMESPACE_BEGIN
51class vtkQuad;
52class vtkDoubleArray;
53
54class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticQuad : public vtkNonLinearCell
55{
56public:
59 void PrintSelf(ostream& os, vtkIndent indent) override;
60
62
66 int GetCellType() override { return VTK_QUADRATIC_QUAD; }
67 int GetCellDimension() override { return 2; }
68 int GetNumberOfEdges() override { return 4; }
69 int GetNumberOfFaces() override { return 0; }
70 vtkCell* GetEdge(int) override;
71 vtkCell* GetFace(int) override { return nullptr; }
73
74 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
75 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
76 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
77 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
78 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
79 double& dist2, double weights[]) override;
80 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
81 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
83 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
84 double* GetParametricCoords() override;
85
90 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
91 vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
92 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
93
98 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
99 double pcoords[3], int& subId) override;
100
104 int GetParametricCenter(double pcoords[3]) override;
105
106 static void InterpolationFunctions(const double pcoords[3], double weights[8]);
107 static void InterpolationDerivs(const double pcoords[3], double derivs[16]);
109
113 void InterpolateFunctions(const double pcoords[3], double weights[8]) override
114 {
116 }
117 void InterpolateDerivs(const double pcoords[3], double derivs[16]) override
118 {
120 }
122
123protected:
126
131
132 // In order to achieve some functionality we introduce a fake center point
133 // which require to have some extra functionalities compare to other non-linar
134 // cells
137 void Subdivide(double* weights);
139 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
140
141private:
142 vtkQuadraticQuad(const vtkQuadraticQuad&) = delete;
143 void operator=(const vtkQuadraticQuad&) = delete;
144};
145//----------------------------------------------------------------------------
146inline int vtkQuadraticQuad::GetParametricCenter(double pcoords[3])
147{
148 pcoords[0] = pcoords[1] = 0.5;
149 pcoords[2] = 0.;
150 return 0;
151}
152
153VTK_ABI_NAMESPACE_END
154#endif
object to represent cell connectivity
represent and manipulate cell attribute data
abstract class to specify cell behavior
Definition vtkCell.h:130
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:133
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
abstract superclass for non-linear cells
represent and manipulate point attribute data
a cell that represents a 2D quadrilateral
Definition vtkQuad.h:87
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 8-node isoparametric quad
vtkQuadraticEdge * Edge
~vtkQuadraticQuad() override
vtkDoubleArray * CellScalars
int GetNumberOfEdges() override
Implement the vtkCell API.
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this quadratic quad using scalar value provided.
void Subdivide(double *weights)
void InterpolateDerivs(const double pcoords[3], double derivs[16]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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 InterpolationDerivs(const double pcoords[3], double derivs[16])
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
vtkCell * GetEdge(int) override
Implement the vtkCell API.
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 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...
vtkDoubleArray * Scalars
int GetNumberOfFaces() override
Implement the vtkCell API.
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.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
static vtkQuadraticQuad * New()
vtkPointData * PointData
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void InterpolateFunctions(const double pcoords[3], double weights[8]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
void InterpolateAttributes(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
vtkCell * GetFace(int) override
Implement the vtkCell API.
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 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.
int GetCellType() override
Implement the vtkCell API.
static void InterpolationFunctions(const double pcoords[3], double weights[8])
int GetCellDimension() override
Implement the vtkCell API.
vtkCellData * CellData
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
@ VTK_QUADRATIC_QUAD
Definition vtkCellType.h:58
int vtkIdType
Definition vtkType.h:315