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 =========================================================================*/
30 #ifndef vtkTetra_h
31 #define vtkTetra_h
32 
33 #include "vtkCell3D.h"
34 #include "vtkCommonDataModelModule.h" // For export macro
35 
36 class vtkLine;
37 class vtkTriangle;
40 
41 class VTKCOMMONDATAMODEL_EXPORT vtkTetra : public vtkCell3D
42 {
43 public:
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  // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
54  VTK_DEPRECATED_IN_9_0_0("Replaced by vtkTetra::GetEdgePoints(vtkIdType, const vtkIdType*&)")
55  void GetEdgePoints(int edgeId, int*& pts) override;
56  vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
57  // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
58  VTK_DEPRECATED_IN_9_0_0("Replaced by vtkTetra::GetFacePoints(vtkIdType, const vtkIdType*&)")
59  void GetFacePoints(int faceId, int*& pts) override;
60  void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
61  vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
62  vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
63  vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
64  vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
65  bool GetCentroid(double centroid[3]) const override;
66  bool IsInsideOut() override;
68 
72  static constexpr vtkIdType NumberOfPoints = 4;
73 
77  static constexpr vtkIdType NumberOfEdges = 6;
78 
82  static constexpr vtkIdType NumberOfFaces = 4;
83 
88  static constexpr vtkIdType MaximumFaceSize = 3;
89 
95  static constexpr vtkIdType MaximumValence = 3;
96 
98 
101  int GetCellType() override { return VTK_TETRA; }
102  int GetNumberOfEdges() override { return 6; }
103  int GetNumberOfFaces() override { return 4; }
104  vtkCell* GetEdge(int edgeId) override;
105  vtkCell* GetFace(int faceId) override;
106  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
107  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
108  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
109  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
110  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
111  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
112  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
113  double& dist2, double weights[]) override;
114  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
115  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
116  double pcoords[3], int& subId) override;
117  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
119  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
120  double* GetParametricCoords() override;
122 
130  static int* GetTriangleCases(int caseId);
131 
137  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
138 
142  int GetParametricCenter(double pcoords[3]) override;
143 
148  double GetParametricDistance(const double pcoords[3]) override;
149 
153  static void TetraCenter(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
154 
160  static double Circumsphere(
161  double x1[3], double x2[3], double x3[3], double x4[3], double center[3]);
162 
168  static double Insphere(double p1[3], double p2[3], double p3[3], double p4[3], double center[3]);
169 
182  static int BarycentricCoords(
183  double x[3], double x1[3], double x2[3], double x3[3], double x4[3], double bcoords[4]);
184 
189  static double ComputeVolume(double p1[3], double p2[3], double p3[3], double p4[3]);
190 
196  int JacobianInverse(double** inverse, double derivs[12]);
197 
198  static void InterpolationFunctions(const double pcoords[3], double weights[4]);
199  static void InterpolationDerivs(const double pcoords[3], double derivs[12]);
201 
205  void InterpolateFunctions(const double pcoords[3], double weights[4]) override
206  {
207  vtkTetra::InterpolationFunctions(pcoords, weights);
208  }
209  void InterpolateDerivs(const double pcoords[3], double derivs[12]) override
210  {
211  vtkTetra::InterpolationDerivs(pcoords, derivs);
212  }
214 
216 
224  static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
225  static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(3);
227 
232 
237 
242 
247 
252 
256  static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
257 
258 protected:
260  ~vtkTetra() override;
261 
264 
265 private:
266  vtkTetra(const vtkTetra&) = delete;
267  void operator=(const vtkTetra&) = delete;
268 };
269 
270 inline int vtkTetra::GetParametricCenter(double pcoords[3])
271 {
272  pcoords[0] = pcoords[1] = pcoords[2] = 0.25;
273  return 0;
274 }
275 
276 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:40
object to represent cell connectivity
Definition: vtkCellArray.h:181
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:58
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
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
Definition: vtkPointData.h:33
represent and manipulate 3D points
Definition: vtkPoints.h:34
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:42
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 * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
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.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
int GetParametricCenter(double pcoords[3]) override
Return the center of the tetrahedron in parametric coordinates.
Definition: vtkTetra.h:270
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
vtkTriangle * Triangle
Definition: vtkTetra.h:263
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
static int * GetTriangleCases(int caseId)
Return the case table for table-based isocontouring (aka marching cubes style implementations).
~vtkTetra() override
void InterpolateFunctions(const double pcoords[3], double weights[4]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkTetra.h:205
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
vtkCell * GetFace(int faceId) override
See the vtkCell API for descriptions of these methods.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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 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.
static vtkTetra * New()
static const vtkIdType * GetPointToIncidentEdgesArray(vtkIdType pointId)
Static method version of GetPointToIncidentEdgesArray.
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkTetra.h:102
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:103
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...
double * GetParametricCoords() override
See the vtkCell API for descriptions of these methods.
vtkLine * Line
Definition: vtkTetra.h:262
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.
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:209
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
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:56
#define VTK_DEPRECATED_IN_9_0_0(reason)
int vtkIdType
Definition: vtkType.h:332
#define VTK_SIZEHINT(...)