VTK  9.3.20240418
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
67 #ifndef vtkQuadraticTetra_h
68 #define vtkQuadraticTetra_h
69 
70 #include "vtkCommonDataModelModule.h" // For export macro
71 #include "vtkNonLinearCell.h"
72 
73 VTK_ABI_NAMESPACE_BEGIN
74 class vtkQuadraticEdge;
76 class vtkTetra;
77 class vtkDoubleArray;
78 
79 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticTetra : public vtkNonLinearCell
80 {
81 public:
84  void PrintSelf(ostream& os, vtkIndent indent) override;
85 
87 
91  int GetCellType() override { return VTK_QUADRATIC_TETRA; }
92  int GetCellDimension() override { return 3; }
93  int GetNumberOfEdges() override { return 6; }
94  int GetNumberOfFaces() override { return 4; }
95  vtkCell* GetEdge(int) override;
96  vtkCell* GetFace(int) override;
98 
99  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
100  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
101  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
102  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
103  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
104  double& dist2, double weights[]) override;
105  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
106  int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
108  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
109  double* GetParametricCoords() override;
110 
115  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
116  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
117  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
118 
127  bool StableClip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
128  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
129  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
130 
135  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
136  double pcoords[3], int& subId) override;
137 
141  int GetParametricCenter(double pcoords[3]) override;
142 
147  double GetParametricDistance(const double pcoords[3]) override;
148 
149  static void InterpolationFunctions(const double pcoords[3], double weights[10]);
150  static void InterpolationDerivs(const double pcoords[3], double derivs[30]);
152 
156  void InterpolateFunctions(const double pcoords[3], double weights[10]) override
157  {
159  }
160  void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
161  {
163  }
166 
173  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
174  static const vtkIdType* GetFaceArray(vtkIdType faceId);
176 
182  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[30]);
183 
184 protected:
186  ~vtkQuadraticTetra() override;
187 
191  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
192 
193 private:
194  vtkQuadraticTetra(const vtkQuadraticTetra&) = delete;
195  void operator=(const vtkQuadraticTetra&) = delete;
196 };
197 
198 VTK_ABI_NAMESPACE_END
199 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:286
represent and manipulate cell attribute data
Definition: vtkCellData.h:141
abstract class to specify cell behavior
Definition: vtkCell.h:130
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:155
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:133
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:108
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:140
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 GetNumberOfFaces() override
Implement the vtkCell API.
vtkQuadraticEdge * Edge
static void InterpolationDerivs(const double pcoords[3], double derivs[30])
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)
vtkCell * GetEdge(int) override
Implement the vtkCell API.
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...
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...
~vtkQuadraticTetra() override
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCell * GetFace(int) override
Implement the vtkCell API.
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 TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
int GetCellDimension() override
Implement the vtkCell API.
static vtkQuadraticTetra * New()
static void InterpolationFunctions(const double pcoords[3], double weights[10])
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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.
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:113
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_QUADRATIC_TETRA
Definition: vtkCellType.h:79
int vtkIdType
Definition: vtkType.h:315