VTK  9.3.20230924
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
81 #ifndef vtkTriQuadraticHexahedron_h
82 #define vtkTriQuadraticHexahedron_h
83 
84 #include "vtkCommonDataModelModule.h" // For export macro
85 #include "vtkNonLinearCell.h"
86 
87 VTK_ABI_NAMESPACE_BEGIN
88 class vtkQuadraticEdge;
89 class vtkBiQuadraticQuad;
90 class vtkHexahedron;
91 class vtkDoubleArray;
92 
93 class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticHexahedron : public vtkNonLinearCell
94 {
95 public:
98  void PrintSelf(ostream& os, vtkIndent indent) override;
99 
101 
105  int GetCellType() override { return VTK_TRIQUADRATIC_HEXAHEDRON; }
106  int GetCellDimension() override { return 3; }
107  int GetNumberOfEdges() override { return 12; }
108  int GetNumberOfFaces() override { return 6; }
109  vtkCell* GetEdge(int) override;
110  vtkCell* GetFace(int) override;
112 
113  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
114  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
115  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
116  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
117  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
118  double& dist2, double* weights) override;
119  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
120  int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
122  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
123  double* GetParametricCoords() override;
124 
130  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
131  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
132  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
133 
138  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
139  double pcoords[3], int& subId) override;
140 
141  static void InterpolationFunctions(const double pcoords[3], double weights[27]);
142  static void InterpolationDerivs(const double pcoords[3], double derivs[81]);
144 
148  void InterpolateFunctions(const double pcoords[3], double weights[27]) override
149  {
151  }
152  void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
153  {
155  }
158 
165  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
166  static const vtkIdType* GetFaceArray(vtkIdType faceId);
168 
174  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[81]);
175 
176 protected:
179 
184 
185 private:
187  void operator=(const vtkTriQuadraticHexahedron&) = delete;
188 };
189 
190 VTK_ABI_NAMESPACE_END
191 #endif
cell represents a parabolic, 9-node isoparametric quad
object to represent cell connectivity
Definition: vtkCellArray.h:285
represent and manipulate cell attribute data
Definition: vtkCellData.h:140
abstract class to specify cell behavior
Definition: vtkCell.h:129
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:154
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
list of point or cell ids
Definition: vtkIdList.h:132
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:108
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:139
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 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:84
int vtkIdType
Definition: vtkType.h:315