41 #ifndef vtkBiQuadraticTriangle_h 42 #define vtkBiQuadraticTriangle_h 44 #include "vtkCommonDataModelModule.h" 78 int& subId,
double pcoords[3],
79 double& dist2,
double weights[])
override;
81 double *weights)
override;
83 void Derivatives(
int subId,
const double pcoords[3],
const double *values,
84 int dim,
double *derivs)
override;
96 int insideOut)
override;
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;
120 static void InterpolationFunctions(
const double pcoords[3],
double weights[7]);
124 static void InterpolationDerivs(
const double pcoords[3],
double derivs[14]);
155 pcoords[0] = pcoords[1] = 1./3;
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
int GetNumberOfFaces() override
Implement the vtkCell API.
represent and manipulate cell attribute data
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.
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
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
static void InterpolationDerivs(const double pcoords[3], double derivs[14])
list of point or cell ids
abstract superclass for arrays of numeric data
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
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
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