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=========================================================================*/
32
33#ifndef vtkTetra_h
34#define vtkTetra_h
35
36#include "vtkCell3D.h"
37#include "vtkCommonDataModelModule.h" // For export macro
38
39class vtkLine;
40class vtkTriangle;
43
44class VTKCOMMONDATAMODEL_EXPORT vtkTetra : public vtkCell3D
45{
46public:
47 static vtkTetra* New();
48 vtkTypeMacro(vtkTetra, vtkCell3D);
49 void PrintSelf(ostream& os, vtkIndent indent) override;
50
52
55 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
56 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
57 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
58 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
59 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
60 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
61 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
62 bool GetCentroid(double centroid[3]) const override;
63 bool IsInsideOut() override;
65
69 static constexpr vtkIdType NumberOfPoints = 4;
70
74 static constexpr vtkIdType NumberOfEdges = 6;
75
79 static constexpr vtkIdType NumberOfFaces = 4;
80
85 static constexpr vtkIdType MaximumFaceSize = 3;
86
92 static constexpr vtkIdType MaximumValence = 3;
93
95
98 int GetCellType() override { return VTK_TETRA; }
99 int GetNumberOfEdges() override { return 6; }
100 int GetNumberOfFaces() override { return 4; }
101 vtkCell* GetEdge(int edgeId) override;
102 vtkCell* GetFace(int faceId) override;
103 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
104 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
105 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
106 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
107 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
108 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
109 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
110 double& dist2, double weights[]) override;
111 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
112 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
113 double pcoords[3], int& subId) 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;
119
127 static int* GetTriangleCases(int caseId);
128
134 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
135
139 int GetParametricCenter(double pcoords[3]) override;
140
145 double GetParametricDistance(const double pcoords[3]) override;
146
150 static void TetraCenter(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
151
157 static double Circumsphere(
158 double x1[3], double x2[3], double x3[3], double x4[3], double center[3]);
159
165 static double Insphere(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
166
180 double x[3], double x1[3], double x2[3], double x3[3], double x4[3], double bcoords[4]);
181
186 static double ComputeVolume(double p1[3], double p2[3], double p3[3], double p4[3]);
187
193 int JacobianInverse(double** inverse, double derivs[12]);
194
195 static void InterpolationFunctions(const double pcoords[3], double weights[4]);
196 static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
198
202 void InterpolateFunctions(const double pcoords[3], double weights[4]) override
203 {
204 vtkTetra::InterpolationFunctions(pcoords, weights);
205 }
206 void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
207 {
208 vtkTetra::InterpolationDerivs(pcoords, derivs);
209 }
210
211
213
224
229
234
239
244
249
253 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
254
255protected:
257 ~vtkTetra() override;
258
261
262private:
263 vtkTetra(const vtkTetra&) = delete;
264 void operator=(const vtkTetra&) = delete;
265};
266
267inline int vtkTetra::GetParametricCenter(double pcoords[3])
268{
269 pcoords[0] = pcoords[1] = pcoords[2] = 0.25;
270 return 0;
271}
272
273#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:42
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
list of point or cell ids
Definition vtkIdList.h:34
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:40
cell represents a 1D line
Definition vtkLine.h:34
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:40
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:267
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:85
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkTriangle * Triangle
Definition vtkTetra.h:260
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:92
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:202
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:74
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:99
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:100
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:79
vtkLine * Line
Definition vtkTetra.h:259
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:98
void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition vtkTetra.h:206
static constexpr vtkIdType NumberOfPoints
static constexpr handle on the number of points.
Definition vtkTetra.h:69
static void InterpolationFunctions(const double pcoords[3], double weights[4])
a cell that represents a triangle
Definition vtkTriangle.h:39
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(...)