VTK  9.2.6
vtkTetra.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkTetra.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=========================================================================*/
29
30#ifndef vtkTetra_h
31#define vtkTetra_h
32
33#include "vtkCell3D.h"
34#include "vtkCommonDataModelModule.h" // For export macro
35
36class vtkLine;
37class vtkTriangle;
40
41class VTKCOMMONDATAMODEL_EXPORT vtkTetra : public vtkCell3D
42{
43public:
44 static vtkTetra* New();
45 vtkTypeMacro(vtkTetra, vtkCell3D);
46 void PrintSelf(ostream& os, vtkIndent indent) override;
47
49
52 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
53 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
54 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
55 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
56 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
57 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
58 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
59 bool GetCentroid(double centroid[3]) const override;
60 bool IsInsideOut() override;
62
66 static constexpr vtkIdType NumberOfPoints = 4;
67
71 static constexpr vtkIdType NumberOfEdges = 6;
72
76 static constexpr vtkIdType NumberOfFaces = 4;
77
82 static constexpr vtkIdType MaximumFaceSize = 3;
83
89 static constexpr vtkIdType MaximumValence = 3;
90
92
95 int GetCellType() override { return VTK_TETRA; }
96 int GetNumberOfEdges() override { return 6; }
97 int GetNumberOfFaces() override { return 4; }
98 vtkCell* GetEdge(int edgeId) override;
99 vtkCell* GetFace(int faceId) 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 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
104 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
105 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
106 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
107 double& dist2, double weights[]) override;
108 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
109 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
110 double pcoords[3], int& subId) override;
111 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
113 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
114 double* GetParametricCoords() override;
116
124 static int* GetTriangleCases(int caseId);
125
131 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
132
136 int GetParametricCenter(double pcoords[3]) override;
137
142 double GetParametricDistance(const double pcoords[3]) override;
143
147 static void TetraCenter(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
148
154 static double Circumsphere(
155 double x1[3], double x2[3], double x3[3], double x4[3], double center[3]);
156
162 static double Insphere(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
163
177 double x[3], double x1[3], double x2[3], double x3[3], double x4[3], double bcoords[4]);
178
183 static double ComputeVolume(double p1[3], double p2[3], double p3[3], double p4[3]);
184
190 int JacobianInverse(double** inverse, double derivs[12]);
191
192 static void InterpolationFunctions(const double pcoords[3], double weights[4]);
193 static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
195
199 void InterpolateFunctions(const double pcoords[3], double weights[4]) override
200 {
201 vtkTetra::InterpolationFunctions(pcoords, weights);
202 }
203 void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
204 {
205 vtkTetra::InterpolationDerivs(pcoords, derivs);
206 }
207
208
210
221
226
231
236
241
246
250 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
251
252protected:
254 ~vtkTetra() override;
255
258
259private:
260 vtkTetra(const vtkTetra&) = delete;
261 void operator=(const vtkTetra&) = delete;
262};
263
264inline int vtkTetra::GetParametricCenter(double pcoords[3])
265{
266 pcoords[0] = pcoords[1] = pcoords[2] = 0.25;
267 return 0;
268}
269
270#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:36
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
list of point or cell ids
Definition vtkIdList.h:31
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:34
cell represents a 1D line
Definition vtkLine.h:31
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:34
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
See the vtkCell API for descriptions of these methods.
bool GetCentroid(double centroid[3]) const override
See vtkCell3D API for description of these methods.
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
static const vtkIdType * GetPointToIncidentEdgesArray(vtkIdType pointId)
Static method version of GetPointToIncidentEdgesArray.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
See the vtkCell API for descriptions of these methods.
static void InterpolationDerivs(const double pcoords[3], double derivs[12])
static void TetraCenter(double p1[3], double p2[3], double p3[3], double p4[3], double center[3])
Compute the center of the tetrahedron,.
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
See the vtkCell API for descriptions of these methods.
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
See vtkCell3D API for description of these methods.
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
int GetParametricCenter(double pcoords[3]) override
Return the center of the tetrahedron in parametric coordinates.
Definition vtkTetra.h:264
static int * GetTriangleCases(int caseId)
Return the case table for table-based isocontouring (aka marching cubes style implementations).
static constexpr vtkIdType MaximumFaceSize
static contexpr handle on the maximum face size.
Definition vtkTetra.h:82
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkTriangle * Triangle
Definition vtkTetra.h:257
double * GetParametricCoords() override
See the vtkCell API for descriptions of these methods.
static constexpr vtkIdType MaximumValence
static constexpr handle on the maximum valence of this cell.
Definition vtkTetra.h:89
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
bool IsInsideOut() override
See vtkCell3D API for description of these methods.
~vtkTetra() override
void InterpolateFunctions(const double pcoords[3], double weights[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition vtkTetra.h:199
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
static int BarycentricCoords(double x[3], double x1[3], double x2[3], double x3[3], double x4[3], double bcoords[4])
Given a 3D point x[3], determine the barycentric coordinates of the point.
static constexpr vtkIdType NumberOfEdges
static contexpr handle on the number of faces.
Definition vtkTetra.h:71
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static double ComputeVolume(double p1[3], double p2[3], double p3[3], double p4[3])
Compute the volume of a tetrahedron defined by the four points p1, p2, p3, and p4.
static bool ComputeCentroid(vtkPoints *points, const vtkIdType *pointIds, double centroid[3])
Static method version of GetCentroid.
int JacobianInverse(double **inverse, double derivs[12])
Given parametric coordinates compute inverse Jacobian transformation matrix.
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition vtkTetra.h:96
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
static double Insphere(double p1[3], double p2[3], double p3[3], double p4[3], double center[3])
Compute the center (center[3]) and radius (method return value) of a sphere that just fits inside the...
void GetEdgePoints(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
See the vtkCell API for descriptions of these methods.
vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
See the vtkCell API for descriptions of these methods.
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition vtkTetra.h:97
vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
static vtkTetra * New()
static double Circumsphere(double x1[3], double x2[3], double x3[3], double x4[3], double center[3])
Compute the circumcenter (center[3]) and radius squared (method return value) of a tetrahedron define...
static constexpr vtkIdType NumberOfFaces
static contexpr handle on the number of edges.
Definition vtkTetra.h:76
vtkLine * Line
Definition vtkTetra.h:256
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
See the vtkCell API for descriptions of these methods.
vtkCell * GetFace(int faceId) override
See the vtkCell API for descriptions of these methods.
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Returns the set of points that are on the boundary of the tetrahedron that are closest parametrically...
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition vtkTetra.h:95
void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition vtkTetra.h:203
static constexpr vtkIdType NumberOfPoints
static constexpr handle on the number of points.
Definition vtkTetra.h:66
static void InterpolationFunctions(const double pcoords[3], double weights[4])
a cell that represents a triangle
Definition vtkTriangle.h:36
dataset represents arbitrary combinations of all possible cell types
@ VTK_TETRA
Definition vtkCellType.h:56
#define vtkDataArray
int vtkIdType
Definition vtkType.h:332
#define VTK_SIZEHINT(...)