VTK  9.3.0
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
39 #ifndef vtkQuadraticPyramid_h
40 #define vtkQuadraticPyramid_h
41 
42 #include "vtkCommonDataModelModule.h" // For export macro
43 #include "vtkNonLinearCell.h"
44 
45 VTK_ABI_NAMESPACE_BEGIN
46 class vtkQuadraticEdge;
47 class vtkQuadraticQuad;
49 class vtkTetra;
50 class vtkPyramid;
51 class vtkDoubleArray;
52 
53 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
54 {
55 public:
58  void PrintSelf(ostream& os, vtkIndent indent) override;
59 
61 
65  int GetCellType() override { return VTK_QUADRATIC_PYRAMID; }
66  int GetCellDimension() override { return 3; }
67  int GetNumberOfEdges() override { return 8; }
68  int GetNumberOfFaces() override { return 5; }
69  vtkCell* GetEdge(int edgeId) override;
70  vtkCell* GetFace(int faceId) override;
72 
73  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
74  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
75  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
76  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
77  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
78  double& dist2, double weights[]) override;
79  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
80  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
82  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
83  double* GetParametricCoords() override;
84 
90  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
91  vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
92  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
93 
98  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
99  double pcoords[3], int& subId) override;
100 
104  int GetParametricCenter(double pcoords[3]) override;
105 
106  static void InterpolationFunctions(const double pcoords[3], double weights[13]);
107  static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
109 
113  void InterpolateFunctions(const double pcoords[3], double weights[13]) override
114  {
116  }
117  void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
118  {
120  }
123 
130  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
131  static const vtkIdType* GetFaceArray(vtkIdType faceId);
133 
139  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[39]);
140 
141 protected:
144 
153  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
154 
156 
162  void Subdivide(
163  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
166 
172  void ResizeArrays(vtkIdType newSize);
174 
175 private:
176  vtkQuadraticPyramid(const vtkQuadraticPyramid&) = delete;
177  void operator=(const vtkQuadraticPyramid&) = delete;
178 };
179 //----------------------------------------------------------------------------
180 // Return the center of the quadratic pyramid in parametric coordinates.
181 //
182 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
183 {
184  pcoords[0] = pcoords[1] = 6.0 / 13.0;
185  pcoords[2] = 3.0 / 13.0;
186  return 0;
187 }
188 
189 VTK_ABI_NAMESPACE_END
190 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:185
represent and manipulate cell attribute data
Definition: vtkCellData.h:40
abstract class to specify cell behavior
Definition: vtkCell.h:59
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:38
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:39
represent and manipulate 3D points
Definition: vtkPoints.h:38
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:45
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 13-node isoparametric pyramid
int GetCellType() override
Implement the vtkCell API.
vtkDoubleArray * Scalars
vtkQuadraticQuad * Face
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCell * GetFace(int faceId) 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’).
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
vtkQuadraticTriangle * TriangleFace
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
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
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])
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 vtkQuadraticPyramid * New()
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.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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:43
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:72
int vtkIdType
Definition: vtkType.h:315