VTK  9.0.20200803
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:
51  static vtkQuadraticTetra* New();
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;
76  void Derivatives(
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
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:33
vtkCell::IntersectWithLine
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
vtkQuadraticTetra::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkQuadraticTetra.h:62
vtkQuadraticTetra::Edge
vtkQuadraticEdge * Edge
Definition: vtkQuadraticTetra.h:145
vtkCell::Contour
virtual void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd)=0
Generate contouring primitives.
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:32
vtkQuadraticTetra::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkQuadraticTetra.h:61
vtkX3D::value
Definition: vtkX3D.h:226
vtkIdType
int vtkIdType
Definition: vtkType.h:330
vtkObject::New
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
vtkQuadraticTetra::Face
vtkQuadraticTriangle * Face
Definition: vtkQuadraticTetra.h:146
vtkQuadraticTetra::Scalars
vtkDoubleArray * Scalars
Definition: vtkQuadraticTetra.h:148
vtkQuadraticTriangle
cell represents a parabolic, isoparametric triangle
Definition: vtkQuadraticTriangle.h:43
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
VTK_QUADRATIC_TETRA
Definition: vtkCellType.h:69
vtkCell::EvaluateLocation
virtual void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkCell::Triangulate
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
vtkQuadraticTetra::InterpolationDerivs
static void InterpolationDerivs(const double pcoords[3], double derivs[30])
vtkQuadraticTetra::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkQuadraticTetra.h:63
vtkQuadraticTetra
cell represents a parabolic, 10-node isoparametric tetrahedron
Definition: vtkQuadraticTetra.h:48
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:57
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
vtkCell::GetFace
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:180
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkQuadraticTetra::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkQuadraticTetra.h:60
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkCell::CellBoundary
virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
vtkQuadraticTetra::InterpolateDerivs
void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
Definition: vtkQuadraticTetra.h:117
vtkNonLinearCell.h
vtkCell::GetParametricDistance
virtual double GetParametricDistance(const double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
vtkCell::EvaluatePosition
virtual int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[])=0
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
vtkCell::GetParametricCoords
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkNonLinearCell
abstract superclass for non-linear cells
Definition: vtkNonLinearCell.h:36
vtkCell::Clip
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
Cut (or clip) the cell based on the input cellScalars and the specified value.
vtkCell::GetEdge
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkNonLinearCell::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkCell::Derivatives
virtual void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
vtkCell::GetParametricCenter
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:35
vtkQuadraticTetra::Tetra
vtkTetra * Tetra
Definition: vtkQuadraticTetra.h:147
vtkTetra
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:41
vtkQuadraticTetra::InterpolationFunctions
static void InterpolationFunctions(const double pcoords[3], double weights[10])
vtkX3D::index
Definition: vtkX3D.h:252
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:40
vtkQuadraticTetra::InterpolateFunctions
void InterpolateFunctions(const double pcoords[3], double weights[10]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkQuadraticTetra.h:113