VTK  9.0.2
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 =========================================================================*/
31 #ifndef vtkPyramid_h
32 #define vtkPyramid_h
33 
34 #include "vtkCell3D.h"
35 #include "vtkCommonDataModelModule.h" // For export macro
36 
37 class vtkLine;
38 class vtkQuad;
39 class vtkTriangle;
42 
43 class VTKCOMMONDATAMODEL_EXPORT vtkPyramid : public vtkCell3D
44 {
45 public:
46  static vtkPyramid* New();
47  vtkTypeMacro(vtkPyramid, vtkCell3D);
48  void PrintSelf(ostream& os, vtkIndent indent) override;
49 
51 
54  void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
55  // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
56  VTK_LEGACY(virtual void GetEdgePoints(int edgeId, int*& pts) override);
57  vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
58  // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
59  VTK_LEGACY(virtual void GetFacePoints(int faceId, int*& pts) override);
60  void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
61  vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
62  vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
63  vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
64  vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
65  bool GetCentroid(double centroid[3]) const override;
66  bool IsInsideOut() override;
68 
72  static constexpr vtkIdType NumberOfPoints = 5;
73 
77  static constexpr vtkIdType NumberOfEdges = 8;
78 
82  static constexpr vtkIdType NumberOfFaces = 5;
83 
88  static constexpr vtkIdType MaximumFaceSize = 4;
89 
95  static constexpr vtkIdType MaximumValence = 4;
97 
100  int GetCellType() override { return VTK_PYRAMID; }
101  int GetCellDimension() override { return 3; }
102  int GetNumberOfEdges() override { return 8; }
103  int GetNumberOfFaces() override { return 5; }
104  vtkCell* GetEdge(int edgeId) override;
105  vtkCell* GetFace(int faceId) override;
106  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
107  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
108  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
109  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
110  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
111  double& dist2, double weights[]) override;
112  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
113  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
114  double pcoords[3], int& subId) override;
115  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
117  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
118  double* GetParametricCoords() override;
120 
128  static int* GetTriangleCases(int caseId);
129 
133  int GetParametricCenter(double pcoords[3]) override;
134 
135  static void InterpolationFunctions(const double pcoords[3], double weights[5]);
136  static void InterpolationDerivs(const double pcoords[3], double derivs[15]);
138 
142  void InterpolateFunctions(const double pcoords[3], double weights[5]) override
143  {
144  vtkPyramid::InterpolationFunctions(pcoords, weights);
145  }
146  void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
147  {
148  vtkPyramid::InterpolationDerivs(pcoords, derivs);
149  }
151 
152  int JacobianInverse(const double pcoords[3], double** inverse, double derivs[15]);
153 
155 
163  static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
164  static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(4);
166 
171 
176 
181 
186 
191 
195  static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
196 
197 protected:
199  ~vtkPyramid() override;
200 
204 
205 private:
206  vtkPyramid(const vtkPyramid&) = delete;
207  void operator=(const vtkPyramid&) = delete;
208 };
209 
210 //----------------------------------------------------------------------------
211 inline int vtkPyramid::GetParametricCenter(double pcoords[3])
212 {
213  pcoords[0] = pcoords[1] = 0.4;
214  pcoords[2] = 0.2;
215  return 0;
216 }
217 
218 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
object to represent cell connectivity
Definition: vtkCellArray.h:180
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:57
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
list of point or cell ids
Definition: vtkIdList.h:31
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:34
cell represents a 1D line
Definition: vtkLine.h:30
represent and manipulate point attribute data
Definition: vtkPointData.h:32
represent and manipulate 3D points
Definition: vtkPoints.h:34
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:44
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.
Definition: vtkPyramid.h:211
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
vtkCell * GetEdge(int edgeId) override
Return the edge cell from the edgeId of the cell.
static vtkPyramid * New()
vtkQuad * Quad
Definition: vtkPyramid.h:203
vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType *&pts) override
Get the ids of a one-ring surrounding point of id pointId.
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.
void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType *&pts) override
Get the ids of the two adjacent faces to edge of id edgeId.
vtkTriangle * Triangle
Definition: vtkPyramid.h:202
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
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
bool IsInsideOut() override
Returns true if the normals of the vtkCell3D point inside the cell.
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
Generate simplices of proper dimension.
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) override
Get the list of vertices that define a face.
static int * GetTriangleCases(int caseId)
Return the case table for table-based isocontouring (aka marching cubes style implementations).
~vtkPyramid() override
virtual void GetFacePoints(int faceId, int *&pts) override
virtual void GetEdgePoints(int edgeId, int *&pts) override
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
void InterpolateDerivs(const double pcoords[3], double derivs[15]) override
Definition: vtkPyramid.h:146
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkPyramid.h:102
vtkCell * GetFace(int faceId) override
Return the face cell from the faceId of the cell.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkLine * Line
Definition: vtkPyramid.h:201
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)
void GetEdgePoints(vtkIdType edgeId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
int GetCellDimension() override
The topological dimension of the cell.
Definition: vtkPyramid.h:101
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
Get the ids of the incident edges to point of id pointId.
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...
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkPyramid.h:103
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkPyramid.h:100
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) override
Get the ids of the adjacent faces to face of id faceId.
bool GetCentroid(double centroid[3]) const override
Computes the centroid of the cell.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Intersect with a ray.
static void InterpolationFunctions(const double pcoords[3], double weights[5])
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
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
Get the ids of the incident faces point of id pointId.
void InterpolateFunctions(const double pcoords[3], double weights[5]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkPyramid.h:142
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:36
a cell that represents a triangle
Definition: vtkTriangle.h:36
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:60
int vtkIdType
Definition: vtkType.h:338
#define VTK_SIZEHINT(...)