VTK  9.2.6
vtkPyramid.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkPyramid.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=========================================================================*/
34#ifndef vtkPyramid_h
35#define vtkPyramid_h
36
37#include "vtkCell3D.h"
38#include "vtkCommonDataModelModule.h" // For export macro
39
40class vtkLine;
41class vtkQuad;
42class vtkTriangle;
45
46class VTKCOMMONDATAMODEL_EXPORT vtkPyramid : public vtkCell3D
47{
48public:
49 static vtkPyramid* New();
50 vtkTypeMacro(vtkPyramid, vtkCell3D);
51 void PrintSelf(ostream& os, vtkIndent indent) override;
52
54
57 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
58 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
59 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
60 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
61 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
62 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
63 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
64 bool GetCentroid(double centroid[3]) const override;
65 bool IsInsideOut() override;
67
71 static constexpr vtkIdType NumberOfPoints = 5;
72
76 static constexpr vtkIdType NumberOfEdges = 8;
77
81 static constexpr vtkIdType NumberOfFaces = 5;
82
87 static constexpr vtkIdType MaximumFaceSize = 4;
88
94 static constexpr vtkIdType MaximumValence = 4;
96
99 int GetCellType() override { return VTK_PYRAMID; }
100 int GetCellDimension() override { return 3; }
101 int GetNumberOfEdges() override { return 8; }
102 int GetNumberOfFaces() override { return 5; }
103 vtkCell* GetEdge(int edgeId) override;
104 vtkCell* GetFace(int faceId) override;
105 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
106 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
107 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
108 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
109 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
110 double& dist2, double weights[]) override;
111 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
112 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
113 double pcoords[3], int& subId) override;
114 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
116 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
117 double* GetParametricCoords() override;
119
127 static int* GetTriangleCases(int caseId);
128
132 int GetParametricCenter(double pcoords[3]) override;
133
134 static void InterpolationFunctions(const double pcoords[3], double weights[5]);
135 static void InterpolationDerivs(const double pcoords[3], double derivs[15]);
137
141 void InterpolateFunctions(const double pcoords[3], double weights[5]) override
142 {
143 vtkPyramid::InterpolationFunctions(pcoords, weights);
144 }
145 void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
146 {
147 vtkPyramid::InterpolationDerivs(pcoords, derivs);
148 }
150
151 int JacobianInverse(const double pcoords[3], double** inverse, double derivs[15]);
152
154
165
170
175
180
185
190
194 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
195
196protected:
198 ~vtkPyramid() override;
199
203
204private:
205 vtkPyramid(const vtkPyramid&) = delete;
206 void operator=(const vtkPyramid&) = delete;
207};
208
209//----------------------------------------------------------------------------
210inline int vtkPyramid::GetParametricCenter(double pcoords[3])
211{
212 pcoords[0] = pcoords[1] = 0.4;
213 pcoords[2] = 0.2;
214 return 0;
215}
216
217#endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
object to represent cell connectivity
Definition: vtkCellArray.h:187
represent and manipulate cell attribute data
Definition: vtkCellData.h:42
abstract class to specify cell behavior
Definition: vtkCell.h:61
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:56
list of point or cell ids
Definition: vtkIdList.h:34
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:40
cell represents a 1D line
Definition: vtkLine.h:34
represent and manipulate point attribute data
Definition: vtkPointData.h:42
represent and manipulate 3D points
Definition: vtkPoints.h:40
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:47
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
Definition: vtkPyramid.h:210
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
static vtkPyramid * New()
vtkQuad * Quad
Definition: vtkPyramid.h:202
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.
vtkTriangle * Triangle
Definition: vtkPyramid.h:201
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 bool ComputeCentroid(vtkPoints *points, const vtkIdType *pointIds, double centroid[3])
Static method version of GetCentroid.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
See the vtkCell API for descriptions of these methods.
static const vtkIdType * GetPointToIncidentEdgesArray(vtkIdType pointId)
Static method version of GetPointToIncidentEdgesArray.
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
~vtkPyramid() override
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkPyramid.h:145
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPyramid.h:101
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.
vtkLine * Line
Definition: vtkPyramid.h:200
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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:100
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 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:102
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPyramid.h:99
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])
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
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.
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
void InterpolateFunctions(const double pcoords[3], double weights[5]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkPyramid.h:141
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:39
a cell that represents a triangle
Definition: vtkTriangle.h:39
dataset represents arbitrary combinations of all possible cell types
@ VTK_PYRAMID
Definition: vtkCellType.h:60
int vtkIdType
Definition: vtkType.h:332
#define VTK_SIZEHINT(...)