VTK  9.3.0
vtkBiQuadraticQuadraticHexahedron.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
59#ifndef vtkBiQuadraticQuadraticHexahedron_h
60#define vtkBiQuadraticQuadraticHexahedron_h
61
62#include "vtkCommonDataModelModule.h" // For export macro
63#include "vtkNonLinearCell.h"
64
65VTK_ABI_NAMESPACE_BEGIN
69class vtkHexahedron;
70class vtkDoubleArray;
71
72class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticHexahedron : public vtkNonLinearCell
73{
74public:
77 void PrintSelf(ostream& os, vtkIndent indent) override;
78
80
85 int GetCellDimension() override { return 3; }
86 int GetNumberOfEdges() override { return 12; }
87 int GetNumberOfFaces() override { return 6; }
88 vtkCell* GetEdge(int) override;
89 vtkCell* GetFace(int) override;
91
92 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
93 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
94 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
95 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
96 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
97 double& dist2, double weights[]) override;
98 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
99 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
101 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
102 double* GetParametricCoords() override;
103
109 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
110 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
111 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
112
117 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
118 double pcoords[3], int& subId) override;
119
120 static void InterpolationFunctions(const double pcoords[3], double weights[24]);
121 static void InterpolationDerivs(const double pcoords[3], double derivs[72]);
123
127 void InterpolateFunctions(const double pcoords[3], double weights[24]) override
128 {
130 }
131 void InterpolateDerivs(const double pcoords[3], double derivs[72]) override
132 {
134 }
137
144 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
145 static const vtkIdType* GetFaceArray(vtkIdType faceId);
147
153 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[72]);
154
155protected:
158
167
169 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
170
171private:
173 void operator=(const vtkBiQuadraticQuadraticHexahedron&) = delete;
174};
175
176VTK_ABI_NAMESPACE_END
177#endif
cell represents a parabolic, 9-node isoparametric quad
cell represents a biquadratic, 24-node isoparametric hexahedron
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[72])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void InterpolateDerivs(const double pcoords[3], double derivs[72]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
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 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...
static vtkBiQuadraticQuadraticHexahedron * New()
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
int GetNumberOfFaces() override
Implement the vtkCell API.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void InterpolateFunctions(const double pcoords[3], double weights[24]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
static void InterpolationFunctions(const double pcoords[3], double weights[24])
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...
vtkCell * GetEdge(int) override
Implement the vtkCell API.
int GetCellDimension() 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.
int GetNumberOfEdges() override
Implement the vtkCell API.
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 void InterpolationDerivs(const double pcoords[3], double derivs[72])
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkCell * GetFace(int) override
Implement the vtkCell API.
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.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
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 biquadratic hexahedron using scalar value provided.
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:31
abstract class to specify cell behavior
Definition vtkCell.h:50
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:23
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:29
abstract superclass for non-linear cells
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:29
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 8-node isoparametric quad
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
Definition vtkCellType.h:70
int vtkIdType
Definition vtkType.h:315