VTK  9.0.3
vtkHigherOrderTetra.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkHigherOrderTetra.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 vtkHigherOrderTetra_h
32 #define vtkHigherOrderTetra_h
33 
34 #include <functional> //For std::function
35 
36 #include "vtkCommonDataModelModule.h" // For export macro
37 #include "vtkNew.h" // For member variable.
38 #include "vtkNonLinearCell.h"
39 #include "vtkSmartPointer.h" // For member variable.
40 
41 #include <vector> //For caching
42 
43 class vtkTetra;
46 class vtkDoubleArray;
47 
48 class VTKCOMMONDATAMODEL_EXPORT vtkHigherOrderTetra : public vtkNonLinearCell
49 {
50 public:
52  void PrintSelf(ostream& os, vtkIndent indent) override;
53 
54  int GetCellType() override = 0;
55  int GetCellDimension() override { return 3; }
56  int RequiresInitialization() override { return 1; }
57  int GetNumberOfEdges() override { return 6; }
58  int GetNumberOfFaces() override { return 4; }
59  vtkCell* GetEdge(int edgeId) override = 0;
60  vtkCell* GetFace(int faceId) override = 0;
61  void SetEdgeIdsAndPoints(int edgeId,
62  const std::function<void(const vtkIdType&)>& set_number_of_ids_and_points,
63  const std::function<void(const vtkIdType&, const vtkIdType&)>& set_ids_and_points);
64  void SetFaceIdsAndPoints(vtkHigherOrderTriangle* result, int edgeId,
65  const std::function<void(const vtkIdType&)>& set_number_of_ids_and_points,
66  const std::function<void(const vtkIdType&, const vtkIdType&)>& set_ids_and_points);
67 
68  void Initialize() override;
69 
70  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
71  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
72  double& dist2, double weights[]) override;
73  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
74  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
75  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
76  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
77  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
78  vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
79  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
80  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
81  double pcoords[3], int& subId) override;
82  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
83  void JacobianInverse(const double pcoords[3], double** inverse, double* derivs);
85  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
87  double* GetParametricCoords() override;
88 
89  int GetParametricCenter(double pcoords[3]) override;
90  double GetParametricDistance(const double pcoords[3]) override;
91 
92  void InterpolateFunctions(const double pcoords[3], double* weights) override = 0;
93  void InterpolateDerivs(const double pcoords[3], double* derivs) override = 0;
94 
95  vtkIdType GetOrder() const { return this->Order; }
97  static vtkIdType ComputeOrder(const vtkIdType nPoints);
98 
100  vtkIdType ToIndex(const vtkIdType* bindex);
101 
103  static vtkIdType Index(const vtkIdType* bindex, vtkIdType order);
106 
107 protected:
110 
111  vtkIdType GetNumberOfSubtetras() const { return this->NumberOfSubtetras; }
113 
114  // Description:
115  // Given the index of the subtriangle, compute the barycentric indices of
116  // the subtriangle's vertices.
117  void SubtetraBarycentricPointIndices(vtkIdType cellIndex, vtkIdType (&pointBIndices)[4][4]);
119  vtkIdType cellIndex, const vtkIdType (&octBIndices)[6][4], vtkIdType (&tetraBIndices)[4][4]);
120 
122  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
126 
127  std::vector<vtkIdType> EdgeIds;
128  std::vector<vtkIdType> BarycentricIndexMap;
129  std::vector<vtkIdType> IndexMap;
130  std::vector<vtkIdType> SubtetraIndexMap;
131 
132 private:
133  vtkHigherOrderTetra(const vtkHigherOrderTetra&) = delete;
134  void operator=(const vtkHigherOrderTetra&) = delete;
135 };
136 
137 #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
A 3D cell that represents an arbitrary order HigherOrder tetrahedron.
vtkIdType GetOrder() const
void InterpolateFunctions(const double pcoords[3], double *weights) override=0
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
~vtkHigherOrderTetra() override
void ToBarycentricIndex(vtkIdType index, vtkIdType *bindex)
int GetNumberOfFaces() override
Return the number of faces in the cell.
void TetraFromOctahedron(vtkIdType cellIndex, const vtkIdType(&octBIndices)[6][4], vtkIdType(&tetraBIndices)[4][4])
void SetFaceIdsAndPoints(vtkHigherOrderTriangle *result, int edgeId, const std::function< void(const vtkIdType &)> &set_number_of_ids_and_points, const std::function< void(const vtkIdType &, const vtkIdType &)> &set_ids_and_points)
vtkIdType ComputeOrder()
void InterpolateDerivs(const double pcoords[3], double *derivs) override=0
static vtkIdType Index(const vtkIdType *bindex, vtkIdType order)
virtual vtkHigherOrderCurve * getEdgeCell()=0
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
static vtkIdType ComputeOrder(const vtkIdType nPoints)
std::vector< vtkIdType > EdgeIds
void SetEdgeIdsAndPoints(int edgeId, const std::function< void(const vtkIdType &)> &set_number_of_ids_and_points, const std::function< void(const vtkIdType &, const vtkIdType &)> &set_ids_and_points)
vtkSmartPointer< vtkPoints > PointParametricCoordinates
vtkDoubleArray * Scalars
void SubtetraBarycentricPointIndices(vtkIdType cellIndex, vtkIdType(&pointBIndices)[4][4])
int GetNumberOfEdges() override
Return the number of edges in 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
Intersect with a ray.
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.
std::vector< vtkIdType > SubtetraIndexMap
vtkCell * GetEdge(int edgeId) override=0
Return the edge cell from the edgeId of the cell.
int RequiresInitialization() override
Some cells require initialization prior to access.
std::vector< vtkIdType > BarycentricIndexMap
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...
vtkIdType GetNumberOfSubtetras() const
vtkIdType ToIndex(const vtkIdType *bindex)
static void BarycentricIndex(vtkIdType index, vtkIdType *bindex, vtkIdType order)
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Cut (or clip) the cell based on the input cellScalars and the specified value.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
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...
virtual vtkHigherOrderTriangle * getFaceCell()=0
int GetCellType() override=0
Return the type of cell.
vtkIdType ComputeNumberOfSubtetras()
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
void Initialize() override
vtkCell * GetFace(int faceId) override=0
Return the face cell from the faceId of the cell.
std::vector< vtkIdType > IndexMap
void JacobianInverse(const double pcoords[3], double **inverse, double *derivs)
int GetParametricCenter(double pcoords[3]) override
Return center of the cell in 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 HigherOrder triangle.
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
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:42
@ order
Definition: vtkX3D.h:446
@ function
Definition: vtkX3D.h:255
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
int vtkIdType
Definition: vtkType.h:338