VTK  9.3.0
vtkQuadraticTetra.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
37#ifndef vtkQuadraticTetra_h
38#define vtkQuadraticTetra_h
39
40#include "vtkCommonDataModelModule.h" // For export macro
41#include "vtkNonLinearCell.h"
42
43VTK_ABI_NAMESPACE_BEGIN
46class vtkTetra;
47class vtkDoubleArray;
48
49class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticTetra : public vtkNonLinearCell
50{
51public:
54 void PrintSelf(ostream& os, vtkIndent indent) override;
55
57
61 int GetCellType() override { return VTK_QUADRATIC_TETRA; }
62 int GetCellDimension() override { return 3; }
63 int GetNumberOfEdges() override { return 6; }
64 int GetNumberOfFaces() override { return 4; }
65 vtkCell* GetEdge(int) override;
66 vtkCell* GetFace(int) override;
68
69 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
70 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
71 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
72 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
73 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
74 double& dist2, double weights[]) override;
75 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
76 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
78 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
79 double* GetParametricCoords() override;
80
85 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
86 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
87 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
88
97 bool StableClip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
98 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
99 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
100
105 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
106 double pcoords[3], int& subId) override;
107
111 int GetParametricCenter(double pcoords[3]) override;
112
117 double GetParametricDistance(const double pcoords[3]) override;
118
119 static void InterpolationFunctions(const double pcoords[3], double weights[10]);
120 static void InterpolationDerivs(const double pcoords[3], double derivs[30]);
122
126 void InterpolateFunctions(const double pcoords[3], double weights[10]) override
127 {
129 }
130 void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
131 {
133 }
136
143 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
144 static const vtkIdType* GetFaceArray(vtkIdType faceId);
146
152 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[30]);
153
154protected:
157
161 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
162
163private:
164 vtkQuadraticTetra(const vtkQuadraticTetra&) = delete;
165 void operator=(const vtkQuadraticTetra&) = delete;
166};
167
168VTK_ABI_NAMESPACE_END
169#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:40
abstract class to specify cell behavior
Definition vtkCell.h:59
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:38
abstract superclass for non-linear cells
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:38
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 10-node isoparametric tetrahedron
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic tetra in parametric coordinates.
int GetNumberOfEdges() override
Implement the vtkCell API.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
int GetNumberOfFaces() override
Implement the vtkCell API.
vtkQuadraticEdge * Edge
static void InterpolationDerivs(const double pcoords[3], double derivs[30])
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
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.
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[30])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
void InterpolateFunctions(const double pcoords[3], double weights[10]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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 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 const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
~vtkQuadraticTetra() override
static vtkQuadraticTetra * New()
vtkQuadraticTriangle * Face
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.
int GetCellDimension() override
Implement the vtkCell API.
static void InterpolationFunctions(const double pcoords[3], double weights[10])
vtkCell * GetFace(int) override
Implement the vtkCell API.
int GetCellType() 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 edge using scalar value provided.
bool StableClip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tetras, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this edge using scalar value provided.
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’).
vtkDoubleArray * Scalars
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.
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:43
@ VTK_QUADRATIC_TETRA
Definition vtkCellType.h:60
int vtkIdType
Definition vtkType.h:315