VTK  9.3.20240328
vtkQuadraticLinearQuad.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
46 #ifndef vtkQuadraticLinearQuad_h
47 #define vtkQuadraticLinearQuad_h
48 
49 #include "vtkCommonDataModelModule.h" // For export macro
50 #include "vtkNonLinearCell.h"
51 
52 VTK_ABI_NAMESPACE_BEGIN
53 class vtkQuadraticEdge;
54 class vtkLine;
55 class vtkQuad;
56 class vtkDoubleArray;
57 
58 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearQuad : public vtkNonLinearCell
59 {
60 public:
63  void PrintSelf(ostream& os, vtkIndent indent) override;
64 
66 
70  int GetCellType() override { return VTK_QUADRATIC_LINEAR_QUAD; }
71  int GetCellDimension() override { return 2; }
72  int GetNumberOfEdges() override { return 4; }
73  int GetNumberOfFaces() override { return 0; }
74  vtkCell* GetEdge(int) override;
75  vtkCell* GetFace(int) override { return nullptr; }
77 
78  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
79  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
80  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
81  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
82  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
83  double& dist2, double* weights) override;
84  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
85  int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
87  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
88  double* GetParametricCoords() override;
89 
94  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
95  vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
96  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
97 
102  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
103  double pcoords[3], int& subId) override;
104 
108  int GetParametricCenter(double pcoords[3]) override;
109 
110  static void InterpolationFunctions(const double pcoords[3], double weights[6]);
111  static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
113 
117  void InterpolateFunctions(const double pcoords[3], double weights[6]) override
118  {
120  }
121  void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
122  {
124  }
126 
136  static int* GetEdgeArray(vtkIdType edgeId);
137 
138 protected:
141 
146 
147 private:
149  void operator=(const vtkQuadraticLinearQuad&) = delete;
150 };
151 //----------------------------------------------------------------------------
152 inline int vtkQuadraticLinearQuad::GetParametricCenter(double pcoords[3])
153 {
154  pcoords[0] = pcoords[1] = 0.5;
155  pcoords[2] = 0.;
156  return 0;
157 }
158 
159 VTK_ABI_NAMESPACE_END
160 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:285
represent and manipulate cell attribute data
Definition: vtkCellData.h:140
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
Definition: vtkDataArray.h:154
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:132
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:108
cell represents a 1D line
Definition: vtkLine.h:132
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:139
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:87
cell represents a parabolic, isoparametric edge
cell represents a quadratic-linear, 6-node isoparametric quad
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
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 linear quad using scalar value provided.
int GetCellDimension() override
Implement the vtkCell API.
int EvaluatePosition(const double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights) override
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
Compute derivatives given cell subId and parametric coordinates.
vtkCell * GetFace(int) override
Implement the vtkCell API.
void InterpolateFunctions(const double pcoords[3], double weights[6]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
static void InterpolationFunctions(const double pcoords[3], double weights[6])
int GetNumberOfEdges() override
Implement the vtkCell API.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
~vtkQuadraticLinearQuad() override
vtkCell * GetEdge(int) override
Implement the vtkCell API.
int GetCellType() override
Implement the vtkCell API.
static vtkQuadraticLinearQuad * New()
static void InterpolationDerivs(const double pcoords[3], double derivs[12])
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
int GetNumberOfFaces() override
Implement the vtkCell API.
static int * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge (edgeId).
void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
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.
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_QUADRATIC_LINEAR_QUAD
Definition: vtkCellType.h:86
int vtkIdType
Definition: vtkType.h:315