VTK  9.3.20240419
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
71 #ifndef vtkWedge_h
72 #define vtkWedge_h
73 
74 #include "vtkCell3D.h"
75 #include "vtkCommonDataModelModule.h" // For export macro
76 
77 VTK_ABI_NAMESPACE_BEGIN
78 class vtkLine;
79 class vtkTriangle;
80 class vtkQuad;
83 
84 class VTKCOMMONDATAMODEL_EXPORT vtkWedge : public vtkCell3D
85 {
86 public:
87  static vtkWedge* New();
88  vtkTypeMacro(vtkWedge, vtkCell3D);
89  void PrintSelf(ostream& os, vtkIndent indent) override;
90 
92 
95  void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
96  vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
97  void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
98  vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
99  vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
100  vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
101  vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
102  bool GetCentroid(double centroid[3]) const override;
103  bool IsInsideOut() override;
105 
109  static constexpr vtkIdType NumberOfPoints = 6;
110 
114  static constexpr vtkIdType NumberOfEdges = 9;
115 
119  static constexpr vtkIdType NumberOfFaces = 5;
120 
125  static constexpr vtkIdType MaximumFaceSize = 4;
126 
132  static constexpr vtkIdType MaximumValence = 3;
133 
135 
138  int GetCellType() override { return VTK_WEDGE; }
139  int GetCellDimension() override { return 3; }
140  int GetNumberOfEdges() override { return static_cast<int>(NumberOfEdges); }
141  int GetNumberOfFaces() override { return static_cast<int>(NumberOfFaces); }
142  vtkCell* GetEdge(int edgeId) override;
143  vtkCell* GetFace(int faceId) override;
144  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
145  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
146  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
147  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
148  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
149  double& dist2, double weights[]) override;
150  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
151  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
152  double pcoords[3], int& subId) override;
153  int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
155  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
156  double* GetParametricCoords() override;
158 
166  static int* GetTriangleCases(int caseId);
167 
171  int GetParametricCenter(double pcoords[3]) override;
172 
173  static void InterpolationFunctions(const double pcoords[3], double weights[6]);
174  static void InterpolationDerivs(const double pcoords[3], double derivs[18]);
176 
180  void InterpolateFunctions(const double pcoords[3], double weights[6]) override
181  {
182  vtkWedge::InterpolationFunctions(pcoords, weights);
183  }
184  void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
185  {
186  vtkWedge::InterpolationDerivs(pcoords, derivs);
187  }
188  int JacobianInverse(const double pcoords[3], double** inverse, double derivs[18]);
190 
192 
200  static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
201  static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(MaximumFaceSize + 1);
203 
208 
213 
218 
223 
228 
232  static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
233 
234 protected:
236  ~vtkWedge() override;
237 
241 
242 private:
243  vtkWedge(const vtkWedge&) = delete;
244  void operator=(const vtkWedge&) = delete;
245 };
246 
247 inline int vtkWedge::GetParametricCenter(double pcoords[3])
248 {
249  pcoords[0] = pcoords[1] = 0.333333;
250  pcoords[2] = 0.5;
251  return 0;
252 }
253 
254 VTK_ABI_NAMESPACE_END
255 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:28
object to represent cell connectivity
Definition: vtkCellArray.h:286
represent and manipulate cell attribute data
Definition: vtkCellData.h:141
abstract class to specify cell behavior
Definition: vtkCell.h:130
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:155
list of point or cell ids
Definition: vtkIdList.h:133
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:108
cell represents a 1D line
Definition: vtkLine.h:132
represent and manipulate point attribute data
Definition: vtkPointData.h:140
represent and manipulate 3D points
Definition: vtkPoints.h:139
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:87
a cell that represents a triangle
Definition: vtkTriangle.h:137
dataset represents arbitrary combinations of all possible cell types
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:85
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:140
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:141
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:238
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:180
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:139
vtkTriangle * Triangle
Definition: vtkWedge.h:239
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.
int TriangulateLocalIds(int index, vtkIdList *ptIds) 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:247
vtkQuad * Quad
Definition: vtkWedge.h:240
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 GetCellType() override
See the vtkCell API for descriptions of these methods.
Definition: vtkWedge.h:138
void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkWedge.h:184
const std::string pointIds
@ points
Definition: vtkX3D.h:446
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_WEDGE
Definition: vtkCellType.h:69
int vtkIdType
Definition: vtkType.h:315
#define VTK_SIZEHINT(...)