VTK  9.3.0
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
24 #ifndef vtkQuadraticQuad_h
25 #define vtkQuadraticQuad_h
26 
27 #include "vtkCommonDataModelModule.h" // For export macro
28 #include "vtkNonLinearCell.h"
29 
30 VTK_ABI_NAMESPACE_BEGIN
31 class vtkQuadraticEdge;
32 class vtkQuad;
33 class vtkDoubleArray;
34 
35 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticQuad : public vtkNonLinearCell
36 {
37 public:
38  static vtkQuadraticQuad* New();
40  void PrintSelf(ostream& os, vtkIndent indent) override;
41 
43 
47  int GetCellType() override { return VTK_QUADRATIC_QUAD; }
48  int GetCellDimension() override { return 2; }
49  int GetNumberOfEdges() override { return 4; }
50  int GetNumberOfFaces() override { return 0; }
51  vtkCell* GetEdge(int) override;
52  vtkCell* GetFace(int) override { return nullptr; }
54 
55  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
56  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
57  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
58  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
59  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
60  double& dist2, double weights[]) override;
61  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
62  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
64  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
65  double* GetParametricCoords() override;
66 
71  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
72  vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
73  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
74 
79  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
80  double pcoords[3], int& subId) override;
81 
85  int GetParametricCenter(double pcoords[3]) override;
86 
87  static void InterpolationFunctions(const double pcoords[3], double weights[8]);
88  static void InterpolationDerivs(const double pcoords[3], double derivs[16]);
90 
94  void InterpolateFunctions(const double pcoords[3], double weights[8]) override
95  {
97  }
98  void InterpolateDerivs(const double pcoords[3], double derivs[16]) override
99  {
100  vtkQuadraticQuad::InterpolationDerivs(pcoords, derivs);
101  }
103 
104 protected:
106  ~vtkQuadraticQuad() override;
107 
112 
113  // In order to achieve some functionality we introduce a fake center point
114  // which require to have some extra functionalities compare to other non-linar
115  // cells
118  void Subdivide(double* weights);
120  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
121 
122 private:
123  vtkQuadraticQuad(const vtkQuadraticQuad&) = delete;
124  void operator=(const vtkQuadraticQuad&) = delete;
125 };
126 //----------------------------------------------------------------------------
127 inline int vtkQuadraticQuad::GetParametricCenter(double pcoords[3])
128 {
129  pcoords[0] = pcoords[1] = 0.5;
130  pcoords[2] = 0.;
131  return 0;
132 }
133 
134 VTK_ABI_NAMESPACE_END
135 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:176
represent and manipulate cell attribute data
Definition: vtkCellData.h:31
abstract class to specify cell behavior
Definition: vtkCell.h:50
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:23
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:29
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:30
represent and manipulate 3D points
Definition: vtkPoints.h:29
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:28
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...
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[16])
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
vtkCell * GetFace(int) override
Implement the vtkCell API.
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 Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
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.
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)
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
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_QUADRATIC_QUAD
Definition: vtkCellType.h:58
int vtkIdType
Definition: vtkType.h:315