VTK  9.3.20240417
vtkQuadraticWedge.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
46 #ifndef vtkQuadraticWedge_h
47 #define vtkQuadraticWedge_h
48 
49 #include "vtkCommonDataModelModule.h" // For export macro
50 #include "vtkNonLinearCell.h"
51 
52 VTK_ABI_NAMESPACE_BEGIN
53 class vtkQuadraticEdge;
54 class vtkQuadraticQuad;
56 class vtkWedge;
57 class vtkDoubleArray;
58 
59 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticWedge : public vtkNonLinearCell
60 {
61 public:
64  void PrintSelf(ostream& os, vtkIndent indent) override;
65 
67 
71  int GetCellType() override { return VTK_QUADRATIC_WEDGE; }
72  int GetCellDimension() override { return 3; }
73  int GetNumberOfEdges() override { return 9; }
74  int GetNumberOfFaces() override { return 5; }
75  vtkCell* GetEdge(int edgeId) override;
76  vtkCell* GetFace(int faceId) override;
78 
79  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
80  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
81  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
82  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
83  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
84  double& dist2, double weights[]) override;
85  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
86  int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
88  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
89  double* GetParametricCoords() override;
90 
96  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
97  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
98  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
99 
104  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
105  double pcoords[3], int& subId) override;
106 
110  int GetParametricCenter(double pcoords[3]) override;
111 
112  static void InterpolationFunctions(const double pcoords[3], double weights[15]);
113  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
115 
119  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
120  {
122  }
123  void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
124  {
126  }
129 
136  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
137  static const vtkIdType* GetFaceArray(vtkIdType faceId);
139 
145  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
146 
147 protected:
149  ~vtkQuadraticWedge() override;
150 
158  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
159 
160  void Subdivide(
161  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
162 
163 private:
164  vtkQuadraticWedge(const vtkQuadraticWedge&) = delete;
165  void operator=(const vtkQuadraticWedge&) = delete;
166 };
167 //----------------------------------------------------------------------------
168 // Return the center of the quadratic wedge in parametric coordinates.
169 inline int vtkQuadraticWedge::GetParametricCenter(double pcoords[3])
170 {
171  pcoords[0] = pcoords[1] = 1. / 3;
172  pcoords[2] = 0.5;
173  return 0;
174 }
175 
176 VTK_ABI_NAMESPACE_END
177 #endif
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
dynamic, self-adjusting array of double
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
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:140
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 8-node isoparametric quad
cell represents a parabolic, isoparametric triangle
cell represents a parabolic, 15-node isoparametric wedge
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
vtkDoubleArray * Scalars
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
vtkQuadraticTriangle * TriangleFace
static vtkQuadraticWedge * New()
~vtkQuadraticWedge() override
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[45])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
vtkQuadraticEdge * Edge
vtkPointData * PointData
int GetNumberOfEdges() override
Implement the vtkCell API.
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
int GetCellType() override
Implement the vtkCell API.
vtkDoubleArray * CellScalars
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Line-edge intersection.
vtkCellData * CellData
vtkQuadraticQuad * Face
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
static void InterpolationFunctions(const double pcoords[3], double weights[15])
void InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetCellDimension() override
Implement the vtkCell API.
int GetNumberOfFaces() override
Implement the vtkCell API.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic wedge in parametric coordinates.
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tetras, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this quadratic hexahedron using scalar value provided.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
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
Generate contouring primitives.
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:85
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_QUADRATIC_WEDGE
Definition: vtkCellType.h:81
int vtkIdType
Definition: vtkType.h:315