VTK  9.6.1
vtkPyramid.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
21
22#ifndef vtkPyramid_h
23#define vtkPyramid_h
24
25#include "vtkCell3D.h"
26#include "vtkCommonDataModelModule.h" // For export macro
27
28VTK_ABI_NAMESPACE_BEGIN
29class vtkLine;
30class vtkQuad;
31class vtkTriangle;
34
35class VTKCOMMONDATAMODEL_EXPORT vtkPyramid : public vtkCell3D
36{
37public:
38 static vtkPyramid* New();
39 vtkTypeMacro(vtkPyramid, vtkCell3D);
40 void PrintSelf(ostream& os, vtkIndent indent) override;
41
43
46 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
47 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
48 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
49 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
50 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
51 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
52 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
53 bool GetCentroid(double centroid[3]) const override;
54 bool IsInsideOut() override;
60 static constexpr vtkIdType NumberOfPoints = 5;
65 static constexpr vtkIdType NumberOfEdges = 8;
70 static constexpr vtkIdType NumberOfFaces = 5;
71
76 static constexpr vtkIdType MaximumFaceSize = 4;
77
83 static constexpr vtkIdType MaximumValence = 4;
85
88 int GetCellType() override { return VTK_PYRAMID; }
89 int GetCellDimension() override { return 3; }
90 int GetNumberOfEdges() override { return 8; }
91 int GetNumberOfFaces() override { return 5; }
92 vtkCell* GetEdge(int edgeId) override;
93 vtkCell* GetFace(int faceId) override;
94 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
95 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
96 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
97 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
98 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
99 double& dist2, double weights[]) override;
100 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
101 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
102 double pcoords[3], int& subId) override;
103 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
105 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
106 double* GetParametricCoords() override;
108
116 static int* GetTriangleCases(int caseId);
117
121 int GetParametricCenter(double pcoords[3]) override;
122
123 static void InterpolationFunctions(const double pcoords[3], double weights[5]);
124 static void InterpolationDerivs(const double pcoords[3], double derivs[15]);
126
130 void InterpolateFunctions(const double pcoords[3], double weights[5]) override
131 {
132 vtkPyramid::InterpolationFunctions(pcoords, weights);
133 }
134 void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
135 {
136 vtkPyramid::InterpolationDerivs(pcoords, derivs);
137 }
138
139
140 int JacobianInverse(const double pcoords[3], double** inverse, double derivs[15]);
141
143
151 static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
152 static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(4);
154
159
164
169
174
179
183 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
184
185protected:
187 ~vtkPyramid() override;
188
192
193private:
194 vtkPyramid(const vtkPyramid&) = delete;
195 void operator=(const vtkPyramid&) = delete;
196};
197
198//----------------------------------------------------------------------------
199inline int vtkPyramid::GetParametricCenter(double pcoords[3])
200{
201 pcoords[0] = pcoords[1] = 0.4;
202 pcoords[2] = 0.2;
203 return 0;
204}
205
206VTK_ABI_NAMESPACE_END
207#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:32
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
list of point or cell ids
Definition vtkIdList.h:24
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:29
cell represents a 1D line
Definition vtkLine.h:23
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:30
a 3D cell that represents a linear pyramid
Definition vtkPyramid.h:36
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
Definition vtkPyramid.h:199
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
static vtkPyramid * New()
vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
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
See the vtkCell API for descriptions of these methods.
void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
See the vtkCell API for descriptions of these methods.
bool IsInsideOut() override
See vtkCell3D API for description of these methods.
static constexpr vtkIdType NumberOfPoints
static constexpr handle on the number of points.
Definition vtkPyramid.h:60
static constexpr vtkIdType NumberOfEdges
static contexpr handle on the number of faces.
Definition vtkPyramid.h:65
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
Definition vtkPyramid.h:134
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition vtkPyramid.h:90
double * GetParametricCoords() override
See the vtkCell API for descriptions of these methods.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static constexpr vtkIdType NumberOfFaces
static contexpr handle on the number of edges.
Definition vtkPyramid.h:70
int JacobianInverse(const double pcoords[3], double **inverse, double derivs[15])
void GetEdgePoints(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition vtkPyramid.h:89
static constexpr vtkIdType MaximumFaceSize
static contexpr handle on the maximum face size.
Definition vtkPyramid.h:76
vtkCell * GetFace(int faceId) override
See the vtkCell API for descriptions of these methods.
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
See vtkCell3D API for description of these methods.
static int * GetTriangleCases(int caseId)
Return the case table for table-based isocontouring (aka marching cubes style implementations).
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
See the vtkCell API for descriptions of these methods.
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
See the vtkCell API for descriptions of these methods.
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition vtkPyramid.h:91
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition vtkPyramid.h:88
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
See the vtkCell API for descriptions of these methods.
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
bool GetCentroid(double centroid[3]) const override
See vtkCell3D API for description of these methods.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
See the vtkCell API for descriptions of these methods.
static void InterpolationFunctions(const double pcoords[3], double weights[5])
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
static constexpr vtkIdType MaximumValence
static constexpr handle on the maximum valence of this cell.
Definition vtkPyramid.h:83
static void InterpolationDerivs(const double pcoords[3], double derivs[15])
vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
void InterpolateFunctions(const double pcoords[3], double weights[5]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives).
Definition vtkPyramid.h:130
a cell that represents a 2D quadrilateral
Definition vtkQuad.h:28
a cell that represents a triangle
Definition vtkTriangle.h:28
dataset represents arbitrary combinations of all possible cell types
@ VTK_PYRAMID
Definition vtkCellType.h:50
#define vtkDataArray
vtkQuad * Quad
vtkLine * Line
vtkTriangle * Triangle
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
static bool ComputeCentroid(vtkPoints *points, const vtkIdType *pointIds, double centroid[3])
Static method version of GetCentroid.
vtkHexagonalPrism vtkCell3D GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
static const vtkIdType * GetPointToIncidentEdgesArray(vtkIdType pointId)
Static method version of GetPointToIncidentEdgesArray.
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
~vtkPyramid() override
vtkPyramid()
int vtkIdType
Definition vtkType.h:368
#define VTK_SIZEHINT(...)