VTK  9.6.20260402
vtkTetra.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
99
100#ifndef vtkTetra_h
101#define vtkTetra_h
102
103#include "vtkCell3D.h"
104#include "vtkCommonDataModelModule.h" // For export macro
105
106VTK_ABI_NAMESPACE_BEGIN
107class vtkLine;
108class vtkTriangle;
111
112class VTKCOMMONDATAMODEL_EXPORT vtkTetra : public vtkCell3D
113{
114public:
115 static vtkTetra* New();
116 vtkTypeMacro(vtkTetra, vtkCell3D);
117 void PrintSelf(ostream& os, vtkIndent indent) override;
118
120
123 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
124 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
125 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
126 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
127 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
128 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
129 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
130 bool GetCentroid(double centroid[3]) const override;
131 bool IsInsideOut() override;
133
137 static constexpr vtkIdType NumberOfPoints = 4;
138
142 static constexpr vtkIdType NumberOfEdges = 6;
143
147 static constexpr vtkIdType NumberOfFaces = 4;
148
153 static constexpr vtkIdType MaximumFaceSize = 3;
154
160 static constexpr vtkIdType MaximumValence = 3;
161
163
166 int GetCellType() override { return VTK_TETRA; }
167 int GetNumberOfEdges() override { return 6; }
168 int GetNumberOfFaces() override { return 4; }
169 vtkCell* GetEdge(int edgeId) override;
170 vtkCell* GetFace(int faceId) override;
171 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
172 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
173 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
174 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
175 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
176 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
177 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
178 double& dist2, double weights[]) override;
179 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
180 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
181 double pcoords[3], int& subId) override;
182 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
184 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
185 double* GetParametricCoords() override;
187
195 VTK_DEPRECATED_IN_9_7_0("Use vtkMarchingCellsContourCases::GetTetraCase instead.")
196 static int* GetTriangleCases(int caseId);
197
203 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
204
208 int GetParametricCenter(double pcoords[3]) override;
209
214 double GetParametricDistance(const double pcoords[3]) override;
215
219 static void TetraCenter(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
220
226 static double Circumsphere(
227 double x1[3], double x2[3], double x3[3], double x4[3], double center[3]);
228
234 static double Insphere(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
235
249 double x[3], double x1[3], double x2[3], double x3[3], double x4[3], double bcoords[4]);
250
255 static double ComputeVolume(double p1[3], double p2[3], double p3[3], double p4[3]);
256
262 int JacobianInverse(double** inverse, double derivs[12]);
263
264 static void InterpolationFunctions(const double pcoords[3], double weights[4]);
265 static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
267
271 void InterpolateFunctions(const double pcoords[3], double weights[4]) override
272 {
273 vtkTetra::InterpolationFunctions(pcoords, weights);
274 }
275 void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
276 {
277 vtkTetra::InterpolationDerivs(pcoords, derivs);
278 }
279
280
282
293
298
303
308
313
318
322 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
323
324protected:
326 ~vtkTetra() override = default;
327
330
331private:
332 vtkTetra(const vtkTetra&) = delete;
333 void operator=(const vtkTetra&) = delete;
334};
335
336inline int vtkTetra::GetParametricCenter(double pcoords[3])
337{
338 pcoords[0] = pcoords[1] = pcoords[2] = 0.25;
339 return 0;
340}
341
342VTK_ABI_NAMESPACE_END
343#endif
object to represent cell connectivity
represent and manipulate cell attribute data
list of point or cell ids
Definition vtkIdList.h:135
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
cell represents a 1D line
Definition vtkLine.h:132
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:140
Hold a reference to a vtkObjectBase instance.
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:336
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:153
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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:160
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=default
void InterpolateFunctions(const double pcoords[3], double weights[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition vtkTetra.h:271
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
vtkSmartPointer< vtkTriangle > Triangle
Definition vtkTetra.h:329
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:142
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:167
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.
vtkSmartPointer< vtkLine > Line
Definition vtkTetra.h:328
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:168
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:147
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
See the vtkCell API for descriptions of these methods.
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:166
void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition vtkTetra.h:275
static constexpr vtkIdType NumberOfPoints
static constexpr handle on the number of points.
Definition vtkTetra.h:137
static void InterpolationFunctions(const double pcoords[3], double weights[4])
a cell that represents a triangle
dataset represents arbitrary combinations of all possible cell types
@ VTK_TETRA
Definition vtkCellType.h:47
#define vtkDataArray
#define VTK_DEPRECATED_IN_9_7_0(reason)
int vtkIdType
Definition vtkType.h:363
#define VTK_SIZEHINT(...)