VTK  9.5.2
vtkTriQuadraticHexahedron.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
53 *
54 * @par Thanks:
55 * Thanks to Soeren Gebbert who developed this class and
56 * integrated it into VTK 5.0.
58 @par Tests:
59 @ref c2_vtk_t_vtkTriQuadraticHexahedron "vtkTriQuadraticHexahedron (Tests)"
60 */
61
62#ifndef vtkTriQuadraticHexahedron_h
63#define vtkTriQuadraticHexahedron_h
64
65#include "vtkCommonDataModelModule.h" // For export macro
66#include "vtkNonLinearCell.h"
67
68VTK_ABI_NAMESPACE_BEGIN
71class vtkHexahedron;
72class vtkDoubleArray;
73
74class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticHexahedron : public vtkNonLinearCell
75{
76public:
79 void PrintSelf(ostream& os, vtkIndent indent) override;
80
82
86 int GetCellType() override { return VTK_TRIQUADRATIC_HEXAHEDRON; }
87 int GetCellDimension() override { return 3; }
88 int GetNumberOfEdges() override { return 12; }
89 int GetNumberOfFaces() override { return 6; }
90 vtkCell* GetEdge(int) override;
91 vtkCell* GetFace(int) override;
93
94 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
95 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
96 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
97 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
98 int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
99 double& dist2, double* weights) override;
100 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
101 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
103 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
104 double* GetParametricCoords() override;
105
111 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
112 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
113 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
114
119 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
120 double pcoords[3], int& subId) override;
121
122 static void InterpolationFunctions(const double pcoords[3], double weights[27]);
123 static void InterpolationDerivs(const double pcoords[3], double derivs[81]);
125
129 void InterpolateFunctions(const double pcoords[3], double weights[27]) override
130 {
132 }
133 void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
134 {
136 }
137
139
146 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
147 static const vtkIdType* GetFaceArray(vtkIdType faceId);
149
155 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[81]);
156
157protected:
160
165
166private:
168 void operator=(const vtkTriQuadraticHexahedron&) = delete;
169};
170
171VTK_ABI_NAMESPACE_END
172#endif
cell represents a parabolic, 9-node isoparametric quad
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:32
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
list of point or cell ids
Definition vtkIdList.h:24
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:29
represent and manipulate point attribute data
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 27-node isoparametric hexahedron
static void InterpolationDerivs(const double pcoords[3], double derivs[81])
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
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 EvaluatePosition(const double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights) override
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 GetCellType() override
Implement the vtkCell API.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
static void InterpolationFunctions(const double pcoords[3], double weights[27])
vtkCell * GetEdge(int) override
Implement the vtkCell API.
vtkCell * GetFace(int) override
Implement the vtkCell API.
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 triquadratic hexahedron using scalar value provided.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetNumberOfEdges() 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 GetCellDimension() 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 vtkTriQuadraticHexahedron * New()
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 InterpolateDerivs(const double pcoords[3], double derivs[81]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
int GetNumberOfFaces() override
Implement the vtkCell API.
void InterpolateFunctions(const double pcoords[3], double weights[27]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
vtkHexahedron * Hex
vtkQuadraticEdge * Edge
vtkBiQuadraticQuadraticHexahedron vtkNonLinearCell JacobianInverse(const double pcoords[3], double **inverse, double derivs[72])
Given parametric coordinates compute inverse Jacobian transformation matrix.
vtkDoubleArray * Scalars
vtkQuadraticQuad * Face
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition vtkCellType.h:65
#define vtkDataArray
~vtkTriQuadraticHexahedron() override
vtkTriQuadraticHexahedron()
int vtkIdType
Definition vtkType.h:332