VTK  9.3.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
30 #ifndef vtkQuadraticPyramid_h
31 #define vtkQuadraticPyramid_h
32 
33 #include "vtkCommonDataModelModule.h" // For export macro
34 #include "vtkNonLinearCell.h"
35 
36 VTK_ABI_NAMESPACE_BEGIN
37 class vtkQuadraticEdge;
38 class vtkQuadraticQuad;
40 class vtkTetra;
41 class vtkPyramid;
42 class vtkDoubleArray;
43 
44 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
45 {
46 public:
47  static vtkQuadraticPyramid* New();
49  void PrintSelf(ostream& os, vtkIndent indent) override;
50 
52 
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;
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 Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
72  void Derivatives(
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  }
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 
132 protected:
134  ~vtkQuadraticPyramid() override;
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 
166 private:
167  vtkQuadraticPyramid(const vtkQuadraticPyramid&) = delete;
168  void operator=(const vtkQuadraticPyramid&) = delete;
169 };
170 //----------------------------------------------------------------------------
171 // Return the center of the quadratic pyramid in parametric coordinates.
172 //
173 inline 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 
180 VTK_ABI_NAMESPACE_END
181 #endif
int GetNumberOfFaces() override
Implement the vtkCell API.
vtkDoubleArray * CellScalars
represent and manipulate point attribute data
Definition: vtkPointData.h:29
cell represents a parabolic, 13-node isoparametric pyramid
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:35
represent and manipulate cell attribute data
Definition: vtkCellData.h:30
int GetCellType() override
Implement the vtkCell API.
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.
vtkQuadraticEdge * Edge
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:315
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic pyramid in parametric coordinates.
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
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:49
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.
cell represents a parabolic, 8-node isoparametric quad
a simple class to control print indentation
Definition: vtkIndent.h:28
vtkQuadraticTriangle * TriangleFace
list of point or cell ids
Definition: vtkIdList.h:22
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:44
void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
static void InterpolationDerivs(const double pcoords[3], double derivs[39])
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.
static void InterpolationFunctions(const double pcoords[3], double weights[13])
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
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
Definition: vtkCellArray.h:175
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkQuadraticQuad * Face
cell represents a parabolic, isoparametric edge
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.
int GetNumberOfEdges() override
Implement the vtkCell API.
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.
cell represents a parabolic, isoparametric triangle
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
vtkDoubleArray * Scalars
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.
void InterpolateFunctions(const double pcoords[3], double weights[13]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
represent and manipulate 3D points
Definition: vtkPoints.h:28