VTK
vtkLagrangeTetra.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkLagrangeTetra.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 =========================================================================*/
31 #ifndef vtkLagrangeTetra_h
32 #define vtkLagrangeTetra_h
33 
34 #include "vtkCommonDataModelModule.h" // For export macro
35 #include "vtkNonLinearCell.h"
36 
37 #define VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER 6
38 
39 #define MAX_POINTS ((VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER + 1) * \
40  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER + 2) * \
41  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER + 3)/6)
42 
43 #define MAX_SUBTETRAHEDRA (((VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER - 2) * \
44  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER - 1) * \
45  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER)) + \
46  4*((VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER - 1) * \
47  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER) * \
48  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER + 1)) + \
49  ((VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER) * \
50  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER + 1) * \
51  (VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER + 2))/6)
52 
53 class vtkTetra;
54 class vtkLagrangeCurve;
56 class vtkDoubleArray;
57 
58 class VTKCOMMONDATAMODEL_EXPORT vtkLagrangeTetra : public vtkNonLinearCell
59 {
60 public:
61  static vtkLagrangeTetra *New();
63  void PrintSelf(ostream& os, vtkIndent indent) override;
64 
65  int GetCellType() override { return VTK_LAGRANGE_TETRAHEDRON; }
66  int GetCellDimension() override { return 3; }
67  int RequiresInitialization() override { return 1; }
68  int GetNumberOfEdges() override { return 6; }
69  int GetNumberOfFaces() override { return 4; }
70  vtkCell *GetEdge(int edgeId) override;
71  vtkCell *GetFace(int faceId) override;
72 
73  void Initialize() override;
74 
76  static int MaximumNumberOfPoints()
77  {
78  return ((vtkLagrangeTetra::MaximumOrder() + 1) *
80  }
81  int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
82  int EvaluatePosition(const double x[3], double closestPoint[3],
83  int& subId, double pcoords[3],
84  double& dist2, double weights[]) override;
85  void EvaluateLocation(int& subId, const double pcoords[3], double x[3],
86  double *weights) override;
87  void Contour(double value, vtkDataArray *cellScalars,
89  vtkCellArray *lines, vtkCellArray *polys,
90  vtkPointData *inPd, vtkPointData *outPd,
91  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override;
92  void Clip(double value, vtkDataArray *cellScalars,
94  vtkPointData *inPd, vtkPointData *outPd,
95  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
96  int insideOut) override;
97  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
98  double x[3], double pcoords[3], int& subId) override;
99  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
100  void JacobianInverse(const double pcoords[3], double** inverse, double* derivs);
101  void Derivatives(int subId, const double pcoords[3], const double *values,
102  int dim, double *derivs) override;
103  double* GetParametricCoords() override;
104 
105  int GetParametricCenter(double pcoords[3]) override;
106  double GetParametricDistance(const double pcoords[3]) override;
107 
108  void InterpolateFunctions(const double pcoords[3], double* weights) override;
109  void InterpolateDerivs(const double pcoords[3], double* derivs) override;
110 
111  vtkIdType GetOrder() const { return this->Order; }
112  vtkIdType ComputeOrder();
113 
114  void ToBarycentricIndex(vtkIdType index, vtkIdType* bindex);
115  vtkIdType ToIndex(const vtkIdType* bindex);
116 
117  static void BarycentricIndex(vtkIdType index, vtkIdType* bindex,
118  vtkIdType order);
119  static vtkIdType Index(const vtkIdType* bindex, vtkIdType order);
120 
121 protected:
123  ~vtkLagrangeTetra() override;
124 
125  vtkIdType GetNumberOfSubtetras() const { return this->NumberOfSubtetras; }
126  vtkIdType ComputeNumberOfSubtetras();
127 
128  // Description:
129  // Given the index of the subtriangle, compute the barycentric indices of
130  // the subtriangle's vertices.
131  void SubtetraBarycentricPointIndices(vtkIdType cellIndex,
132  vtkIdType (&pointBIndices)[4][4]);
133  void TetraFromOctahedron(vtkIdType cellIndex,
134  const vtkIdType (&octBIndices)[6][4],
135  vtkIdType (&tetraBIndices)[4][4]);
136 
140  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
144 
146  vtkIdType BarycentricIndexMap[4*MAX_POINTS];
150  vtkIdType SubtetraIndexMap[16*MAX_SUBTETRAHEDRA];
151 
152 private:
153  vtkLagrangeTetra(const vtkLagrangeTetra&) = delete;
154  void operator=(const vtkLagrangeTetra&) = delete;
155 };
156 
157 #undef MAX_POINTS
158 #undef MAX_SUBTETRAHEDRA
159 
160 #endif
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkIdType GetNumberOfSubtetras() const
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this leve...
Definition: vtkCell.h:354
abstract superclass for non-linear cells
int GetCellType() override
Return the type of cell.
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.
A 3D cell that represents an arbitrary order Lagrange tetrahedron.
int RequiresInitialization() override
Some cells require initialization prior to access.
vtkLagrangeTriangle * 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...
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:357
dynamic, self-adjusting array of double
double * ParametricCoordinates
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:41
virtual double GetParametricDistance(const double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
abstract class to specify cell behavior
Definition: vtkCell.h:56
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.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
A 2D cell that represents an arbitrary order Lagrange triangle.
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkDoubleArray * Scalars
list of point or cell ids
Definition: vtkIdList.h:30
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
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.
int GetNumberOfFaces() override
Return the number of faces in the cell.
static int MaximumOrder()
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:44
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkIdType NumberOfSubtetras
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.
vtkLagrangeCurve * Edge
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.
#define VTK_LAGRANGE_TETRAHEDRON_MAX_ORDER
#define MAX_SUBTETRAHEDRA
static int MaximumNumberOfPoints()
#define MAX_POINTS
virtual void Initialize()
Definition: vtkCell.h:111
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.
int GetNumberOfEdges() override
Return the number of edges in the cell.
vtkIdType GetOrder() const
represent and manipulate 3D points
Definition: vtkPoints.h:33