forked from Kitware/VTK
-
Notifications
You must be signed in to change notification settings - Fork 0
/
vtkQuadraticHexahedron.h
160 lines (138 loc) · 5.77 KB
/
vtkQuadraticHexahedron.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
/*=========================================================================
Program: Visualization Toolkit
Module: vtkQuadraticHexahedron.h
Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
All rights reserved.
See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
This software is distributed WITHOUT ANY WARRANTY; without even
the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
PURPOSE. See the above copyright notice for more information.
=========================================================================*/
/**
* @class vtkQuadraticHexahedron
* @brief cell represents a parabolic, 20-node isoparametric hexahedron
*
* vtkQuadraticHexahedron is a concrete implementation of vtkNonLinearCell to
* represent a three-dimensional, 20-node isoparametric parabolic
* hexahedron. The interpolation is the standard finite element, quadratic
* isoparametric shape function. The cell includes a mid-edge node. The
* ordering of the twenty points defining the cell is point ids (0-7,8-19)
* where point ids 0-7 are the eight corner vertices of the cube; followed by
* twelve midedge nodes (8-19). Note that these midedge nodes correspond lie
* on the edges defined by (0,1), (1,2), (2,3), (3,0), (4,5), (5,6), (6,7),
* (7,4), (0,4), (1,5), (2,6), (3,7).
*
* @sa
* vtkQuadraticEdge vtkQuadraticTriangle vtkQuadraticTetra
* vtkQuadraticQuad vtkQuadraticPyramid vtkQuadraticWedge
*/
#ifndef vtkQuadraticHexahedron_h
#define vtkQuadraticHexahedron_h
#include "vtkCommonDataModelModule.h" // For export macro
#include "vtkNonLinearCell.h"
class vtkQuadraticEdge;
class vtkQuadraticQuad;
class vtkHexahedron;
class vtkDoubleArray;
class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticHexahedron : public vtkNonLinearCell
{
public:
static vtkQuadraticHexahedron *New();
vtkTypeMacro(vtkQuadraticHexahedron,vtkNonLinearCell);
void PrintSelf(ostream& os, vtkIndent indent) override;
//@{
/**
* Implement the vtkCell API. See the vtkCell API for descriptions
* of these methods.
*/
int GetCellType() override {return VTK_QUADRATIC_HEXAHEDRON;}
int GetCellDimension() override {return 3;}
int GetNumberOfEdges() override {return 12;}
int GetNumberOfFaces() override {return 6;}
vtkCell *GetEdge(int) override;
vtkCell *GetFace(int) override;
//@}
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
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;
int EvaluatePosition(const double x[3], double closestPoint[3],
int& subId, double pcoords[3],
double& dist2, double weights[]) override;
void EvaluateLocation(int& subId, const double pcoords[3], double x[3],
double *weights) override;
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
void Derivatives(int subId, const double pcoords[3], const double *values,
int dim, double *derivs) override;
double *GetParametricCoords() override;
/**
* Clip this quadratic hexahedron using scalar value provided. Like
* contouring, except that it cuts the hex to produce linear
* tetrahedron.
*/
void Clip(double value, vtkDataArray *cellScalars,
vtkIncrementalPointLocator *locator, vtkCellArray *tetras,
vtkPointData *inPd, vtkPointData *outPd,
vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
int insideOut) override;
/**
* Line-edge intersection. Intersection has to occur within [0,1] parametric
* coordinates and with specified tolerance.
*/
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
double x[3], double pcoords[3], int& subId) override;
/**
* @deprecated Replaced by vtkQuadraticHexahedron::InterpolateFunctions as of VTK 5.2
*/
static void InterpolationFunctions(const double pcoords[3], double weights[20]);
/**
* @deprecated Replaced by vtkQuadraticHexahedron::InterpolateDerivs as of VTK 5.2
*/
static void InterpolationDerivs(const double pcoords[3], double derivs[60]);
//@{
/**
* Compute the interpolation functions/derivatives
* (aka shape functions/derivatives)
*/
void InterpolateFunctions(const double pcoords[3], double weights[20]) override
{
vtkQuadraticHexahedron::InterpolationFunctions(pcoords,weights);
}
void InterpolateDerivs(const double pcoords[3], double derivs[60]) override
{
vtkQuadraticHexahedron::InterpolationDerivs(pcoords,derivs);
}
//@}
//@{
/**
* Return the ids of the vertices defining edge/face (`edgeId`/`faceId').
* Ids are related to the cell, not to the dataset.
*/
static int *GetEdgeArray(int edgeId);
static int *GetFaceArray(int faceId);
//@}
/**
* Given parametric coordinates compute inverse Jacobian transformation
* matrix. Returns 9 elements of 3x3 inverse Jacobian plus interpolation
* function derivatives.
*/
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[60]);
protected:
vtkQuadraticHexahedron();
~vtkQuadraticHexahedron() override;
vtkQuadraticEdge *Edge;
vtkQuadraticQuad *Face;
vtkHexahedron *Hex;
vtkPointData *PointData;
vtkCellData *CellData;
vtkDoubleArray *CellScalars;
vtkDoubleArray *Scalars;
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId,
vtkDataArray *cellScalars);
private:
vtkQuadraticHexahedron(const vtkQuadraticHexahedron&) = delete;
void operator=(const vtkQuadraticHexahedron&) = delete;
};
#endif