VTK  9.1.0
vtkQuadraticTetra.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkQuadraticTetra.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=========================================================================*/
76#ifndef vtkQuadraticTetra_h
77#define vtkQuadraticTetra_h
78
79#include "vtkCommonDataModelModule.h" // For export macro
80#include "vtkNonLinearCell.h"
81
84class vtkTetra;
85class vtkDoubleArray;
86
87class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticTetra : public vtkNonLinearCell
88{
89public:
92 void PrintSelf(ostream& os, vtkIndent indent) override;
93
95
99 int GetCellType() override { return VTK_QUADRATIC_TETRA; }
100 int GetCellDimension() override { return 3; }
101 int GetNumberOfEdges() override { return 6; }
102 int GetNumberOfFaces() override { return 4; }
103 vtkCell* GetEdge(int) override;
104 vtkCell* GetFace(int) override;
106
107 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
108 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
109 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
110 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
111 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
112 double& dist2, double weights[]) override;
113 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
114 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
116 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
117 double* GetParametricCoords() override;
118
123 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
124 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
125 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
126
131 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
132 double pcoords[3], int& subId) override;
133
137 int GetParametricCenter(double pcoords[3]) override;
138
143 double GetParametricDistance(const double pcoords[3]) override;
144
145 static void InterpolationFunctions(const double pcoords[3], double weights[10]);
146 static void InterpolationDerivs(const double pcoords[3], double derivs[30]);
148
152 void InterpolateFunctions(const double pcoords[3], double weights[10]) override
153 {
155 }
156 void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
157 {
159 }
162
169 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
170 static const vtkIdType* GetFaceArray(vtkIdType faceId);
172
178 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[30]);
179
180protected:
183
187 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
188
189private:
190 vtkQuadraticTetra(const vtkQuadraticTetra&) = delete;
191 void operator=(const vtkQuadraticTetra&) = delete;
192};
193
194#endif
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
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, 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.
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:111
@ VTK_QUADRATIC_TETRA
Definition: vtkCellType.h:108
int vtkIdType
Definition: vtkType.h:332