VTK  9.1.0
vtkQuadraticQuad.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkQuadraticQuad.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=========================================================================*/
42#ifndef vtkQuadraticQuad_h
43#define vtkQuadraticQuad_h
44
45#include "vtkCommonDataModelModule.h" // For export macro
46#include "vtkNonLinearCell.h"
47
49class vtkQuad;
50class vtkDoubleArray;
51
52class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticQuad : public vtkNonLinearCell
53{
54public:
57 void PrintSelf(ostream& os, vtkIndent indent) override;
58
60
64 int GetCellType() override { return VTK_QUADRATIC_QUAD; }
65 int GetCellDimension() override { return 2; }
66 int GetNumberOfEdges() override { return 4; }
67 int GetNumberOfFaces() override { return 0; }
68 vtkCell* GetEdge(int) override;
69 vtkCell* GetFace(int) override { return nullptr; }
71
72 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
73 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
74 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
75 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
76 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
77 double& dist2, double weights[]) override;
78 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
79 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
81 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
82 double* GetParametricCoords() override;
83
88 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
89 vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
90 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
91
96 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
97 double pcoords[3], int& subId) override;
98
102 int GetParametricCenter(double pcoords[3]) override;
103
104 static void InterpolationFunctions(const double pcoords[3], double weights[8]);
105 static void InterpolationDerivs(const double pcoords[3], double derivs[16]);
107
111 void InterpolateFunctions(const double pcoords[3], double weights[8]) override
112 {
114 }
115 void InterpolateDerivs(const double pcoords[3], double derivs[16]) override
116 {
118 }
120
121protected:
124
129
130 // In order to achieve some functionality we introduce a fake center point
131 // which require to have some extra functionalities compare to other non-linar
132 // cells
135 void Subdivide(double* weights);
137 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
138
139private:
140 vtkQuadraticQuad(const vtkQuadraticQuad&) = delete;
141 void operator=(const vtkQuadraticQuad&) = delete;
142};
143//----------------------------------------------------------------------------
144inline int vtkQuadraticQuad::GetParametricCenter(double pcoords[3])
145{
146 pcoords[0] = pcoords[1] = 0.5;
147 pcoords[2] = 0.;
148 return 0;
149}
150
151#endif
object to represent cell connectivity
Definition: vtkCellArray.h:190
represent and manipulate cell attribute data
Definition: vtkCellData.h:42
abstract class to specify cell behavior
Definition: vtkCell.h:67
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:59
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:40
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:43
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:42
represent and manipulate 3D points
Definition: vtkPoints.h:43
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:45
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 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.
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
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_QUAD
Definition: vtkCellType.h:67
int vtkIdType
Definition: vtkType.h:332