VTK
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 =========================================================================*/
40 #ifndef vtkQuadraticTetra_h
41 #define vtkQuadraticTetra_h
42 
43 #include "vtkCommonDataModelModule.h" // For export macro
44 #include "vtkNonLinearCell.h"
45 
46 class vtkQuadraticEdge;
48 class vtkTetra;
49 class vtkDoubleArray;
50 
51 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticTetra : public vtkNonLinearCell
52 {
53 public:
54  static vtkQuadraticTetra *New();
56  void PrintSelf(ostream& os, vtkIndent indent) override;
57 
59 
63  int GetCellType() override {return VTK_QUADRATIC_TETRA;}
64  int GetCellDimension() override {return 3;}
65  int GetNumberOfEdges() override {return 6;}
66  int GetNumberOfFaces() override {return 4;}
67  vtkCell *GetEdge(int) override;
68  vtkCell *GetFace(int) override;
70 
71  int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
72  void Contour(double value, vtkDataArray *cellScalars,
74  vtkCellArray *lines, vtkCellArray *polys,
75  vtkPointData *inPd, vtkPointData *outPd,
76  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override;
77  int EvaluatePosition(const double x[3], double closestPoint[3],
78  int& subId, double pcoords[3],
79  double& dist2, double weights[]) override;
80  void EvaluateLocation(int& subId, const double pcoords[3], double x[3],
81  double *weights) override;
82  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
83  void Derivatives(int subId, const double pcoords[3], const double *values,
84  int dim, double *derivs) override;
85  double *GetParametricCoords() override;
86 
91  void Clip(double value, vtkDataArray *cellScalars,
92  vtkIncrementalPointLocator *locator, vtkCellArray *tetras,
93  vtkPointData *inPd, vtkPointData *outPd,
94  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
95  int insideOut) override;
96 
101  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
102  double x[3], double pcoords[3], int& subId) override;
103 
104 
108  int GetParametricCenter(double pcoords[3]) override;
109 
114  double GetParametricDistance(const double pcoords[3]) override;
115 
119  static void InterpolationFunctions(const double pcoords[3], double weights[10]);
123  static void InterpolationDerivs(const double pcoords[3], double derivs[30]);
125 
129  void InterpolateFunctions(const double pcoords[3], double weights[10]) override
130  {
132  }
133  void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
134  {
136  }
138 
139 
143  static int *GetEdgeArray(int edgeId);
144  static int *GetFaceArray(int faceId);
146 
152  void JacobianInverse(const double pcoords[3], double **inverse, double derivs[30]);
153 
154 protected:
156  ~vtkQuadraticTetra() override;
157 
161  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
162 
163 private:
164  vtkQuadraticTetra(const vtkQuadraticTetra&) = delete;
165  void operator=(const vtkQuadraticTetra&) = delete;
166 };
167 
168 #endif
169 
170 
int GetNumberOfEdges() override
Implement the vtkCell API.
static void InterpolationFunctions(const double pcoords[3], double weights[10])
represent and manipulate point attribute data
Definition: vtkPointData.h:37
vtkDoubleArray * Scalars
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
Abstract class in support of both point location and point insertion.
vtkQuadraticEdge * Edge
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:347
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkQuadraticTriangle * Face
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...
dynamic, self-adjusting array of double
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:47
virtual double GetParametricDistance(const double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
int GetCellType() override
Implement the vtkCell API.
abstract class to specify cell behavior
Definition: vtkCell.h:59
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.
a simple class to control print indentation
Definition: vtkIndent.h:39
static void InterpolationDerivs(const double pcoords[3], double derivs[30])
list of point or cell ids
Definition: vtkIdList.h:36
int GetNumberOfFaces() override
Implement the vtkCell API.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
cell represents a parabolic, 10-node isoparametric tetrahedron
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.
void InterpolateDerivs(const double pcoords[3], double derivs[30]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
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...
object to represent cell connectivity
Definition: vtkCellArray.h:50
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
cell represents a parabolic, isoparametric edge
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.
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.
cell represents a parabolic, isoparametric triangle
void InterpolateFunctions(const double pcoords[3], double weights[10]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
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.
represent and manipulate 3D points
Definition: vtkPoints.h:39
int GetCellDimension() override
Implement the vtkCell API.