VTK  9.3.0
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
71 #ifndef vtkTriQuadraticHexahedron_h
72 #define vtkTriQuadraticHexahedron_h
73 
74 #include "vtkCommonDataModelModule.h" // For export macro
75 #include "vtkNonLinearCell.h"
76 
77 VTK_ABI_NAMESPACE_BEGIN
78 class vtkQuadraticEdge;
79 class vtkBiQuadraticQuad;
80 class vtkHexahedron;
81 class vtkDoubleArray;
82 
83 class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticHexahedron : public vtkNonLinearCell
84 {
85 public:
88  void PrintSelf(ostream& os, vtkIndent indent) override;
89 
91 
95  int GetCellType() override { return VTK_TRIQUADRATIC_HEXAHEDRON; }
96  int GetCellDimension() override { return 3; }
97  int GetNumberOfEdges() override { return 12; }
98  int GetNumberOfFaces() override { return 6; }
99  vtkCell* GetEdge(int) override;
100  vtkCell* GetFace(int) override;
102 
103  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
104  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
105  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
106  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
107  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
108  double& dist2, double* weights) override;
109  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
110  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
112  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
113  double* GetParametricCoords() override;
114 
120  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
121  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
122  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
123 
128  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
129  double pcoords[3], int& subId) override;
130 
131  static void InterpolationFunctions(const double pcoords[3], double weights[27]);
132  static void InterpolationDerivs(const double pcoords[3], double derivs[81]);
134 
138  void InterpolateFunctions(const double pcoords[3], double weights[27]) override
139  {
141  }
142  void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
143  {
145  }
148 
155  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
156  static const vtkIdType* GetFaceArray(vtkIdType faceId);
158 
164  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[81]);
165 
166 protected:
169 
174 
175 private:
177  void operator=(const vtkTriQuadraticHexahedron&) = delete;
178 };
179 
180 VTK_ABI_NAMESPACE_END
181 #endif
cell represents a parabolic, 9-node isoparametric quad
object to represent cell connectivity
Definition: vtkCellArray.h:185
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
Definition: vtkDataArray.h:54
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:43
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
Definition: vtkPointData.h:39
represent and manipulate 3D points
Definition: vtkPoints.h:38
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 27-node isoparametric hexahedron
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
static void InterpolationDerivs(const double pcoords[3], double derivs[81])
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 JacobianInverse(const double pcoords[3], double **inverse, double derivs[81])
Given parametric coordinates compute inverse Jacobian transformation matrix.
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.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static void InterpolationFunctions(const double pcoords[3], double weights[27])
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.
static vtkTriQuadraticHexahedron * New()
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.
vtkCell * GetFace(int) 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...
void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
int GetNumberOfFaces() override
Implement the vtkCell API.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void InterpolateFunctions(const double pcoords[3], double weights[27]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
~vtkTriQuadraticHexahedron() override
vtkCell * GetEdge(int) override
Implement the vtkCell API.
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:74
int vtkIdType
Definition: vtkType.h:315