VTK
vtkQuadraticTriangle.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticTriangle.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 vtkQuadraticTriangle_h
34 #define vtkQuadraticTriangle_h
35 
36 #include "vtkCommonDataModelModule.h" // For export macro
37 #include "vtkNonLinearCell.h"
38 
39 class vtkQuadraticEdge;
40 class vtkTriangle;
41 class vtkDoubleArray;
42 
43 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticTriangle : public vtkNonLinearCell
44 {
45 public:
46  static vtkQuadraticTriangle *New();
48  void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
49 
51 
55  int GetCellType() VTK_OVERRIDE {return VTK_QUADRATIC_TRIANGLE;};
56  int GetCellDimension() VTK_OVERRIDE {return 2;}
57  int GetNumberOfEdges() VTK_OVERRIDE {return 3;}
58  int GetNumberOfFaces() VTK_OVERRIDE {return 0;}
59  vtkCell *GetEdge(int edgeId) VTK_OVERRIDE;
60  vtkCell *GetFace(int) VTK_OVERRIDE {return 0;}
62 
63  int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) VTK_OVERRIDE;
64  void Contour(double value, vtkDataArray *cellScalars,
66  vtkCellArray *lines, vtkCellArray *polys,
67  vtkPointData *inPd, vtkPointData *outPd,
68  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) VTK_OVERRIDE;
69  int EvaluatePosition(double x[3], double* closestPoint,
70  int& subId, double pcoords[3],
71  double& dist2, double *weights) VTK_OVERRIDE;
72  void EvaluateLocation(int& subId, double pcoords[3], double x[3],
73  double *weights) VTK_OVERRIDE;
74  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) VTK_OVERRIDE;
75  void Derivatives(int subId, double pcoords[3], double *values,
76  int dim, double *derivs) VTK_OVERRIDE;
77  double *GetParametricCoords() VTK_OVERRIDE;
78 
84  void Clip(double value, vtkDataArray *cellScalars,
86  vtkPointData *inPd, vtkPointData *outPd,
87  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
88  int insideOut) VTK_OVERRIDE;
89 
94  int IntersectWithLine(double p1[3], double p2[3], double tol, double& t,
95  double x[3], double pcoords[3], int& subId) VTK_OVERRIDE;
96 
97 
101  int GetParametricCenter(double pcoords[3]) VTK_OVERRIDE;
102 
107  double GetParametricDistance(double pcoords[3]) VTK_OVERRIDE;
108 
112  static void InterpolationFunctions(double pcoords[3], double weights[6]);
116  static void InterpolationDerivs(double pcoords[3], double derivs[12]);
118 
122  void InterpolateFunctions(double pcoords[3], double weights[6]) VTK_OVERRIDE
123  {
125  }
126  void InterpolateDerivs(double pcoords[3], double derivs[12]) VTK_OVERRIDE
127  {
129  }
131 
132 protected:
134  ~vtkQuadraticTriangle() VTK_OVERRIDE;
135 
138  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
139 
140 private:
141  vtkQuadraticTriangle(const vtkQuadraticTriangle&) VTK_DELETE_FUNCTION;
142  void operator=(const vtkQuadraticTriangle&) VTK_DELETE_FUNCTION;
143 };
144 //----------------------------------------------------------------------------
145 inline int vtkQuadraticTriangle::GetParametricCenter(double pcoords[3])
146 {
147  pcoords[0] = pcoords[1] = 1./3;
148  pcoords[2] = 0.0;
149  return 0;
150 }
151 
152 
153 #endif
154 
155 
vtkQuadraticEdge * Edge
represent and manipulate point attribute data
Definition: vtkPointData.h:31
void InterpolateDerivs(double pcoords[3], double derivs[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
vtkCell * GetFace(int) override
Implement the vtkCell API.
virtual double * GetParametricCoords()
Return a contiguous array of parametric coordinates of the points defining this cell.
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
int GetNumberOfEdges() override
Implement the vtkCell API.
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 EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
virtual int EvaluatePosition(double x[3], double *closestPoint, 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...
abstract superclass for non-linear cells
int GetCellDimension() override
Implement the vtkCell API.
int vtkIdType
Definition: vtkType.h:345
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
dynamic, self-adjusting array of double
abstract class to specify cell behavior
Definition: vtkCell.h:56
virtual double GetParametricDistance(double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
a simple class to control print indentation
Definition: vtkIndent.h:33
list of point or cell ids
Definition: vtkIdList.h:30
int GetNumberOfFaces() override
Implement the vtkCell API.
virtual void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
virtual int IntersectWithLine(double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
void InterpolateFunctions(double pcoords[3], double weights[6]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
static void InterpolationFunctions(double pcoords[3], double weights[6])
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.
object to represent cell connectivity
Definition: vtkCellArray.h:44
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
static void InterpolationDerivs(double pcoords[3], double derivs[12])
cell represents a parabolic, isoparametric edge
a cell that represents a triangle
Definition: vtkTriangle.h:35
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.
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic triangle in parametric coordinates.
cell represents a parabolic, isoparametric triangle
int GetCellType() override
Implement the vtkCell API.
virtual int CellBoundary(int subId, double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
represent and manipulate 3D points
Definition: vtkPoints.h:33