VTK  9.0.1
vtkQuadraticPyramid.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticPyramid.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
39 #ifndef vtkQuadraticPyramid_h
40 #define vtkQuadraticPyramid_h
41 
42 #include "vtkCommonDataModelModule.h" // For export macro
43 #include "vtkNonLinearCell.h"
44 
45 class vtkQuadraticEdge;
46 class vtkQuadraticQuad;
48 class vtkTetra;
49 class vtkPyramid;
50 class vtkDoubleArray;
51 
52 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
53 {
54 public:
57  void PrintSelf(ostream& os, vtkIndent indent) override;
58 
60 
64  int GetCellType() override { return VTK_QUADRATIC_PYRAMID; }
65  int GetCellDimension() override { return 3; }
66  int GetNumberOfEdges() override { return 8; }
67  int GetNumberOfFaces() override { return 5; }
68  vtkCell* GetEdge(int edgeId) override;
69  vtkCell* GetFace(int faceId) override;
71 
72  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
73  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
74  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
75  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
76  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
77  double& dist2, double weights[]) override;
78  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
79  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
81  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
82  double* GetParametricCoords() override;
83 
89  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
90  vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
91  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
92 
97  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
98  double pcoords[3], int& subId) override;
99 
103  int GetParametricCenter(double pcoords[3]) override;
104 
105  static void InterpolationFunctions(const double pcoords[3], double weights[13]);
106  static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
108 
112  void InterpolateFunctions(const double pcoords[3], double weights[13]) override
113  {
115  }
116  void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
117  {
119  }
121 
122 
129  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
130  static const vtkIdType* GetFaceArray(vtkIdType faceId);
132 
138  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[39]);
139 
140 protected:
143 
152  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
153 
155 
161  void Subdivide(
162  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
164 
165 
171  void ResizeArrays(vtkIdType newSize);
173 
174 private:
175  vtkQuadraticPyramid(const vtkQuadraticPyramid&) = delete;
176  void operator=(const vtkQuadraticPyramid&) = delete;
177 };
178 //----------------------------------------------------------------------------
179 // Return the center of the quadratic pyramid in parametric coordinates.
180 //
181 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
182 {
183  pcoords[0] = pcoords[1] = 6.0 / 13.0;
184  pcoords[2] = 3.0 / 13.0;
185  return 0;
186 }
187 
188 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:34
vtkQuadraticPyramid::InterpolationDerivs
static void InterpolationDerivs(const double pcoords[3], double derivs[39])
vtkQuadraticPyramid::Triangulate
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
vtkQuadraticPyramid::Edge
vtkQuadraticEdge * Edge
Definition: vtkQuadraticPyramid.h:144
vtkQuadraticPyramid::InterpolationFunctions
static void InterpolationFunctions(const double pcoords[3], double weights[13])
vtkQuadraticPyramid::Derivatives
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:32
vtkX3D::value
@ value
Definition: vtkX3D.h:226
vtkIdType
int vtkIdType
Definition: vtkType.h:338
vtkQuadraticPyramid::CellData
vtkCellData * CellData
Definition: vtkQuadraticPyramid.h:150
vtkQuadraticPyramid::Clip
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.
vtkQuadraticPyramid::GetParametricCenter
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic pyramid in parametric coordinates.
Definition: vtkQuadraticPyramid.h:181
vtkQuadraticTriangle
cell represents a parabolic, isoparametric triangle
Definition: vtkQuadraticTriangle.h:44
vtkPyramid
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:44
vtkQuadraticPyramid::Scalars
vtkDoubleArray * Scalars
Definition: vtkQuadraticPyramid.h:152
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
vtkQuadraticPyramid::Subdivide
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...
vtkQuadraticQuad
cell represents a parabolic, 8-node isoparametric quad
Definition: vtkQuadraticQuad.h:44
vtkQuadraticPyramid::InterpolateFunctions
void InterpolateFunctions(const double pcoords[3], double weights[13]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkQuadraticPyramid.h:112
vtkQuadraticPyramid::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkQuadraticPyramid.h:66
VTK_QUADRATIC_PYRAMID
@ VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:72
vtkQuadraticPyramid
cell represents a parabolic, 13-node isoparametric pyramid
Definition: vtkQuadraticPyramid.h:53
vtkQuadraticPyramid::IntersectWithLine
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.
vtkQuadraticPyramid::GetEdgeArray
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkQuadraticPyramid::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkQuadraticPyramid.h:67
vtkQuadraticPyramid::InterpolateDerivs
void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
Definition: vtkQuadraticPyramid.h:116
vtkQuadraticPyramid::GetFaceArray
static const vtkIdType * GetFaceArray(vtkIdType faceId)
vtkQuadraticPyramid::EvaluateLocation
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:57
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
vtkQuadraticPyramid::GetParametricCoords
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkQuadraticPyramid::Face
vtkQuadraticQuad * Face
Definition: vtkQuadraticPyramid.h:146
vtkQuadraticPyramid::New
static vtkQuadraticPyramid * New()
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:34
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:180
vtkQuadraticPyramid::EvaluatePosition
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...
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:31
vtkQuadraticPyramid::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkQuadraticPyramid::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkQuadraticPyramid.h:65
vtkQuadraticPyramid::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkQuadraticPyramid.h:64
vtkQuadraticPyramid::PointData
vtkPointData * PointData
Definition: vtkQuadraticPyramid.h:149
vtkNonLinearCell.h
vtkNonLinearCell
abstract superclass for non-linear cells
Definition: vtkNonLinearCell.h:37
vtkQuadraticPyramid::Pyramid
vtkPyramid * Pyramid
Definition: vtkQuadraticPyramid.h:148
vtkQuadraticPyramid::vtkQuadraticPyramid
vtkQuadraticPyramid()
vtkCell::GetParametricCenter
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
vtkQuadraticPyramid::~vtkQuadraticPyramid
~vtkQuadraticPyramid() override
vtkQuadraticPyramid::Tetra
vtkTetra * Tetra
Definition: vtkQuadraticPyramid.h:147
vtkQuadraticPyramid::JacobianInverse
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[39])
Given parametric coordinates compute inverse Jacobian transformation matrix.
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:36
vtkQuadraticPyramid::CellScalars
vtkDoubleArray * CellScalars
Definition: vtkQuadraticPyramid.h:151
vtkQuadraticPyramid::Contour
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::GetFace
vtkCell * GetFace(int faceId) override
Return the face cell from the faceId of the cell.
vtkTetra
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:42
vtkQuadraticPyramid::CellBoundary
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...
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:41
vtkQuadraticPyramid::ResizeArrays
void ResizeArrays(vtkIdType newSize)
Resize the superclasses' member arrays to newSize where newSize should either be 13 or 14.
vtkQuadraticPyramid::GetEdge
vtkCell * GetEdge(int edgeId) override
Return the edge cell from the edgeId of the cell.
vtkQuadraticPyramid::TriangleFace
vtkQuadraticTriangle * TriangleFace
Definition: vtkQuadraticPyramid.h:145