VTK  9.3.0
vtkWedge.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
22 #ifndef vtkWedge_h
23 #define vtkWedge_h
24 
25 #include "vtkCell3D.h"
26 #include "vtkCommonDataModelModule.h" // For export macro
27 
28 VTK_ABI_NAMESPACE_BEGIN
29 class vtkLine;
30 class vtkTriangle;
31 class vtkQuad;
34 
35 class VTKCOMMONDATAMODEL_EXPORT vtkWedge : public vtkCell3D
36 {
37 public:
38  static vtkWedge* New();
39  vtkTypeMacro(vtkWedge, 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;
56 
60  static constexpr vtkIdType NumberOfPoints = 6;
61 
65  static constexpr vtkIdType NumberOfEdges = 9;
66 
70  static constexpr vtkIdType NumberOfFaces = 5;
71 
76  static constexpr vtkIdType MaximumFaceSize = 4;
77 
83  static constexpr vtkIdType MaximumValence = 3;
84 
86 
89  int GetCellType() override { return VTK_WEDGE; }
90  int GetCellDimension() override { return 3; }
91  int GetNumberOfEdges() override { return static_cast<int>(NumberOfEdges); }
92  int GetNumberOfFaces() override { return static_cast<int>(NumberOfFaces); }
93  vtkCell* GetEdge(int edgeId) override;
94  vtkCell* GetFace(int faceId) override;
95  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
96  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
97  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
98  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
99  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
100  double& dist2, double weights[]) override;
101  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
102  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
103  double pcoords[3], int& subId) override;
104  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
106  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
107  double* GetParametricCoords() override;
109 
117  static int* GetTriangleCases(int caseId);
118 
122  int GetParametricCenter(double pcoords[3]) override;
123 
124  static void InterpolationFunctions(const double pcoords[3], double weights[6]);
125  static void InterpolationDerivs(const double pcoords[3], double derivs[18]);
127 
131  void InterpolateFunctions(const double pcoords[3], double weights[6]) override
132  {
133  vtkWedge::InterpolationFunctions(pcoords, weights);
134  }
135  void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
136  {
137  vtkWedge::InterpolationDerivs(pcoords, derivs);
138  }
139  int JacobianInverse(const double pcoords[3], double** inverse, double derivs[18]);
141 
143 
151  static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
152  static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(MaximumFaceSize + 1);
154 
159 
164 
169 
174 
179 
183  static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
184 
185 protected:
187  ~vtkWedge() override;
188 
192 
193 private:
194  vtkWedge(const vtkWedge&) = delete;
195  void operator=(const vtkWedge&) = delete;
196 };
197 
198 inline int vtkWedge::GetParametricCenter(double pcoords[3])
199 {
200  pcoords[0] = pcoords[1] = 0.333333;
201  pcoords[2] = 0.5;
202  return 0;
203 }
204 
205 VTK_ABI_NAMESPACE_END
206 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:28
object to represent cell connectivity
Definition: vtkCellArray.h:176
represent and manipulate cell attribute data
Definition: vtkCellData.h:31
abstract class to specify cell behavior
Definition: vtkCell.h:50
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
list of point or cell ids
Definition: vtkIdList.h:23
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
Definition: vtkPointData.h:30
represent and manipulate 3D points
Definition: vtkPoints.h:29
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
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:36
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.
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.
vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkWedge.h:91
static vtkWedge * New()
static bool ComputeCentroid(vtkPoints *points, const vtkIdType *pointIds, double centroid[3])
Static method version of GetCentroid.
int GetNumberOfFaces() override
See the vtkCell API for descriptions of these methods.
Definition: vtkWedge.h:92
vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
~vtkWedge() override
vtkCell * GetFace(int faceId) override
See the vtkCell API for descriptions of these methods.
double * GetParametricCoords() override
See the vtkCell API for descriptions of these methods.
static void InterpolationDerivs(const double pcoords[3], double derivs[18])
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static const vtkIdType * GetPointToOneRingPointsArray(vtkIdType pointId)
Static method version of GetPointToOneRingPoints.
static const vtkIdType * GetPointToIncidentEdgesArray(vtkIdType pointId)
Static method version of GetPointToIncidentEdgesArray.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
vtkLine * Line
Definition: vtkWedge.h:189
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
static void InterpolationFunctions(const double pcoords[3], double weights[6])
void InterpolateFunctions(const double pcoords[3], double weights[6]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkWedge.h:131
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition: vtkWedge.h:90
vtkTriangle * Triangle
Definition: vtkWedge.h:190
vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType *&edgeIds) override
See vtkCell3D API for description of these methods.
void GetEdgePoints(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.
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) 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).
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
bool GetCentroid(double centroid[3]) const 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 GetParametricCenter(double pcoords[3]) override
Return the center of the wedge in parametric coordinates.
Definition: vtkWedge.h:198
vtkQuad * Quad
Definition: vtkWedge.h:191
int JacobianInverse(const double pcoords[3], double **inverse, double derivs[18])
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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.
bool IsInsideOut() override
See vtkCell3D API for description of these methods.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
See the vtkCell API for descriptions of these methods.
int GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkWedge.h:89
void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkWedge.h:135
@ points
Definition: vtkX3D.h:446
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_WEDGE
Definition: vtkCellType.h:50
int vtkIdType
Definition: vtkType.h:315
#define VTK_SIZEHINT(...)