VTK  9.0.2
vtkQuadraticTetra.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticTetra.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 =========================================================================*/
37 #ifndef vtkQuadraticTetra_h
38 #define vtkQuadraticTetra_h
39 
40 #include "vtkCommonDataModelModule.h" // For export macro
41 #include "vtkNonLinearCell.h"
42 
43 class vtkQuadraticEdge;
45 class vtkTetra;
46 class vtkDoubleArray;
47 
48 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticTetra : public vtkNonLinearCell
49 {
50 public:
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
56 
60  int GetCellType() override { return VTK_QUADRATIC_TETRA; }
61  int GetCellDimension() override { return 3; }
62  int GetNumberOfEdges() override { return 6; }
63  int GetNumberOfFaces() override { return 4; }
64  vtkCell* GetEdge(int) override;
65  vtkCell* GetFace(int) override;
67 
68  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
69  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
70  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
71  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
72  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
73  double& dist2, double weights[]) override;
74  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
75  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
77  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
78  double* GetParametricCoords() override;
79 
84  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
85  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
86  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
87 
92  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
93  double pcoords[3], int& subId) override;
94 
98  int GetParametricCenter(double pcoords[3]) override;
99 
104  double GetParametricDistance(const double pcoords[3]) override;
105 
106  static void InterpolationFunctions(const double pcoords[3], double weights[10]);
107  static void InterpolationDerivs(const double pcoords[3], double derivs[30]);
109 
113  void InterpolateFunctions(const double pcoords[3], double weights[10]) override
114  {
116  }
117  void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
118  {
120  }
122 
123 
130  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
131  static const vtkIdType* GetFaceArray(vtkIdType faceId);
133 
139  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[30]);
140 
141 protected:
143  ~vtkQuadraticTetra() override;
144 
148  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
149 
150 private:
151  vtkQuadraticTetra(const vtkQuadraticTetra&) = delete;
152  void operator=(const vtkQuadraticTetra&) = delete;
153 };
154 
155 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:180
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:57
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
dynamic, self-adjusting array of double
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
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:32
represent and manipulate 3D points
Definition: vtkPoints.h:34
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 10-node isoparametric tetrahedron
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic tetra in parametric coordinates.
int GetNumberOfEdges() override
Return the number of edges in the cell.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
int GetNumberOfFaces() override
Return the number of faces in the cell.
vtkQuadraticEdge * Edge
static void InterpolationDerivs(const double pcoords[3], double derivs[30])
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Line-edge intersection.
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[30])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
void InterpolateFunctions(const double pcoords[3], double weights[10]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
~vtkQuadraticTetra() override
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
vtkQuadraticTriangle * Face
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
static vtkQuadraticTetra * New()
static void InterpolationFunctions(const double pcoords[3], double weights[10])
static const vtkIdType * GetFaceArray(vtkIdType faceId)
int GetCellType() override
Implement the vtkCell API.
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tetras, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this edge using scalar value provided.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkDoubleArray * Scalars
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
Generate contouring primitives.
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:42
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_TETRA
Definition: vtkCellType.h:69
int vtkIdType
Definition: vtkType.h:338