VTK  9.1.20211115
vtkWedge.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkWedge.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 =========================================================================*/
80 #ifndef vtkWedge_h
81 #define vtkWedge_h
82 
83 #include "vtkCell3D.h"
84 #include "vtkCommonDataModelModule.h" // For export macro
85 
86 class vtkLine;
87 class vtkTriangle;
88 class vtkQuad;
91 
92 class VTKCOMMONDATAMODEL_EXPORT vtkWedge : public vtkCell3D
93 {
94 public:
95  static vtkWedge* New();
96  vtkTypeMacro(vtkWedge, vtkCell3D);
97  void PrintSelf(ostream& os, vtkIndent indent) override;
98 
100 
103  void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
104  vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
105  void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
106  vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
107  vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
108  vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
109  vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
110  bool GetCentroid(double centroid[3]) const override;
111  bool IsInsideOut() override;
113 
117  static constexpr vtkIdType NumberOfPoints = 6;
118 
122  static constexpr vtkIdType NumberOfEdges = 9;
123 
127  static constexpr vtkIdType NumberOfFaces = 5;
128 
133  static constexpr vtkIdType MaximumFaceSize = 4;
134 
140  static constexpr vtkIdType MaximumValence = 3;
141 
143 
146  int GetCellType() override { return VTK_WEDGE; }
147  int GetCellDimension() override { return 3; }
148  int GetNumberOfEdges() override { return static_cast<int>(NumberOfEdges); }
149  int GetNumberOfFaces() override { return static_cast<int>(NumberOfFaces); }
150  vtkCell* GetEdge(int edgeId) override;
151  vtkCell* GetFace(int faceId) override;
152  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
153  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
154  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
155  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
156  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
157  double& dist2, double weights[]) override;
158  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
159  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
160  double pcoords[3], int& subId) override;
161  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
163  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
164  double* GetParametricCoords() override;
166 
174  static int* GetTriangleCases(int caseId);
175 
179  int GetParametricCenter(double pcoords[3]) override;
180 
181  static void InterpolationFunctions(const double pcoords[3], double weights[6]);
182  static void InterpolationDerivs(const double pcoords[3], double derivs[18]);
184 
188  void InterpolateFunctions(const double pcoords[3], double weights[6]) override
189  {
190  vtkWedge::InterpolationFunctions(pcoords, weights);
191  }
192  void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
193  {
194  vtkWedge::InterpolationDerivs(pcoords, derivs);
195  }
196  int JacobianInverse(const double pcoords[3], double** inverse, double derivs[18]);
198 
200 
208  static const vtkIdType* GetEdgeArray(vtkIdType edgeId) VTK_SIZEHINT(2);
209  static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(MaximumFaceSize + 1);
211 
216 
221 
226 
231 
236 
240  static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
241 
242 protected:
244  ~vtkWedge() override;
245 
249 
250 private:
251  vtkWedge(const vtkWedge&) = delete;
252  void operator=(const vtkWedge&) = delete;
253 };
254 
255 inline int vtkWedge::GetParametricCenter(double pcoords[3])
256 {
257  pcoords[0] = pcoords[1] = 0.333333;
258  pcoords[2] = 0.5;
259  return 0;
260 }
261 
262 #endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:39
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
list of point or cell ids
Definition: vtkIdList.h:140
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:113
cell represents a 1D line
Definition: vtkLine.h:140
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:95
a cell that represents a triangle
Definition: vtkTriangle.h:145
dataset represents arbitrary combinations of all possible cell types
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:93
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:148
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:149
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:246
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:188
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:147
vtkTriangle * Triangle
Definition: vtkWedge.h:247
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:255
vtkQuad * Quad
Definition: vtkWedge.h:248
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:146
void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkWedge.h:192
@ points
Definition: vtkX3D.h:452
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_WEDGE
Definition: vtkCellType.h:98
int vtkIdType
Definition: vtkType.h:332
#define VTK_SIZEHINT(...)