VTK  9.2.6
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=========================================================================*/
30
31#ifndef vtkWedge_h
32#define vtkWedge_h
33
34#include "vtkCell3D.h"
35#include "vtkCommonDataModelModule.h" // For export macro
36
37class vtkLine;
38class vtkTriangle;
39class vtkQuad;
42
43class VTKCOMMONDATAMODEL_EXPORT vtkWedge : public vtkCell3D
44{
45public:
46 static vtkWedge* New();
47 vtkTypeMacro(vtkWedge, vtkCell3D);
48 void PrintSelf(ostream& os, vtkIndent indent) override;
49
51
54 void GetEdgePoints(vtkIdType edgeId, const vtkIdType*& pts) override;
55 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
56 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
57 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
58 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
59 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
60 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
61 bool GetCentroid(double centroid[3]) const override;
62 bool IsInsideOut() override;
64
68 static constexpr vtkIdType NumberOfPoints = 6;
69
73 static constexpr vtkIdType NumberOfEdges = 9;
74
78 static constexpr vtkIdType NumberOfFaces = 5;
79
84 static constexpr vtkIdType MaximumFaceSize = 4;
85
91 static constexpr vtkIdType MaximumValence = 3;
92
94
97 int GetCellType() override { return VTK_WEDGE; }
98 int GetCellDimension() override { return 3; }
99 int GetNumberOfEdges() override { return static_cast<int>(NumberOfEdges); }
100 int GetNumberOfFaces() override { return static_cast<int>(NumberOfFaces); }
101 vtkCell* GetEdge(int edgeId) override;
102 vtkCell* GetFace(int faceId) override;
103 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
104 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
105 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
106 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
107 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
108 double& dist2, double weights[]) override;
109 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
110 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
111 double pcoords[3], int& subId) override;
112 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
114 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
115 double* GetParametricCoords() override;
117
125 static int* GetTriangleCases(int caseId);
126
130 int GetParametricCenter(double pcoords[3]) override;
131
132 static void InterpolationFunctions(const double pcoords[3], double weights[6]);
133 static void InterpolationDerivs(const double pcoords[3], double derivs[18]);
135
139 void InterpolateFunctions(const double pcoords[3], double weights[6]) override
140 {
141 vtkWedge::InterpolationFunctions(pcoords, weights);
142 }
143 void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
144 {
145 vtkWedge::InterpolationDerivs(pcoords, derivs);
146 }
147 int JacobianInverse(const double pcoords[3], double** inverse, double derivs[18]);
149
151
162
167
172
177
182
187
191 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
192
193protected:
195 ~vtkWedge() override;
196
200
201private:
202 vtkWedge(const vtkWedge&) = delete;
203 void operator=(const vtkWedge&) = delete;
204};
205
206inline int vtkWedge::GetParametricCenter(double pcoords[3])
207{
208 pcoords[0] = pcoords[1] = 0.333333;
209 pcoords[2] = 0.5;
210 return 0;
211}
212
213#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:36
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
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:31
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:34
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
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.
static vtkWedge * New()
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:99
static const vtkIdType * GetEdgeToAdjacentFacesArray(vtkIdType edgeId)
Static method version of GetEdgeToAdjacentFaces.
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:100
vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType *&pts) override
See vtkCell3D API for description of these methods.
~vtkWedge() override
static constexpr vtkIdType MaximumFaceSize
static contexpr handle on the maximum face size.
Definition vtkWedge.h:84
static void InterpolationDerivs(const double pcoords[3], double derivs[18])
double * GetParametricCoords() override
See the vtkCell API for descriptions of these methods.
static const vtkIdType * GetFaceToAdjacentFacesArray(vtkIdType faceId)
Static method version of GetFaceToAdjacentFaces.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static const vtkIdType * GetPointToIncidentFacesArray(vtkIdType pointId)
Static method version of GetPointToIncidentFacesArray.
vtkLine * Line
Definition vtkWedge.h:197
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
static constexpr vtkIdType MaximumValence
static constexpr handle on the maximum valence of this cell.
Definition vtkWedge.h:91
static int * GetTriangleCases(int caseId)
Return the case table for table-based isocontouring (aka marching cubes style implementations).
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:139
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.
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition vtkWedge.h:98
vtkTriangle * Triangle
Definition vtkWedge.h:198
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.
static constexpr vtkIdType NumberOfEdges
static contexpr handle on the number of faces.
Definition vtkWedge.h:73
vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType *&faceIds) override
See vtkCell3D API for description of these methods.
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 * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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:206
vtkQuad * Quad
Definition vtkWedge.h:199
static constexpr vtkIdType NumberOfPoints
static constexpr handle on the number of points.
Definition vtkWedge.h:68
static constexpr vtkIdType NumberOfFaces
static contexpr handle on the number of edges.
Definition vtkWedge.h:78
int JacobianInverse(const double pcoords[3], double **inverse, double derivs[18])
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
vtkCell * GetFace(int faceId) override
See the vtkCell API for descriptions of these methods.
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:97
void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition vtkWedge.h:143
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
@ VTK_WEDGE
Definition vtkCellType.h:59
#define vtkDataArray
int vtkIdType
Definition vtkType.h:332
#define VTK_SIZEHINT(...)