VTK  9.1.0
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=========================================================================*/
99#ifndef vtkTetra_h
100#define vtkTetra_h
101
102#include "vtkCell3D.h"
103#include "vtkCommonDataModelModule.h" // For export macro
104
105class vtkLine;
106class vtkTriangle;
109
110class VTKCOMMONDATAMODEL_EXPORT vtkTetra : public vtkCell3D
111{
112public:
113 static vtkTetra* New();
114 vtkTypeMacro(vtkTetra, vtkCell3D);
115 void PrintSelf(ostream& os, vtkIndent indent) override;
116
118
121 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
122 // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
123 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkTetra::GetEdgePoints(vtkIdType, const vtkIdType*&)")
124 void GetEdgePoints(int edgeId, int*& pts) override;
125 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
126 // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
127 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkTetra::GetFacePoints(vtkIdType, const vtkIdType*&)")
128 void GetFacePoints(int faceId, int*& pts) override;
129 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
130 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
131 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
132 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
133 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
134 bool GetCentroid(double centroid[3]) const override;
135 bool IsInsideOut() override;
137
141 static constexpr vtkIdType NumberOfPoints = 4;
142
146 static constexpr vtkIdType NumberOfEdges = 6;
147
151 static constexpr vtkIdType NumberOfFaces = 4;
152
157 static constexpr vtkIdType MaximumFaceSize = 3;
158
164 static constexpr vtkIdType MaximumValence = 3;
165
167
170 int GetCellType() override { return VTK_TETRA; }
171 int GetNumberOfEdges() override { return 6; }
172 int GetNumberOfFaces() override { return 4; }
173 vtkCell* GetEdge(int edgeId) override;
174 vtkCell* GetFace(int faceId) override;
175 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
176 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
177 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
178 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
179 vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
180 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
181 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
182 double& dist2, double weights[]) override;
183 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
184 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
185 double pcoords[3], int& subId) override;
186 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
188 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
189 double* GetParametricCoords() override;
191
199 static int* GetTriangleCases(int caseId);
200
206 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
207
211 int GetParametricCenter(double pcoords[3]) override;
212
217 double GetParametricDistance(const double pcoords[3]) override;
218
222 static void TetraCenter(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
223
229 static double Circumsphere(
230 double x1[3], double x2[3], double x3[3], double x4[3], double center[3]);
231
237 static double Insphere(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
238
252 double x[3], double x1[3], double x2[3], double x3[3], double x4[3], double bcoords[4]);
253
258 static double ComputeVolume(double p1[3], double p2[3], double p3[3], double p4[3]);
259
265 int JacobianInverse(double** inverse, double derivs[12]);
266
267 static void InterpolationFunctions(const double pcoords[3], double weights[4]);
268 static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
270
274 void InterpolateFunctions(const double pcoords[3], double weights[4]) override
275 {
276 vtkTetra::InterpolationFunctions(pcoords, weights);
277 }
278 void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
279 {
280 vtkTetra::InterpolationDerivs(pcoords, derivs);
281 }
283
285
296
301
306
311
316
321
325 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
326
327protected:
329 ~vtkTetra() override;
330
333
334private:
335 vtkTetra(const vtkTetra&) = delete;
336 void operator=(const vtkTetra&) = delete;
337};
338
339inline int vtkTetra::GetParametricCenter(double pcoords[3])
340{
341 pcoords[0] = pcoords[1] = pcoords[2] = 0.25;
342 return 0;
343}
344
345#endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:40
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
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
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
cell represents a 1D line
Definition: vtkLine.h:140
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:111
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.
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
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.
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:339
static int * GetTriangleCases(int caseId)
Return the case table for table-based isocontouring (aka marching cubes style implementations).
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkTriangle * Triangle
Definition: vtkTetra.h:332
double * GetParametricCoords() override
See the vtkCell API for descriptions of these methods.
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.
~vtkTetra() override
void InterpolateFunctions(const double pcoords[3], double weights[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkTetra.h:274
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 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:171
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.
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.
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition: vtkTetra.h:172
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...
vtkLine * Line
Definition: vtkTetra.h:331
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...
void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkTetra.h:278
static void InterpolationFunctions(const double pcoords[3], double weights[4])
a cell that represents a triangle
Definition: vtkTriangle.h:145
dataset represents arbitrary combinations of all possible cell types
int GetCellType(const Ioss::ElementTopology *topology)
Returns VTK celltype for a Ioss topology element.
@ points
Definition: vtkX3D.h:452
@ value
Definition: vtkX3D.h:226
@ center
Definition: vtkX3D.h:236
@ index
Definition: vtkX3D.h:252
@ VTK_TETRA
Definition: vtkCellType.h:95
#define VTK_DEPRECATED_IN_9_0_0(reason)
int vtkIdType
Definition: vtkType.h:332
#define VTK_SIZEHINT(...)