VTK  9.2.6
vtkQuadraticLinearQuad.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticLinearQuad.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 =========================================================================*/
58 #ifndef vtkQuadraticLinearQuad_h
59 #define vtkQuadraticLinearQuad_h
60 
61 #include "vtkCommonDataModelModule.h" // For export macro
62 #include "vtkNonLinearCell.h"
63 
64 class vtkQuadraticEdge;
65 class vtkLine;
66 class vtkQuad;
67 class vtkDoubleArray;
68 
69 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearQuad : public vtkNonLinearCell
70 {
71 public:
74  void PrintSelf(ostream& os, vtkIndent indent) override;
75 
77 
81  int GetCellType() override { return VTK_QUADRATIC_LINEAR_QUAD; }
82  int GetCellDimension() override { return 2; }
83  int GetNumberOfEdges() override { return 4; }
84  int GetNumberOfFaces() override { return 0; }
85  vtkCell* GetEdge(int) override;
86  vtkCell* GetFace(int) override { return nullptr; }
88 
89  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
90  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
91  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
92  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
93  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
94  double& dist2, double* weights) override;
95  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
96  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
98  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
99  double* GetParametricCoords() override;
100 
105  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
106  vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
107  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
108 
113  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
114  double pcoords[3], int& subId) override;
115 
119  int GetParametricCenter(double pcoords[3]) override;
120 
121  static void InterpolationFunctions(const double pcoords[3], double weights[6]);
122  static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
124 
128  void InterpolateFunctions(const double pcoords[3], double weights[6]) override
129  {
131  }
132  void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
133  {
135  }
137 
147  static int* GetEdgeArray(vtkIdType edgeId);
148 
149 protected:
152 
157 
158 private:
160  void operator=(const vtkQuadraticLinearQuad&) = delete;
161 };
162 //----------------------------------------------------------------------------
163 inline int vtkQuadraticLinearQuad::GetParametricCenter(double pcoords[3])
164 {
165  pcoords[0] = pcoords[1] = 0.5;
166  pcoords[2] = 0.;
167  return 0;
168 }
169 
170 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:296
represent and manipulate cell attribute data
Definition: vtkCellData.h:151
abstract class to specify cell behavior
Definition: vtkCell.h:150
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:143
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:119
cell represents a 1D line
Definition: vtkLine.h:143
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:151
represent and manipulate 3D points
Definition: vtkPoints.h:149
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:98
cell represents a parabolic, isoparametric edge
cell represents a quadratic-linear, 6-node isoparametric quad
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 Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
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:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_LINEAR_QUAD
Definition: vtkCellType.h:115
int vtkIdType
Definition: vtkType.h:332