VTK  9.1.0
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=========================================================================*/
58#ifndef vtkQuadraticPyramid_h
59#define vtkQuadraticPyramid_h
60
61#include "vtkCommonDataModelModule.h" // For export macro
62#include "vtkNonLinearCell.h"
63
67class vtkTetra;
68class vtkPyramid;
69class vtkDoubleArray;
70
71class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
72{
73public:
76 void PrintSelf(ostream& os, vtkIndent indent) override;
77
79
83 int GetCellType() override { return VTK_QUADRATIC_PYRAMID; }
84 int GetCellDimension() override { return 3; }
85 int GetNumberOfEdges() override { return 8; }
86 int GetNumberOfFaces() override { return 5; }
87 vtkCell* GetEdge(int edgeId) override;
88 vtkCell* GetFace(int faceId) override;
90
91 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
92 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
93 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
94 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
95 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
96 double& dist2, double weights[]) override;
97 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
98 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
100 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
101 double* GetParametricCoords() override;
102
108 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
109 vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
110 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
111
116 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
117 double pcoords[3], int& subId) override;
118
122 int GetParametricCenter(double pcoords[3]) override;
123
124 static void InterpolationFunctions(const double pcoords[3], double weights[13]);
125 static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
127
131 void InterpolateFunctions(const double pcoords[3], double weights[13]) override
132 {
134 }
135 void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
136 {
138 }
141
148 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
149 static const vtkIdType* GetFaceArray(vtkIdType faceId);
151
157 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[39]);
158
159protected:
162
171 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
172
174
181 vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
184
190 void ResizeArrays(vtkIdType newSize);
192
193private:
195 void operator=(const vtkQuadraticPyramid&) = delete;
196};
197//----------------------------------------------------------------------------
198// Return the center of the quadratic pyramid in parametric coordinates.
199//
200inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
201{
202 pcoords[0] = pcoords[1] = 6.0 / 13.0;
203 pcoords[2] = 3.0 / 13.0;
204 return 0;
205}
206
207#endif
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:140
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:113
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:103
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
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’).
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
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
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:111
@ VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:111
int vtkIdType
Definition: vtkType.h:332