VTK
vtkBiQuadraticTriangle.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkBiQuadraticTriangle.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 =========================================================================*/
41 #ifndef vtkBiQuadraticTriangle_h
42 #define vtkBiQuadraticTriangle_h
43 
44 #include "vtkCommonDataModelModule.h" // For export macro
45 #include "vtkNonLinearCell.h"
46 
47 class vtkQuadraticEdge;
48 class vtkTriangle;
49 class vtkDoubleArray;
50 
51 class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticTriangle : public vtkNonLinearCell
52 {
53 public:
54  static vtkBiQuadraticTriangle *New();
56  void PrintSelf(ostream& os, vtkIndent indent) override;
57 
59 
63  int GetCellType() override {return VTK_BIQUADRATIC_TRIANGLE;};
64  int GetCellDimension() override {return 2;}
65  int GetNumberOfEdges() override {return 3;}
66  int GetNumberOfFaces() override {return 0;}
67  vtkCell *GetEdge(int edgeId) override;
68  vtkCell *GetFace(int) override {return nullptr;}
70 
71  int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
72  void Contour(double value, vtkDataArray *cellScalars,
74  vtkCellArray *lines, vtkCellArray *polys,
75  vtkPointData *inPd, vtkPointData *outPd,
76  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override;
77  int EvaluatePosition(const double x[3], double closestPoint[3],
78  int& subId, double pcoords[3],
79  double& dist2, double weights[]) override;
80  void EvaluateLocation(int& subId, const double pcoords[3], double x[3],
81  double *weights) override;
82  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
83  void Derivatives(int subId, const double pcoords[3], const double *values,
84  int dim, double *derivs) override;
85  double *GetParametricCoords() override;
86 
92  void Clip(double value, vtkDataArray *cellScalars,
94  vtkPointData *inPd, vtkPointData *outPd,
95  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
96  int insideOut) override;
97 
102  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
103  double x[3], double pcoords[3], int& subId) override;
104 
105 
109  int GetParametricCenter(double pcoords[3]) override;
110 
115  double GetParametricDistance(const double pcoords[3]) override;
116 
120  static void InterpolationFunctions(const double pcoords[3], double weights[7]);
124  static void InterpolationDerivs(const double pcoords[3], double derivs[14]);
126 
130  void InterpolateFunctions(const double pcoords[3], double weights[7]) override
131  {
133  }
134  void InterpolateDerivs(const double pcoords[3], double derivs[14]) override
135  {
137  }
139 
140 protected:
142  ~vtkBiQuadraticTriangle() override;
143 
146  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
147 
148 private:
150  void operator=(const vtkBiQuadraticTriangle&) = delete;
151 };
152 //----------------------------------------------------------------------------
153 inline int vtkBiQuadraticTriangle::GetParametricCenter(double pcoords[3])
154 {
155  pcoords[0] = pcoords[1] = 1./3;
156  pcoords[2] = 0.0;
157  return 0;
158 }
159 
160 
161 #endif
162 
163 
void InterpolateFunctions(const double pcoords[3], double weights[7]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
represent and manipulate point attribute data
Definition: vtkPointData.h:37
int GetNumberOfFaces() override
Implement the vtkCell API.
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
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.
abstract superclass for non-linear cells
vtkCell * GetFace(int) override
Implement the vtkCell API.
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.
int GetCellType() override
Implement the vtkCell API.
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...
dynamic, self-adjusting array of double
static void InterpolationFunctions(const double pcoords[3], double weights[7])
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:59
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
Implement the vtkCell API.
a simple class to control print indentation
Definition: vtkIndent.h:39
static void InterpolationDerivs(const double pcoords[3], double derivs[14])
list of point or cell ids
Definition: vtkIdList.h:36
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
cell represents a parabolic, isoparametric triangle
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.
void InterpolateDerivs(const double pcoords[3], double derivs[14]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic triangle in parametric coordinates.
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:50
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
int GetNumberOfEdges() override
Implement the vtkCell API.
cell represents a parabolic, isoparametric edge
a cell that represents a triangle
Definition: vtkTriangle.h:41
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.
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.
represent and manipulate 3D points
Definition: vtkPoints.h:39