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
31 #ifndef vtkWedge_h
32 #define vtkWedge_h
33 
34 #include "vtkCell3D.h"
35 #include "vtkCommonDataModelModule.h" // For export macro
36 
37 VTK_ABI_NAMESPACE_BEGIN
38 class vtkLine;
39 class vtkTriangle;
40 class vtkQuad;
43 
44 class VTKCOMMONDATAMODEL_EXPORT vtkWedge : public vtkCell3D
45 {
46 public:
47  static vtkWedge* New();
48  vtkTypeMacro(vtkWedge, vtkCell3D);
49  void PrintSelf(ostream& os, vtkIndent indent) override;
50 
52 
55  void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
56  vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
57  void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
58  vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
59  vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
60  vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
61  vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
62  bool GetCentroid(double centroid[3]) const override;
63  bool IsInsideOut() override;
65 
69  static constexpr vtkIdType NumberOfPoints = 6;
70 
74  static constexpr vtkIdType NumberOfEdges = 9;
75 
79  static constexpr vtkIdType NumberOfFaces = 5;
80 
85  static constexpr vtkIdType MaximumFaceSize = 4;
86 
92  static constexpr vtkIdType MaximumValence = 3;
93 
95 
98  int GetCellType() override { return VTK_WEDGE; }
99  int GetCellDimension() override { return 3; }
100  int GetNumberOfEdges() override { return static_cast<int>(NumberOfEdges); }
101  int GetNumberOfFaces() override { return static_cast<int>(NumberOfFaces); }
102  vtkCell* GetEdge(int edgeId) override;
103  vtkCell* GetFace(int faceId) override;
104  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
105  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
106  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
107  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
108  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
109  double& dist2, double weights[]) override;
110  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
111  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
112  double pcoords[3], int& subId) override;
113  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
115  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
116  double* GetParametricCoords() override;
118 
126  static int* GetTriangleCases(int caseId);
127 
131  int GetParametricCenter(double pcoords[3]) override;
132 
133  static void InterpolationFunctions(const double pcoords[3], double weights[6]);
134  static void InterpolationDerivs(const double pcoords[3], double derivs[18]);
136 
140  void InterpolateFunctions(const double pcoords[3], double weights[6]) override
141  {
142  vtkWedge::InterpolationFunctions(pcoords, weights);
143  }
144  void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
145  {
146  vtkWedge::InterpolationDerivs(pcoords, derivs);
147  }
148  int JacobianInverse(const double pcoords[3], double** inverse, double derivs[18]);
150 
152 
160  static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
161  static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(MaximumFaceSize + 1);
163 
168 
173 
178 
183 
188 
192  static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
193 
194 protected:
196  ~vtkWedge() override;
197 
201 
202 private:
203  vtkWedge(const vtkWedge&) = delete;
204  void operator=(const vtkWedge&) = delete;
205 };
206 
207 inline int vtkWedge::GetParametricCenter(double pcoords[3])
208 {
209  pcoords[0] = pcoords[1] = 0.333333;
210  pcoords[2] = 0.5;
211  return 0;
212 }
213 
214 VTK_ABI_NAMESPACE_END
215 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:28
object to represent cell connectivity
Definition: vtkCellArray.h:185
represent and manipulate cell attribute data
Definition: vtkCellData.h:40
abstract class to specify cell behavior
Definition: vtkCell.h:59
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
list of point or cell ids
Definition: vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:38
cell represents a 1D line
Definition: vtkLine.h:32
represent and manipulate point attribute data
Definition: vtkPointData.h:39
represent and manipulate 3D points
Definition: vtkPoints.h:38
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:37
a cell that represents a triangle
Definition: vtkTriangle.h:37
dataset represents arbitrary combinations of all possible cell types
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:45
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:100
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:101
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:198
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:140
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:99
vtkTriangle * Triangle
Definition: vtkWedge.h:199
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:207
vtkQuad * Quad
Definition: vtkWedge.h:200
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:98
void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkWedge.h:144
@ points
Definition: vtkX3D.h:446
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_WEDGE
Definition: vtkCellType.h:59
int vtkIdType
Definition: vtkType.h:315
#define VTK_SIZEHINT(...)