VTK  9.2.6
vtkQuadraticPolygon.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticPolygon.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 =========================================================================*/
36 #ifndef vtkQuadraticPolygon_h
37 #define vtkQuadraticPolygon_h
38 
39 #include "vtkCommonDataModelModule.h" // For export macro
40 #include "vtkNonLinearCell.h"
41 
42 class vtkQuadraticEdge;
43 class vtkPolygon;
44 class vtkIdTypeArray;
45 
46 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPolygon : public vtkNonLinearCell
47 {
48 public:
49  static vtkQuadraticPolygon* New();
51  void PrintSelf(ostream& os, vtkIndent indent) override;
52 
57  int GetCellType() override { return VTK_QUADRATIC_POLYGON; }
58  int GetCellDimension() override { return 2; }
59  int GetNumberOfEdges() override { return this->GetNumberOfPoints() / 2; }
60  int GetNumberOfFaces() override { return 0; }
61  vtkCell* GetEdge(int) override;
62  vtkCell* GetFace(int) override { return nullptr; }
63  int IsPrimaryCell() override { return 0; }
64 
66 
72  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
73  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
74  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
75  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
76  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
77  vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
78  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
79  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
80  double& dist2, double weights[]) override;
81  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
82  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
83  double pcoords[3], int& subId) override;
84  void InterpolateFunctions(const double x[3], double* weights) override;
85  static void ComputeCentroid(vtkIdTypeArray* ids, vtkPoints* pts, double centroid[3]);
86  int ParameterizePolygon(
87  double p0[3], double p10[3], double& l10, double p20[3], double& l20, double n[3]);
88  static int PointInPolygon(double x[3], int numPts, double* pts, double bounds[6], double n[3]);
89  int Triangulate(vtkIdList* outTris);
90  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
91  int NonDegenerateTriangulate(vtkIdList* outTris);
92  static double DistanceToPolygon(
93  double x[3], int numPts, double* pts, double bounds[6], double closest[3]);
94  static int IntersectPolygonWithPolygon(int npts, double* pts, double bounds[6], int npts2,
95  double* pts2, double bounds2[6], double tol, double x[3]);
96  static int IntersectConvex2DCells(
97  vtkCell* cell1, vtkCell* cell2, double tol, double p0[3], double p1[3]);
99 
100  // Not implemented
101  void Derivatives(
102  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
103 
105 
111  vtkGetMacro(UseMVCInterpolation, bool);
112  vtkSetMacro(UseMVCInterpolation, bool);
114 
115 protected:
117  ~vtkQuadraticPolygon() override;
118 
119  // variables used by instances of this class
122 
123  // Parameter indicating whether to use Mean Value Coordinate algorithm
124  // for interpolation. The parameter is true by default.
126 
128 
132  static void GetPermutationFromPolygon(vtkIdType nb, vtkIdList* permutation);
133  static void PermuteToPolygon(vtkIdType nbPoints, double* inPoints, double* outPoints);
134  static void PermuteToPolygon(vtkCell* inCell, vtkCell* outCell);
135  static void PermuteToPolygon(vtkPoints* inPoints, vtkPoints* outPoints);
136  static void PermuteToPolygon(vtkIdTypeArray* inIds, vtkIdTypeArray* outIds);
137  static void PermuteToPolygon(vtkDataArray* inDataArray, vtkDataArray* outDataArray);
138  void InitializePolygon();
140 
142 
146  static void GetPermutationToPolygon(vtkIdType nb, vtkIdList* permutation);
147  static void PermuteFromPolygon(vtkIdType nb, double* values);
148  static void ConvertFromPolygon(vtkIdList* ids);
150 
151 private:
152  vtkQuadraticPolygon(const vtkQuadraticPolygon&) = delete;
153  void operator=(const vtkQuadraticPolygon&) = delete;
154 };
155 
156 #endif
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
vtkQuadraticEdge * Edge
represent and manipulate point attribute data
Definition: vtkPointData.h:41
vtkIdType GetNumberOfPoints() const
Return the number of points in the cell.
Definition: vtkCell.h:141
represent and manipulate cell attribute data
Definition: vtkCellData.h:41
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:391
abstract superclass for non-linear cells
dynamic, self-adjusting array of vtkIdType
int vtkIdType
Definition: vtkType.h:332
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
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...
abstract class to specify cell behavior
Definition: vtkCell.h:60
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
list of point or cell ids
Definition: vtkIdList.h:33
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:42
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.
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:186
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
a cell that represents a parabolic n-sided polygon
int GetCellType() override
Implement the vtkCell API.
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.
int GetNumberOfFaces() override
Return the number of faces in the cell.
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
int IsPrimaryCell() override
Return whether this cell type has a fixed topology or whether the topology varies depending on the da...
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 GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
represent and manipulate 3D points
Definition: vtkPoints.h:39
int GetNumberOfEdges() override
Return the number of edges in the cell.