VTK  9.3.0
vtkQuadraticHexahedron.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
34#ifndef vtkQuadraticHexahedron_h
35#define vtkQuadraticHexahedron_h
36
37#include "vtkCommonDataModelModule.h" // For export macro
38#include "vtkNonLinearCell.h"
39
40VTK_ABI_NAMESPACE_BEGIN
43class vtkHexahedron;
44class vtkDoubleArray;
45
46class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticHexahedron : public vtkNonLinearCell
47{
48public:
51 void PrintSelf(ostream& os, vtkIndent indent) override;
52
54
58 int GetCellType() override { return VTK_QUADRATIC_HEXAHEDRON; }
59 int GetCellDimension() override { return 3; }
60 int GetNumberOfEdges() override { return 12; }
61 int GetNumberOfFaces() override { return 6; }
62 vtkCell* GetEdge(int) override;
63 vtkCell* GetFace(int) override;
65
66 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
67 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
68 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
69 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
70 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
71 double& dist2, double weights[]) override;
72 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
73 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
75 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
76 double* GetParametricCoords() override;
77
83 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
84 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
85 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
86
91 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
92 double pcoords[3], int& subId) override;
93
94 static void InterpolationFunctions(const double pcoords[3], double weights[20]);
95 static void InterpolationDerivs(const double pcoords[3], double derivs[60]);
97
101 void InterpolateFunctions(const double pcoords[3], double weights[20]) override
102 {
104 }
105 void InterpolateDerivs(const double pcoords[3], double derivs[60]) override
106 {
108 }
111
118 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
119 static const vtkIdType* GetFaceArray(vtkIdType faceId);
121
127 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[60]);
128
129protected:
132
140
142 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
143
144private:
146 void operator=(const vtkQuadraticHexahedron&) = delete;
147};
148
149VTK_ABI_NAMESPACE_END
150#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:40
abstract class to specify cell behavior
Definition vtkCell.h:59
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
list of point or cell ids
Definition vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:38
abstract superclass for non-linear cells
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:38
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 20-node isoparametric hexahedron
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.
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[60])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
vtkCell * GetEdge(int) override
Implement the vtkCell API.
int GetNumberOfEdges() override
Implement the vtkCell API.
static void InterpolationFunctions(const double pcoords[3], double weights[20])
void InterpolateDerivs(const double pcoords[3], double derivs[60]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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 InterpolateFunctions(const double pcoords[3], double weights[20]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetCellType() override
Implement the vtkCell API.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
int GetCellDimension() override
Implement the vtkCell API.
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...
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.
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 const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static vtkQuadraticHexahedron * New()
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
static void InterpolationDerivs(const double pcoords[3], double derivs[60])
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 PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tetras, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this quadratic hexahedron using scalar value provided.
~vtkQuadraticHexahedron() override
vtkCell * GetFace(int) override
Implement the vtkCell API.
int GetNumberOfFaces() override
Implement the vtkCell API.
cell represents a parabolic, 8-node isoparametric quad
@ VTK_QUADRATIC_HEXAHEDRON
Definition vtkCellType.h:61
int vtkIdType
Definition vtkType.h:315