VTK  9.1.0
vtkQuadraticHexahedron.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticHexahedron.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 =========================================================================*/
34 #ifndef vtkQuadraticHexahedron_h
35 #define vtkQuadraticHexahedron_h
36 
37 #include "vtkCommonDataModelModule.h" // For export macro
38 #include "vtkNonLinearCell.h"
39 
40 class vtkQuadraticEdge;
41 class vtkQuadraticQuad;
42 class vtkHexahedron;
43 class vtkDoubleArray;
44 
45 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticHexahedron : public vtkNonLinearCell
46 {
47 public:
50  void PrintSelf(ostream& os, vtkIndent indent) override;
51 
53 
57  int GetCellType() override { return VTK_QUADRATIC_HEXAHEDRON; }
58  int GetCellDimension() override { return 3; }
59  int GetNumberOfEdges() override { return 12; }
60  int GetNumberOfFaces() override { return 6; }
61  vtkCell* GetEdge(int) override;
62  vtkCell* GetFace(int) override;
64 
65  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
66  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
67  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
68  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
69  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
70  double& dist2, double weights[]) override;
71  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
72  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
74  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
75  double* GetParametricCoords() override;
76 
82  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
83  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
84  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
85 
90  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
91  double pcoords[3], int& subId) override;
92 
93  static void InterpolationFunctions(const double pcoords[3], double weights[20]);
94  static void InterpolationDerivs(const double pcoords[3], double derivs[60]);
96 
100  void InterpolateFunctions(const double pcoords[3], double weights[20]) override
101  {
103  }
104  void InterpolateDerivs(const double pcoords[3], double derivs[60]) override
105  {
107  }
110 
117  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
118  static const vtkIdType* GetFaceArray(vtkIdType faceId);
120 
126  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[60]);
127 
128 protected:
131 
139 
140  void Subdivide(
141  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
142 
143 private:
145  void operator=(const vtkQuadraticHexahedron&) = delete;
146 };
147 
148 #endif
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
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:42
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
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:33
represent and manipulate 3D points
Definition: vtkPoints.h:34
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 20-node isoparametric hexahedron
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.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[60])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
int GetNumberOfEdges() override
Implement the vtkCell API.
static void InterpolationFunctions(const double pcoords[3], double weights[20])
vtkCell * GetEdge(int) override
Implement the vtkCell API.
void InterpolateDerivs(const double pcoords[3], double derivs[60]) 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.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void InterpolateFunctions(const double pcoords[3], double weights[20]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetCellType() override
Implement the vtkCell API.
vtkCell * GetFace(int) override
Implement the vtkCell API.
int GetCellDimension() override
Implement the vtkCell API.
static vtkQuadraticHexahedron * 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...
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 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.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
static void InterpolationDerivs(const double pcoords[3], double derivs[60])
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
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 quadratic hexahedron using scalar value provided.
~vtkQuadraticHexahedron() override
int GetNumberOfFaces() override
Implement the vtkCell API.
cell represents a parabolic, 8-node isoparametric quad
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:70
int vtkIdType
Definition: vtkType.h:332