VTK  9.0.1
vtkHigherOrderTriangle.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkHigherOrderTriangle.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 vtkHigherOrderTriangle_h
32 #define vtkHigherOrderTriangle_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 vtkDoubleArray;
45 class vtkTriangle;
46 
47 class VTKCOMMONDATAMODEL_EXPORT vtkHigherOrderTriangle : public vtkNonLinearCell
48 {
49 public:
51  void PrintSelf(ostream& os, vtkIndent indent) override;
52 
53  int GetCellType() override = 0;
54  int GetCellDimension() override { return 2; }
55  int RequiresInitialization() override { return 1; }
56  int GetNumberOfEdges() override { return 3; }
57  int GetNumberOfFaces() override { return 0; }
58  vtkCell* GetEdge(int edgeId) override = 0;
59  void SetEdgeIdsAndPoints(int edgeId,
60  const std::function<void(const vtkIdType&)>& set_number_of_ids_and_points,
61  const std::function<void(const vtkIdType&, const vtkIdType&)>& set_ids_and_points);
62  vtkCell* GetFace(int) override { return nullptr; }
63 
64  void Initialize() override;
65 
66  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
67  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
68  double& dist2, double weights[]) override;
69  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
70  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
71  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
72  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
73  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
74  vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
75  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
76  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
77  double pcoords[3], int& subId) override;
78  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
79  void JacobianInverse(const double pcoords[3], double** inverse, double* derivs);
80  void Derivatives(
81  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
82  void SetParametricCoords();
83  double* GetParametricCoords() override;
84 
85  int GetParametricCenter(double pcoords[3]) override;
86  double GetParametricDistance(const double pcoords[3]) override;
87 
88  void InterpolateFunctions(const double pcoords[3], double* weights) override = 0;
89  void InterpolateDerivs(const double pcoords[3], double* derivs) override = 0;
90 
91  vtkIdType GetOrder() const { return this->Order; }
92  vtkIdType ComputeOrder();
93 
94  void ToBarycentricIndex(vtkIdType index, vtkIdType* bindex);
95  vtkIdType ToIndex(const vtkIdType* bindex);
96 
97  static void BarycentricIndex(vtkIdType index, vtkIdType* bindex, vtkIdType order);
98  static vtkIdType Index(const vtkIdType* bindex, vtkIdType order);
99 
100  static double eta(vtkIdType n, vtkIdType chi, double sigma);
101  static double d_eta(vtkIdType n, vtkIdType chi, double sigma);
102  virtual vtkHigherOrderCurve* getEdgeCell() = 0;
103 
104 protected:
106  ~vtkHigherOrderTriangle() override;
107 
108  vtkIdType GetNumberOfSubtriangles() const { return this->NumberOfSubtriangles; }
109  vtkIdType ComputeNumberOfSubtriangles();
110 
111  // Description:
112  // Given the index of the subtriangle, compute the barycentric indices of
113  // the subtriangle's vertices.
114  void SubtriangleBarycentricPointIndices(vtkIdType cellIndex, vtkIdType (&pointBIndices)[3][3]);
115 
117  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
121 
122  std::vector<vtkIdType> BarycentricIndexMap;
123  std::vector<vtkIdType> IndexMap;
124  std::vector<vtkIdType> SubtriangleIndexMap;
125 
126 private:
128  void operator=(const vtkHigherOrderTriangle&) = delete;
129 };
130 
131 #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.
vtkHigherOrderTriangle::SubtriangleIndexMap
std::vector< vtkIdType > SubtriangleIndexMap
Definition: vtkHigherOrderTriangle.h:124
vtkX3D::function
@ function
Definition: vtkX3D.h:255
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.
vtkHigherOrderTriangle::IndexMap
std::vector< vtkIdType > IndexMap
Definition: vtkHigherOrderTriangle.h:123
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkHigherOrderTriangle::NumberOfSubtriangles
vtkIdType NumberOfSubtriangles
Definition: vtkHigherOrderTriangle.h:119
vtkHigherOrderTriangle::GetNumberOfSubtriangles
vtkIdType GetNumberOfSubtriangles() const
Definition: vtkHigherOrderTriangle.h:108
vtkX3D::value
@ value
Definition: vtkX3D.h:226
vtkHigherOrderTriangle::Scalars
vtkDoubleArray * Scalars
Definition: vtkHigherOrderTriangle.h:117
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkHigherOrderTriangle::PointParametricCoordinates
vtkSmartPointer< vtkPoints > PointParametricCoordinates
Definition: vtkHigherOrderTriangle.h:120
vtkCell::Initialize
virtual void Initialize()
Definition: vtkCell.h:111
vtkSmartPointer< vtkPoints >
vtkHigherOrderTriangle::GetFace
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
Definition: vtkHigherOrderTriangle.h:62
vtkHigherOrderTriangle::Order
vtkIdType Order
Definition: vtkHigherOrderTriangle.h:118
vtkHigherOrderTriangle::GetOrder
vtkIdType GetOrder() const
Definition: vtkHigherOrderTriangle.h:91
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
vtkCell::GetCellType
virtual int GetCellType()=0
Return the type of cell.
vtkHigherOrderCurve
Definition: vtkHigherOrderCurve.h:37
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.
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:56
vtkHigherOrderTriangle::Face
vtkTriangle * Face
Definition: vtkHigherOrderTriangle.h:116
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:179
vtkSmartPointer.h
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkHigherOrderTriangle::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkHigherOrderTriangle.h:54
vtkTriangle
a cell that represents a triangle
Definition: vtkTriangle.h:35
vtkX3D::order
@ order
Definition: vtkX3D.h:446
vtkHigherOrderTriangle::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkHigherOrderTriangle.h:57
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...
vtkNonLinearCell.h
vtkCell::InterpolateDerivs
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:344
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...
vtkHigherOrderTriangle::RequiresInitialization
int RequiresInitialization() override
Some cells require initialization prior to access.
Definition: vtkHigherOrderTriangle.h:55
vtkCell::GetParametricCoords
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkNew.h
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::InterpolateFunctions
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:341
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
vtkHigherOrderTriangle::BarycentricIndexMap
std::vector< vtkIdType > BarycentricIndexMap
Definition: vtkHigherOrderTriangle.h:122
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkHigherOrderTriangle
A 2D cell that represents an arbitrary order HigherOrder triangle.
Definition: vtkHigherOrderTriangle.h:47
vtkHigherOrderTriangle::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkHigherOrderTriangle.h:56