VTK  9.2.5
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 =========================================================================*/
93 #ifndef vtkPyramid_h
94 #define vtkPyramid_h
95 
96 #include "vtkCell3D.h"
97 #include "vtkCommonDataModelModule.h" // For export macro
98 
99 class vtkLine;
100 class vtkQuad;
101 class vtkTriangle;
102 class vtkUnstructuredGrid;
104 
105 class VTKCOMMONDATAMODEL_EXPORT vtkPyramid : public vtkCell3D
106 {
107 public:
108  static vtkPyramid* New();
109  vtkTypeMacro(vtkPyramid, vtkCell3D);
110  void PrintSelf(ostream& os, vtkIndent indent) override;
111 
113 
116  void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
117  vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
118  void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
119  vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
120  vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
121  vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
122  vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
123  bool GetCentroid(double centroid[3]) const override;
124  bool IsInsideOut() override;
126 
130  static constexpr vtkIdType NumberOfPoints = 5;
131 
135  static constexpr vtkIdType NumberOfEdges = 8;
136 
140  static constexpr vtkIdType NumberOfFaces = 5;
141 
146  static constexpr vtkIdType MaximumFaceSize = 4;
147 
153  static constexpr vtkIdType MaximumValence = 4;
155 
158  int GetCellType() override { return VTK_PYRAMID; }
159  int GetCellDimension() override { return 3; }
160  int GetNumberOfEdges() override { return 8; }
161  int GetNumberOfFaces() override { return 5; }
162  vtkCell* GetEdge(int edgeId) override;
163  vtkCell* GetFace(int faceId) override;
164  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
165  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
166  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
167  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
168  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
169  double& dist2, double weights[]) override;
170  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
171  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
172  double pcoords[3], int& subId) override;
173  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
175  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
176  double* GetParametricCoords() override;
178 
186  static int* GetTriangleCases(int caseId);
187 
191  int GetParametricCenter(double pcoords[3]) override;
192 
193  static void InterpolationFunctions(const double pcoords[3], double weights[5]);
194  static void InterpolationDerivs(const double pcoords[3], double derivs[15]);
196 
200  void InterpolateFunctions(const double pcoords[3], double weights[5]) override
201  {
202  vtkPyramid::InterpolationFunctions(pcoords, weights);
203  }
204  void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
205  {
206  vtkPyramid::InterpolationDerivs(pcoords, derivs);
207  }
209 
210  int JacobianInverse(const double pcoords[3], double** inverse, double derivs[15]);
211 
213 
221  static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
222  static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(4);
224 
229 
234 
239 
244 
249 
253  static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
254 
255 protected:
257  ~vtkPyramid() override;
258 
262 
263 private:
264  vtkPyramid(const vtkPyramid&) = delete;
265  void operator=(const vtkPyramid&) = delete;
266 };
267 
268 //----------------------------------------------------------------------------
269 inline int vtkPyramid::GetParametricCenter(double pcoords[3])
270 {
271  pcoords[0] = pcoords[1] = 0.4;
272  pcoords[2] = 0.2;
273  return 0;
274 }
275 
276 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
object to represent cell connectivity
Definition: vtkCellArray.h:296
represent and manipulate cell attribute data
Definition: vtkCellData.h:151
abstract class to specify cell behavior
Definition: vtkCell.h:150
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
list of point or cell ids
Definition: vtkIdList.h:143
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:119
cell represents a 1D line
Definition: vtkLine.h:143
represent and manipulate point attribute data
Definition: vtkPointData.h:151
represent and manipulate 3D points
Definition: vtkPoints.h:149
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:106
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
Definition: vtkPyramid.h:269
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
static vtkPyramid * New()
vtkQuad * Quad
Definition: vtkPyramid.h:261
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:260
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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.
double * GetParametricCoords() 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.
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) 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).
~vtkPyramid() override
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkPyramid.h:204
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPyramid.h:160
vtkCell * GetFace(int faceId) 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:259
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
int JacobianInverse(const double pcoords[3], double **inverse, double derivs[15])
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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:159
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
See vtkCell3D API for description 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:161
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPyramid.h:158
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 void InterpolationDerivs(const double pcoords[3], double derivs[15])
static const vtkIdType * GetPointToIncidentEdgesArray(vtkIdType pointId)
Static method version of GetPointToIncidentEdgesArray.
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
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:200
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:98
a cell that represents a triangle
Definition: vtkTriangle.h:148
dataset represents arbitrary combinations of all possible cell types
@ points
Definition: vtkX3D.h:452
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_PYRAMID
Definition: vtkCellType.h:99
int vtkIdType
Definition: vtkType.h:332
#define VTK_SIZEHINT(...)