VTK
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 =========================================================================*/
33 #ifndef vtkQuadraticPolygon_h
34 #define vtkQuadraticPolygon_h
35 
36 #include "vtkCommonDataModelModule.h" // For export macro
37 #include "vtkNonLinearCell.h"
38 
39 class vtkQuadraticEdge;
40 class vtkPolygon;
41 class vtkIdTypeArray;
42 
43 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPolygon : public vtkNonLinearCell
44 {
45 public:
46  static vtkQuadraticPolygon *New();
48  void PrintSelf(ostream& os, vtkIndent indent) override;
49 
54  int GetCellType() override { return VTK_QUADRATIC_POLYGON; }
55  int GetCellDimension() override { return 2;}
56  int GetNumberOfEdges() override { return this->GetNumberOfPoints() / 2; }
57  int GetNumberOfFaces() override { return 0; }
58  vtkCell *GetEdge(int) override;
59  vtkCell *GetFace(int) override { return nullptr; }
60  int IsPrimaryCell() override { return 0; }
61 
63 
69  int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
70  void Contour(double value, vtkDataArray *cellScalars,
72  vtkCellArray *lines, vtkCellArray *polys,
73  vtkPointData *inPd, vtkPointData *outPd,
74  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override;
75  void Clip(double value, vtkDataArray *cellScalars,
77  vtkPointData *inPd, vtkPointData *outPd,
78  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
79  int insideOut) override;
80  int EvaluatePosition(const double x[3], double closestPoint[3],
81  int& subId, double pcoords[3],
82  double& dist2, double weights[]) override;
83  void EvaluateLocation(int& subId, const double pcoords[3], double x[3],
84  double *weights) override;
85  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
86  double x[3], double pcoords[3], int& subId) override;
87  void InterpolateFunctions(const double x[3], double *weights) override;
88  static void ComputeCentroid(vtkIdTypeArray *ids, vtkPoints *pts,
89  double centroid[3]);
90  int ParameterizePolygon(double p0[3], double p10[3], double &l10,
91  double p20[3], double &l20, double n[3]);
92  static int PointInPolygon(double x[3], int numPts, double *pts,
93  double bounds[6], double n[3]);
94  int Triangulate(vtkIdList *outTris);
95  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
96  int NonDegenerateTriangulate(vtkIdList *outTris);
97  static double DistanceToPolygon(double x[3], int numPts, double *pts,
98  double bounds[6], double closest[3]);
99  static int IntersectPolygonWithPolygon(int npts, double *pts, double bounds[6],
100  int npts2, double *pts2,
101  double bounds2[3], double tol,
102  double x[3]);
103  static int IntersectConvex2DCells(vtkCell *cell1, vtkCell *cell2,
104  double tol, double p0[3], double p1[3]);
106 
107  // Not implemented
108  void Derivatives(int subId, const double pcoords[3], const double *values,
109  int dim, double *derivs) override;
110 
112 
118  vtkGetMacro(UseMVCInterpolation, bool);
119  vtkSetMacro(UseMVCInterpolation, bool);
121 
122 protected:
124  ~vtkQuadraticPolygon() override;
125 
126  // variables used by instances of this class
129 
130  // Parameter indicating whether to use Mean Value Coordinate algorithm
131  // for interpolation. The parameter is true by default.
133 
135 
139  static void GetPermutationFromPolygon(vtkIdType nb, vtkIdList *permutation);
140  static void PermuteToPolygon(vtkIdType nbPoints, double *inPoints, double *outPoints);
141  static void PermuteToPolygon(vtkCell *inCell, vtkCell *outCell);
142  static void PermuteToPolygon(vtkPoints *inPoints, vtkPoints *outPoints);
143  static void PermuteToPolygon(vtkIdTypeArray *inIds, vtkIdTypeArray *outIds);
144  static void PermuteToPolygon(vtkDataArray *inDataArray, vtkDataArray *outDataArray);
145  void InitializePolygon();
147 
149 
153  static void GetPermutationToPolygon(vtkIdType nb, vtkIdList *permutation);
154  static void PermuteFromPolygon(vtkIdType nb, double *values);
155  static void ConvertFromPolygon(vtkIdList *ids);
157 
158 private:
159  vtkQuadraticPolygon(const vtkQuadraticPolygon&) = delete;
160  void operator=(const vtkQuadraticPolygon&) = delete;
161 };
162 
163 #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.
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.
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkX3D::value
Definition: vtkX3D.h:220
vtkIdType
int vtkIdType
Definition: vtkType.h:347
vtkPolygon
a cell that represents an n-sided polygon
Definition: vtkPolygon.h:39
vtkObject::New
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
vtkQuadraticPolygon::IsPrimaryCell
int IsPrimaryCell() override
Return whether this cell type has a fixed topology or whether the topology varies depending on the da...
Definition: vtkQuadraticPolygon.h:60
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
vtkQuadraticPolygon
a cell that represents a parabolic n-sided polygon
Definition: vtkQuadraticPolygon.h:43
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
VTK_QUADRATIC_POLYGON
Definition: vtkCellType.h:67
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:44
vtkQuadraticPolygon::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkQuadraticPolygon.h:56
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
vtkCell::GetNumberOfPoints
vtkIdType GetNumberOfPoints()
Return the number of points in the cell.
Definition: vtkCell.h:137
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
vtkQuadraticPolygon::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkQuadraticPolygon.h:54
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...
vtkIdTypeArray
dynamic, self-adjusting array of vtkIdType
Definition: vtkIdTypeArray.h:35
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.
vtkQuadraticPolygon::Polygon
vtkPolygon * Polygon
Definition: vtkQuadraticPolygon.h:127
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:354
vtkQuadraticPolygon::GetFace
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
Definition: vtkQuadraticPolygon.h:59
vtkQuadraticPolygon::Edge
vtkQuadraticEdge * Edge
Definition: vtkQuadraticPolygon.h:128
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.
vtkQuadraticPolygon::UseMVCInterpolation
bool UseMVCInterpolation
Definition: vtkQuadraticPolygon.h:132
vtkQuadraticPolygon::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkQuadraticPolygon.h:57
vtkQuadraticPolygon::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkQuadraticPolygon.h:55
vtkX3D::index
Definition: vtkX3D.h:246
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:40