VTK  9.3.20240424
vtkHexahedron.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
91#ifndef vtkHexahedron_h
92#define vtkHexahedron_h
93
94#include "vtkCell3D.h"
95#include "vtkCommonDataModelModule.h" // For export macro
96
97VTK_ABI_NAMESPACE_BEGIN
98class vtkLine;
99class vtkQuad;
101
102class VTKCOMMONDATAMODEL_EXPORT vtkHexahedron : public vtkCell3D
103{
104public:
106 vtkTypeMacro(vtkHexahedron, vtkCell3D);
107 void PrintSelf(ostream& os, vtkIndent indent) override;
108
110
113 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
114 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
115 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
116 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
117 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
118 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
119 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
120 bool GetCentroid(double centroid[3]) const override;
122
126 static constexpr vtkIdType NumberOfPoints = 8;
127
131 static constexpr vtkIdType NumberOfEdges = 12;
132
136 static constexpr vtkIdType NumberOfFaces = 6;
137
142 static constexpr vtkIdType MaximumFaceSize = 4;
143
149 static constexpr vtkIdType MaximumValence = 3;
150
152
155 int GetCellType() override { return VTK_HEXAHEDRON; }
156 int GetNumberOfEdges() override { return 12; }
157 int GetNumberOfFaces() override { return 6; }
158 vtkCell* GetEdge(int edgeId) override;
159 vtkCell* GetFace(int faceId) override;
160 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
161 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
162 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
163 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
165
166 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
167 double& dist2, double weights[]) override;
168 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
169 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
170 double pcoords[3], int& subId) override;
171 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
173 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
174 double* GetParametricCoords() override;
175
183 static int* GetTriangleCases(int caseId);
184
185 static void InterpolationFunctions(const double pcoords[3], double weights[8]);
186 static void InterpolationDerivs(const double pcoords[3], double derivs[24]);
188
192 void InterpolateFunctions(const double pcoords[3], double weights[8]) override
193 {
195 }
196 void InterpolateDerivs(const double pcoords[3], double derivs[24]) override
197 {
198 vtkHexahedron::InterpolationDerivs(pcoords, derivs);
199 }
201
203
214
219
224
229
234
239
243 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
244
250 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[24]);
251
252protected:
254 ~vtkHexahedron() override;
255
258
259private:
260 vtkHexahedron(const vtkHexahedron&) = delete;
261 void operator=(const vtkHexahedron&) = delete;
262};
263
264VTK_ABI_NAMESPACE_END
265#endif
abstract class to specify 3D cell interface
Definition vtkCell3D.h:28
object to represent cell connectivity
represent and manipulate cell attribute data
abstract class to specify cell behavior
Definition vtkCell.h:130
abstract superclass for arrays of numeric data
a cell that represents a linear 3D hexahedron
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
void InterpolateDerivs(const double pcoords[3], double derivs[24]) 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.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Intersect with a ray.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
void InterpolateFunctions(const double pcoords[3], double weights[8]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
bool GetCentroid(double centroid[3]) const override
See vtkCell3D API for description of these methods.
static int * GetTriangleCases(int caseId)
Return the case table for table-based isocontouring (aka marching cubes style implementations).
static const vtkIdType * GetPointToIncidentEdgesArray(vtkIdType pointId)
Static method version of GetPointToIncidentEdgesArray.
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void GetEdgePoints(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
static void InterpolationFunctions(const double pcoords[3], double weights[8])
vtkCell * GetFace(int faceId) override
See the vtkCell API for descriptions of these methods.
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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
See the vtkCell API for descriptions of these methods.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static bool ComputeCentroid(vtkPoints *points, const vtkIdType *pointIds, double centroid[3])
Static method version of GetCentroid.
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
See vtkCell3D API for description of these methods.
~vtkHexahedron() override
static void InterpolationDerivs(const double pcoords[3], double derivs[24])
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
static vtkHexahedron * New()
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...
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
See the vtkCell API for descriptions of these methods.
int GetCellType() override
See the vtkCell API for descriptions of these methods.
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[24])
Given parametric coordinates compute inverse Jacobian transformation matrix.
vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
list of point or cell ids
Definition vtkIdList.h:133
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
cell represents a 1D line
Definition vtkLine.h:132
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:139
a cell that represents a 2D quadrilateral
Definition vtkQuad.h:87
@ VTK_HEXAHEDRON
Definition vtkCellType.h:68
int vtkIdType
Definition vtkType.h:315
#define VTK_SIZEHINT(...)