VTK  9.3.20240424
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
49#ifndef vtkQuadraticPyramid_h
50#define vtkQuadraticPyramid_h
51
52#include "vtkCommonDataModelModule.h" // For export macro
53#include "vtkNonLinearCell.h"
54
55VTK_ABI_NAMESPACE_BEGIN
59class vtkTetra;
60class vtkPyramid;
61class vtkDoubleArray;
62
63class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
64{
65public:
68 void PrintSelf(ostream& os, vtkIndent indent) override;
69
71
75 int GetCellType() override { return VTK_QUADRATIC_PYRAMID; }
76 int GetCellDimension() override { return 3; }
77 int GetNumberOfEdges() override { return 8; }
78 int GetNumberOfFaces() override { return 5; }
79 vtkCell* GetEdge(int edgeId) override;
80 vtkCell* GetFace(int faceId) override;
82
83 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
84 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
85 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
86 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
87 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
88 double& dist2, double weights[]) override;
89 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
90 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
92 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
93 double* GetParametricCoords() override;
94
100 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
101 vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
102 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
103
108 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
109 double pcoords[3], int& subId) override;
110
114 int GetParametricCenter(double pcoords[3]) override;
115
116 static void InterpolationFunctions(const double pcoords[3], double weights[13]);
117 static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
119
123 void InterpolateFunctions(const double pcoords[3], double weights[13]) override
124 {
126 }
127 void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
128 {
130 }
133
140 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
141 static const vtkIdType* GetFaceArray(vtkIdType faceId);
143
149 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[39]);
150
151protected:
154
163 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
164
166
173 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
176
182 void ResizeArrays(vtkIdType newSize);
184
185private:
187 void operator=(const vtkQuadraticPyramid&) = delete;
188};
189//----------------------------------------------------------------------------
190// Return the center of the quadratic pyramid in parametric coordinates.
191//
192inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
193{
194 pcoords[0] = pcoords[1] = 6.0 / 13.0;
195 pcoords[2] = 3.0 / 13.0;
196 return 0;
197}
198
199VTK_ABI_NAMESPACE_END
200#endif
object to represent cell connectivity
represent and manipulate cell attribute data
abstract class to specify cell behavior
Definition vtkCell.h:130
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:133
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
abstract superclass for non-linear cells
represent and manipulate point attribute data
a 3D cell that represents a linear pyramid
Definition vtkPyramid.h:95
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 13-node isoparametric pyramid
int GetCellType() override
Implement the vtkCell API.
vtkQuadraticQuad * Face
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
void ResizeArrays(vtkIdType newSize)
Resize the superclasses' member arrays to newSize where newSize should either be 13 or 14.
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 const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static vtkQuadraticPyramid * New()
vtkQuadraticTriangle * TriangleFace
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.
~vtkQuadraticPyramid() override
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)
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[39])
Given parametric coordinates compute inverse Jacobian transformation matrix.
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.
vtkQuadraticEdge * Edge
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)
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
This method adds in a point at the center of the quadrilateral face and then interpolates values to t...
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.
vtkDoubleArray * CellScalars
cell represents a parabolic, 8-node isoparametric quad
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:113
@ VTK_QUADRATIC_PYRAMID
Definition vtkCellType.h:82
int vtkIdType
Definition vtkType.h:315