VTK  9.2.6
vtkTriQuadraticHexahedron.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkTriQuadraticHexahedron.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=========================================================================*/
74#ifndef vtkTriQuadraticHexahedron_h
75#define vtkTriQuadraticHexahedron_h
76
77#include "vtkCommonDataModelModule.h" // For export macro
78#include "vtkNonLinearCell.h"
79
82class vtkHexahedron;
83class vtkDoubleArray;
84
85class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticHexahedron : public vtkNonLinearCell
86{
87public:
90 void PrintSelf(ostream& os, vtkIndent indent) override;
91
93
97 int GetCellType() override { return VTK_TRIQUADRATIC_HEXAHEDRON; }
98 int GetCellDimension() override { return 3; }
99 int GetNumberOfEdges() override { return 12; }
100 int GetNumberOfFaces() override { return 6; }
101 vtkCell* GetEdge(int) override;
102 vtkCell* GetFace(int) override;
104
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;
109 int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
110 double& dist2, double* weights) override;
111 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
112 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
114 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
115 double* GetParametricCoords() override;
116
122 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
123 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
124 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
125
130 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
131 double pcoords[3], int& subId) override;
132
133 static void InterpolationFunctions(const double pcoords[3], double weights[27]);
134 static void InterpolationDerivs(const double pcoords[3], double derivs[81]);
136
140 void InterpolateFunctions(const double pcoords[3], double weights[27]) override
141 {
143 }
144 void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
145 {
147 }
150
157 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
158 static const vtkIdType* GetFaceArray(vtkIdType faceId);
160
166 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[81]);
167
168protected:
171
176
177private:
179 void operator=(const vtkTriQuadraticHexahedron&) = delete;
180};
181
182#endif
cell represents a parabolic, 9-node isoparametric quad
object to represent cell connectivity
Definition: vtkCellArray.h:187
represent and manipulate cell attribute data
Definition: vtkCellData.h:42
abstract class to specify cell behavior
Definition: vtkCell.h:61
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:56
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:45
list of point or cell ids
Definition: vtkIdList.h:34
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:40
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:42
represent and manipulate 3D points
Definition: vtkPoints.h:40
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’).
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
static void InterpolationFunctions(const double pcoords[3], double weights[27])
vtkCell * GetEdge(int) override
Implement the vtkCell API.
vtkCell * GetFace(int) override
Implement the vtkCell API.
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.
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.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static vtkTriQuadraticHexahedron * 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...
void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetNumberOfFaces() override
Implement the vtkCell API.
void InterpolateFunctions(const double pcoords[3], double weights[27]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
~vtkTriQuadraticHexahedron() override
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:74
int vtkIdType
Definition: vtkType.h:332