VTK  9.1.0
vtkBiQuadraticQuadraticHexahedron.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkBiQuadraticQuadraticHexahedron.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=========================================================================*/
87#ifndef vtkBiQuadraticQuadraticHexahedron_h
88#define vtkBiQuadraticQuadraticHexahedron_h
89
90#include "vtkCommonDataModelModule.h" // For export macro
91#include "vtkNonLinearCell.h"
92
96class vtkHexahedron;
97class vtkDoubleArray;
98
99class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticHexahedron : public vtkNonLinearCell
100{
101public:
104 void PrintSelf(ostream& os, vtkIndent indent) override;
105
107
112 int GetCellDimension() override { return 3; }
113 int GetNumberOfEdges() override { return 12; }
114 int GetNumberOfFaces() override { return 6; }
115 vtkCell* GetEdge(int) override;
116 vtkCell* GetFace(int) override;
118
119 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
120 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
121 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
122 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
123 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
124 double& dist2, double weights[]) override;
125 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
126 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
128 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
129 double* GetParametricCoords() override;
130
136 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
137 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
138 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
139
144 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
145 double pcoords[3], int& subId) override;
146
147 static void InterpolationFunctions(const double pcoords[3], double weights[24]);
148 static void InterpolationDerivs(const double pcoords[3], double derivs[72]);
150
154 void InterpolateFunctions(const double pcoords[3], double weights[24]) override
155 {
157 }
158 void InterpolateDerivs(const double pcoords[3], double derivs[72]) override
159 {
161 }
164
171 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
172 static const vtkIdType* GetFaceArray(vtkIdType faceId);
174
180 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[72]);
181
182protected:
185
194
196 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
197
198private:
200 void operator=(const vtkBiQuadraticQuadraticHexahedron&) = delete;
201};
202
203#endif
cell represents a parabolic, 9-node isoparametric quad
cell represents a biquadratic, 24-node isoparametric hexahedron
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[72])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void InterpolateDerivs(const double pcoords[3], double derivs[72]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetCellType() override
Implement the vtkCell API.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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...
static vtkBiQuadraticQuadraticHexahedron * New()
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
int GetNumberOfFaces() override
Implement the vtkCell API.
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[24]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
static void InterpolationFunctions(const double pcoords[3], double weights[24])
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...
vtkCell * GetEdge(int) override
Implement the vtkCell API.
int GetCellDimension() override
Implement the vtkCell API.
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 GetNumberOfEdges() 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 void InterpolationDerivs(const double pcoords[3], double derivs[72])
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkCell * GetFace(int) override
Implement the vtkCell API.
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.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
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 biquadratic hexahedron using scalar value provided.
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
list of point or cell ids
Definition: vtkIdList.h:140
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:113
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 8-node isoparametric quad
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:118
int vtkIdType
Definition: vtkType.h:332