VTK  9.6.1
vtkQuadraticPyramid.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
29
30#ifndef vtkQuadraticPyramid_h
31#define vtkQuadraticPyramid_h
32
33#include "vtkCommonDataModelModule.h" // For export macro
34#include "vtkNonLinearCell.h"
36VTK_ABI_NAMESPACE_BEGIN
40class vtkTetra;
41class vtkPyramid;
42class vtkDoubleArray;
44class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
45{
46public:
49 void PrintSelf(ostream& os, vtkIndent indent) override;
50
56 int GetCellType() override { return VTK_QUADRATIC_PYRAMID; }
57 int GetCellDimension() override { return 3; }
58 int GetNumberOfEdges() override { return 8; }
59 int GetNumberOfFaces() override { return 5; }
60 vtkCell* GetEdge(int edgeId) override;
61 vtkCell* GetFace(int faceId) override;
62 ///@}
63
64 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
65 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
66 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
67 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
68 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
69 double& dist2, double weights[]) override;
70 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
71 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
73 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
74 double* GetParametricCoords() override;
75
81 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
82 vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
83 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
84
89 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
90 double pcoords[3], int& subId) override;
91
95 int GetParametricCenter(double pcoords[3]) override;
96
97 static void InterpolationFunctions(const double pcoords[3], double weights[13]);
98 static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
100
104 void InterpolateFunctions(const double pcoords[3], double weights[13]) override
105 {
107 }
108 void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
109 {
111 }
112
114
121 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
122 static const vtkIdType* GetFaceArray(vtkIdType faceId);
124
130 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[39]);
131
132protected:
135
144 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
145
147
153 void Subdivide(
154 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
157
163 void ResizeArrays(vtkIdType newSize);
165
166private:
168 void operator=(const vtkQuadraticPyramid&) = delete;
169};
170//----------------------------------------------------------------------------
171// Return the center of the quadratic pyramid in parametric coordinates.
172//
173inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
174{
175 pcoords[0] = pcoords[1] = 6.0 / 13.0;
176 pcoords[2] = 3.0 / 13.0;
177 return 0;
178}
179
180VTK_ABI_NAMESPACE_END
181#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:32
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:24
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:29
represent and manipulate point attribute data
a 3D cell that represents a linear pyramid
Definition vtkPyramid.h:36
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 13-node isoparametric pyramid
int GetCellType() override
Implement the vtkCell API.
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic pyramid in parametric coordinates.
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.
int GetCellDimension() override
Implement the vtkCell API.
static vtkQuadraticPyramid * New()
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.
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
static void InterpolationDerivs(const double pcoords[3], double derivs[39])
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
void InterpolateFunctions(const double pcoords[3], double weights[13]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
int GetNumberOfFaces() override
Implement the vtkCell API.
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...
static void InterpolationFunctions(const double pcoords[3], double weights[13])
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
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 PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
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 GetNumberOfEdges() override
Implement the vtkCell API.
cell represents a parabolic, 8-node isoparametric quad
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:34
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
vtkQuadraticEdge * Edge
vtkPointData * PointData
vtkCellData * CellData
vtkBiQuadraticQuadraticHexahedron vtkNonLinearCell JacobianInverse(const double pcoords[3], double **inverse, double derivs[72])
Given parametric coordinates compute inverse Jacobian transformation matrix.
vtkDoubleArray * Scalars
vtkQuadraticQuad * Face
vtkDoubleArray * CellScalars
vtkQuadraticTriangle * TriangleFace
@ VTK_QUADRATIC_PYRAMID
Definition vtkCellType.h:62
#define vtkDataArray
vtkPyramid * Pyramid
vtkTetra * Tetra
vtkQuadraticPyramid()
~vtkQuadraticPyramid() override
int vtkIdType
Definition vtkType.h:368