VTK  9.1.0
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=========================================================================*/
55#ifndef vtkQuadraticLinearQuad_h
56#define vtkQuadraticLinearQuad_h
57
58#include "vtkCommonDataModelModule.h" // For export macro
59#include "vtkNonLinearCell.h"
60
62class vtkLine;
63class vtkQuad;
64class vtkDoubleArray;
65
66class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearQuad : public vtkNonLinearCell
67{
68public:
71 void PrintSelf(ostream& os, vtkIndent indent) override;
72
74
78 int GetCellType() override { return VTK_QUADRATIC_LINEAR_QUAD; }
79 int GetCellDimension() override { return 2; }
80 int GetNumberOfEdges() override { return 4; }
81 int GetNumberOfFaces() override { return 0; }
82 vtkCell* GetEdge(int) override;
83 vtkCell* GetFace(int) override { return nullptr; }
85
86 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
87 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
88 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
89 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
90 int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
91 double& dist2, double* weights) override;
92 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
93 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
95 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
96 double* GetParametricCoords() override;
97
102 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
103 vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
104 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
105
110 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
111 double pcoords[3], int& subId) override;
112
116 int GetParametricCenter(double pcoords[3]) override;
117
118 static void InterpolationFunctions(const double pcoords[3], double weights[6]);
119 static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
121
125 void InterpolateFunctions(const double pcoords[3], double weights[6]) override
126 {
128 }
129 void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
130 {
132 }
134
144 static int* GetEdgeArray(vtkIdType edgeId);
145
146protected:
149
154
155private:
157 void operator=(const vtkQuadraticLinearQuad&) = delete;
158};
159//----------------------------------------------------------------------------
161{
162 pcoords[0] = pcoords[1] = 0.5;
163 pcoords[2] = 0.;
164 return 0;
165}
166
167#endif
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:140
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:113
cell represents a 1D line
Definition: vtkLine.h:140
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:95
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.
static vtkQuadraticLinearQuad * New()
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.
static int * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge (edgeId).
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.
~vtkQuadraticLinearQuad() override
int GetCellType() override
Implement the vtkCell API.
vtkCell * GetEdge(int) override
Implement the vtkCell API.
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[12])
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
int GetNumberOfFaces() override
Implement the vtkCell API.
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.
@ VTK_QUADRATIC_LINEAR_QUAD
Definition: vtkCellType.h:115
int vtkIdType
Definition: vtkType.h:332