VTK  9.1.0
vtkHexahedron.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkHexahedron.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
31#ifndef vtkHexahedron_h
32#define vtkHexahedron_h
33
34#include "vtkCell3D.h"
35#include "vtkCommonDataModelModule.h" // For export macro
36
37class vtkLine;
38class vtkQuad;
40
41class VTKCOMMONDATAMODEL_EXPORT vtkHexahedron : public vtkCell3D
42{
43public:
44 static vtkHexahedron* New();
45 vtkTypeMacro(vtkHexahedron, vtkCell3D);
46 void PrintSelf(ostream& os, vtkIndent indent) override;
47
49
52 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
53 // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
54 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkHexahedron::GetEdgePoints(vtkIdType, const vtkIdType*&)")
55 void GetEdgePoints(int edgeId, int*& pts) override;
56 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
57 // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
58 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkHexahedron::GetFacePoints(vtkIdType, const vtkIdType*&)")
59 void GetFacePoints(int faceId, int*& pts) override;
60 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
61 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
62 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
63 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
64 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
65 bool GetCentroid(double centroid[3]) const override;
67
71 static constexpr vtkIdType NumberOfPoints = 8;
72
76 static constexpr vtkIdType NumberOfEdges = 12;
77
81 static constexpr vtkIdType NumberOfFaces = 6;
82
87 static constexpr vtkIdType MaximumFaceSize = 4;
88
94 static constexpr vtkIdType MaximumValence = 3;
95
97
100 int GetCellType() override { return VTK_HEXAHEDRON; }
101 int GetNumberOfEdges() override { return 12; }
102 int GetNumberOfFaces() override { return 6; }
103 vtkCell* GetEdge(int edgeId) override;
104 vtkCell* GetFace(int faceId) override;
105 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
106 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
107 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
108 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
110
111 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
112 double& dist2, double weights[]) override;
113 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
114 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
115 double pcoords[3], int& subId) override;
116 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
118 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
119 double* GetParametricCoords() override;
120
128 static int* GetTriangleCases(int caseId);
129
130 static void InterpolationFunctions(const double pcoords[3], double weights[8]);
131 static void InterpolationDerivs(const double pcoords[3], double derivs[24]);
133
137 void InterpolateFunctions(const double pcoords[3], double weights[8]) override
138 {
140 }
141 void InterpolateDerivs(const double pcoords[3], double derivs[24]) override
142 {
143 vtkHexahedron::InterpolationDerivs(pcoords, derivs);
144 }
146
148
159
164
169
174
179
184
188 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
189
195 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[24]);
196
197protected:
199 ~vtkHexahedron() override;
200
203
204private:
205 vtkHexahedron(const vtkHexahedron&) = delete;
206 void operator=(const vtkHexahedron&) = delete;
207};
208
209#endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:40
object to represent cell connectivity
Definition: vtkCellArray.h:181
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:58
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:42
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
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)
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
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.
~vtkHexahedron() override
static void InterpolationDerivs(const double pcoords[3], double derivs[24])
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
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.
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
vtkLine * Line
vtkQuad * Quad
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.
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[24])
Given parametric coordinates compute inverse Jacobian transformation matrix.
list of point or cell ids
Definition: vtkIdList.h:31
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:34
cell represents a 1D line
Definition: vtkLine.h:31
represent and manipulate point attribute data
Definition: vtkPointData.h:33
represent and manipulate 3D points
Definition: vtkPoints.h:34
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:36
@ VTK_HEXAHEDRON
Definition: vtkCellType.h:58
#define VTK_DEPRECATED_IN_9_0_0(reason)
int vtkIdType
Definition: vtkType.h:332
#define VTK_SIZEHINT(...)