VTK  9.1.0
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=========================================================================*/
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 // @deprecated Replaced by GetEdgePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
56 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkWedge::GetEdgePoints(vtkIdType, const vtkIdType*&)")
57 void GetEdgePoints(int edgeId, int*& pts) override;
58 vtkIdType GetFacePoints(vtkIdType faceId, const vtkIdType*& pts) override;
59 // @deprecated Replaced by GetFacePoints(vtkIdType, const vtkIdType*&) as of VTK 9.0
60 VTK_DEPRECATED_IN_9_0_0("Replaced by vtkWedge::GetFacePoints(vtkIdType, const vtkIdType*&)")
61 void GetFacePoints(int faceId, int*& pts) override;
62 void GetEdgeToAdjacentFaces(vtkIdType edgeId, const vtkIdType*& pts) override;
63 vtkIdType GetFaceToAdjacentFaces(vtkIdType faceId, const vtkIdType*& faceIds) override;
64 vtkIdType GetPointToIncidentEdges(vtkIdType pointId, const vtkIdType*& edgeIds) override;
65 vtkIdType GetPointToIncidentFaces(vtkIdType pointId, const vtkIdType*& faceIds) override;
66 vtkIdType GetPointToOneRingPoints(vtkIdType pointId, const vtkIdType*& pts) override;
67 bool GetCentroid(double centroid[3]) const override;
68 bool IsInsideOut() override;
70
74 static constexpr vtkIdType NumberOfPoints = 6;
75
79 static constexpr vtkIdType NumberOfEdges = 9;
80
84 static constexpr vtkIdType NumberOfFaces = 5;
85
90 static constexpr vtkIdType MaximumFaceSize = 4;
91
97 static constexpr vtkIdType MaximumValence = 3;
98
100
103 int GetCellType() override { return VTK_WEDGE; }
104 int GetCellDimension() override { return 3; }
105 int GetNumberOfEdges() override { return static_cast<int>(NumberOfEdges); }
106 int GetNumberOfFaces() override { return static_cast<int>(NumberOfFaces); }
107 vtkCell* GetEdge(int edgeId) override;
108 vtkCell* GetFace(int faceId) override;
109 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
110 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
111 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
112 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
113 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
114 double& dist2, double weights[]) override;
115 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
116 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
117 double pcoords[3], int& subId) override;
118 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
120 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
121 double* GetParametricCoords() override;
123
131 static int* GetTriangleCases(int caseId);
132
136 int GetParametricCenter(double pcoords[3]) override;
137
138 static void InterpolationFunctions(const double pcoords[3], double weights[6]);
139 static void InterpolationDerivs(const double pcoords[3], double derivs[18]);
141
145 void InterpolateFunctions(const double pcoords[3], double weights[6]) override
146 {
147 vtkWedge::InterpolationFunctions(pcoords, weights);
148 }
149 void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
150 {
151 vtkWedge::InterpolationDerivs(pcoords, derivs);
152 }
153 int JacobianInverse(const double pcoords[3], double** inverse, double derivs[18]);
155
157
166 static const vtkIdType* GetFaceArray(vtkIdType faceId) VTK_SIZEHINT(MaximumFaceSize + 1);
168
173
178
183
188
193
197 static bool ComputeCentroid(vtkPoints* points, const vtkIdType* pointIds, double centroid[3]);
198
199protected:
201 ~vtkWedge() override;
202
206
207private:
208 vtkWedge(const vtkWedge&) = delete;
209 void operator=(const vtkWedge&) = delete;
210};
211
212inline int vtkWedge::GetParametricCenter(double pcoords[3])
213{
214 pcoords[0] = pcoords[1] = 0.333333;
215 pcoords[2] = 0.5;
216 return 0;
217}
218
219#endif
abstract class to specify 3D cell interface
Definition: vtkCell3D.h:40
object to represent cell connectivity
Definition: vtkCellArray.h:181
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:58
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
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
Definition: vtkPointData.h:33
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
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:44
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()
int GetNumberOfEdges() override
See the vtkCell API for descriptions of these methods.
Definition: vtkWedge.h:105
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:106
~vtkWedge() override
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:203
vtkCell * GetEdge(int edgeId) override
See the vtkCell API for descriptions of these methods.
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:145
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
See the vtkCell API for descriptions of these methods.
int GetCellDimension() override
See the vtkCell API for descriptions of these methods.
Definition: vtkWedge.h:104
vtkTriangle * Triangle
Definition: vtkWedge.h:204
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 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’).
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:212
vtkQuad * Quad
Definition: vtkWedge.h:205
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.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
See the vtkCell API for descriptions of these methods.
void InterpolateDerivs(const double pcoords[3], double derivs[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkWedge.h:149
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
int GetCellType(const Ioss::ElementTopology *topology)
Returns VTK celltype for a Ioss topology element.
@ points
Definition: vtkX3D.h:452
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_WEDGE
Definition: vtkCellType.h:59
#define VTK_DEPRECATED_IN_9_0_0(reason)
int vtkIdType
Definition: vtkType.h:332
#define VTK_SIZEHINT(...)