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 
93 class vtkQuadraticEdge;
94 class vtkQuadraticQuad;
95 class vtkBiQuadraticQuad;
96 class vtkHexahedron;
97 class vtkDoubleArray;
98 
99 class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticHexahedron : public vtkNonLinearCell
100 {
101 public:
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 
182 protected:
185 
194 
195  void Subdivide(
196  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
197 
198 private:
200  void operator=(const vtkBiQuadraticQuadraticHexahedron&) = delete;
201 };
202 
203 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:143
vtkBiQuadraticQuadraticHexahedron::New
static vtkBiQuadraticQuadraticHexahedron * New()
vtkBiQuadraticQuadraticHexahedron::CellData
vtkCellData * CellData
Definition: vtkBiQuadraticQuadraticHexahedron.h:191
vtkBiQuadraticQuadraticHexahedron::GetNumberOfEdges
int GetNumberOfEdges() override
Implement the vtkCell API.
Definition: vtkBiQuadraticQuadraticHexahedron.h:113
vtkBiQuadraticQuadraticHexahedron::InterpolateFunctions
void InterpolateFunctions(const double pcoords[3], double weights[24]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkBiQuadraticQuadraticHexahedron.h:154
vtkHexahedron
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:111
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:142
vtkX3D::value
@ value
Definition: vtkX3D.h:226
vtkBiQuadraticQuadraticHexahedron::Triangulate
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
vtkIdType
int vtkIdType
Definition: vtkType.h:332
vtkBiQuadraticQuadraticHexahedron::Scalars
vtkDoubleArray * Scalars
Definition: vtkBiQuadraticQuadraticHexahedron.h:193
vtkBiQuadraticQuadraticHexahedron::PointData
vtkPointData * PointData
Definition: vtkBiQuadraticQuadraticHexahedron.h:190
vtkBiQuadraticQuadraticHexahedron::Subdivide
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
vtkBiQuadraticQuadraticHexahedron::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkBiQuadraticQuadraticHexahedron::GetEdge
vtkCell * GetEdge(int) override
Implement the vtkCell API.
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
vtkBiQuadraticQuadraticHexahedron
cell represents a biquadratic, 24-node isoparametric hexahedron
Definition: vtkBiQuadraticQuadraticHexahedron.h:100
vtkQuadraticQuad
cell represents a parabolic, 8-node isoparametric quad
Definition: vtkQuadraticQuad.h:63
vtkBiQuadraticQuadraticHexahedron::GetNumberOfFaces
int GetNumberOfFaces() override
Implement the vtkCell API.
Definition: vtkBiQuadraticQuadraticHexahedron.h:114
vtkBiQuadraticQuadraticHexahedron::IntersectWithLine
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.
vtkBiQuadraticQuadraticHexahedron::CellScalars
vtkDoubleArray * CellScalars
Definition: vtkBiQuadraticQuadraticHexahedron.h:192
vtkBiQuadraticQuad
cell represents a parabolic, 9-node isoparametric quad
Definition: vtkBiQuadraticQuad.h:71
vtkBiQuadraticQuadraticHexahedron::EvaluateLocation
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkBiQuadraticQuadraticHexahedron::Hex
vtkHexahedron * Hex
Definition: vtkBiQuadraticQuadraticHexahedron.h:189
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:147
vtkBiQuadraticQuadraticHexahedron::GetParametricCoords
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:113
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:290
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:140
VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
@ VTK_BIQUADRATIC_QUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:118
vtkBiQuadraticQuadraticHexahedron::GetCellDimension
int GetCellDimension() override
Implement the vtkCell API.
Definition: vtkBiQuadraticQuadraticHexahedron.h:112
vtkBiQuadraticQuadraticHexahedron::InterpolateDerivs
void InterpolateDerivs(const double pcoords[3], double derivs[72]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkBiQuadraticQuadraticHexahedron.h:158
vtkBiQuadraticQuadraticHexahedron::CellBoundary
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...
vtkBiQuadraticQuadraticHexahedron::JacobianInverse
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[72])
Given parametric coordinates compute inverse Jacobian transformation matrix.
vtkBiQuadraticQuadraticHexahedron::Contour
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.
vtkBiQuadraticQuadraticHexahedron::Clip
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.
vtkBiQuadraticQuadraticHexahedron::~vtkBiQuadraticQuadraticHexahedron
~vtkBiQuadraticQuadraticHexahedron() override
vtkBiQuadraticQuadraticHexahedron::Edge
vtkQuadraticEdge * Edge
Definition: vtkBiQuadraticQuadraticHexahedron.h:186
vtkNonLinearCell.h
vtkBiQuadraticQuadraticHexahedron::Derivatives
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkBiQuadraticQuadraticHexahedron::GetFace
vtkCell * GetFace(int) override
Implement the vtkCell API.
vtkBiQuadraticQuadraticHexahedron::GetEdgeArray
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkNonLinearCell
abstract superclass for non-linear cells
Definition: vtkNonLinearCell.h:37
vtkBiQuadraticQuadraticHexahedron::BiQuadFace
vtkBiQuadraticQuad * BiQuadFace
Definition: vtkBiQuadraticQuadraticHexahedron.h:188
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:145
vtkBiQuadraticQuadraticHexahedron::InterpolationDerivs
static void InterpolationDerivs(const double pcoords[3], double derivs[72])
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:60
vtkBiQuadraticQuadraticHexahedron::Face
vtkQuadraticQuad * Face
Definition: vtkBiQuadraticQuadraticHexahedron.h:187
vtkBiQuadraticQuadraticHexahedron::GetFaceArray
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkBiQuadraticQuadraticHexahedron::vtkBiQuadraticQuadraticHexahedron
vtkBiQuadraticQuadraticHexahedron()
vtkBiQuadraticQuadraticHexahedron::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkBiQuadraticQuadraticHexahedron.h:111
vtkBiQuadraticQuadraticHexahedron::EvaluatePosition
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...
vtkBiQuadraticQuadraticHexahedron::InterpolationFunctions
static void InterpolationFunctions(const double pcoords[3], double weights[24])