 |
VTK
9.1.0
|
Go to the documentation of this file.
90 #ifndef vtkTriQuadraticPyramid_h
91 #define vtkTriQuadraticPyramid_h
93 #include "vtkCommonDataModelModule.h"
128 int EvaluatePosition(
const double x[3],
double closestPoint[3],
int& subId,
double pcoords[3],
129 double& dist2,
double weights[])
override;
130 void EvaluateLocation(
int& subId,
const double pcoords[3],
double x[3],
double* weights)
override;
136 int IntersectWithLine(
const double p1[3],
const double p2[3],
double tol,
double& t,
double x[3],
137 double pcoords[3],
int& subId)
override;
141 int subId,
const double pcoords[3],
const double* values,
int dim,
double* derivs)
override;
221 pcoords[0] = pcoords[1] = 0.5;
225 pcoords[2] = 283.0 / 456.0;
represent and manipulate 3D points
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
represent and manipulate point attribute data
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[57])
Given parametric coordinates compute inverse Jacobian transformation matrix.
int GetParametricCenter(double pcoords[3]) override
Return the center of the tri-quadratic pyramid in parametric coordinates.
cell represents a parabolic, isoparametric triangle
void InterpolateFunctions(const double pcoords[3], double weights[19]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
a 3D cell that represents a linear pyramid
abstract superclass for arrays of numeric data
static void InterpolationDerivs(const double pcoords[3], double derivs[57])
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tets, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this quadratic triangle using scalar value provided.
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
void InterpolateDerivs(const double pcoords[3], double derivs[57]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
~vtkTriQuadraticPyramid() override
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
Generate contouring primitives.
cell represents a parabolic, 9-node isoparametric quad
static vtkTriQuadraticPyramid * New()
abstract class to specify cell behavior
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
represent and manipulate cell attribute data
vtkNew< vtkBiQuadraticTriangle > TriangleFace2
a simple class to control print indentation
object to represent cell connectivity
Abstract class in support of both point location and point insertion.
list of point or cell ids
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Line-edge intersection.
int GetCellDimension() override
Implement the vtkCell API.
int GetCellType() override
Implement the vtkCell API.
vtkNew< vtkBiQuadraticTriangle > TriangleFace
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
abstract superclass for non-linear cells
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetNumberOfFaces() override
Implement the vtkCell API.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkNew< vtkDoubleArray > Scalars
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
cell represents a parabolic, 13-node isoparametric pyramid
dynamic, self-adjusting array of double
@ VTK_TRIQUADRATIC_PYRAMID
vtkNew< vtkQuadraticEdge > Edge
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
a 3D cell that represents a tetrahedron
int GetNumberOfEdges() override
Implement the vtkCell API.
cell represents a parabolic, isoparametric edge
static void InterpolationFunctions(const double pcoords[3], double weights[19])
vtkNew< vtkBiQuadraticQuad > QuadFace
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
vtkNew< vtkPyramid > Pyramid