VTK  9.3.0
vtkQuadraticLinearWedge.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
30 #ifndef vtkQuadraticLinearWedge_h
31 #define vtkQuadraticLinearWedge_h
32 
33 #include "vtkCommonDataModelModule.h" // For export macro
34 #include "vtkNonLinearCell.h"
35 
36 VTK_ABI_NAMESPACE_BEGIN
37 class vtkQuadraticEdge;
38 class vtkLine;
41 class vtkWedge;
42 class vtkDoubleArray;
43 
44 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearWedge : public vtkNonLinearCell
45 {
46 public:
49  void PrintSelf(ostream& os, vtkIndent indent) override;
50 
52 
56  int GetCellType() override { return VTK_QUADRATIC_LINEAR_WEDGE; }
57  int GetCellDimension() override { return 3; }
58  int GetNumberOfEdges() override { return 9; }
59  int GetNumberOfFaces() override { return 5; }
60  vtkCell* GetEdge(int edgeId) override;
61  vtkCell* GetFace(int faceId) override;
63 
64  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
65 
67 
71  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
72  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
73  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
74  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
75  double& dist2, double* weights) override;
76  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
77  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
79  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
80  double* GetParametricCoords() override;
82 
88  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
89  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
90  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
91 
96  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
97  double pcoords[3], int& subId) override;
98 
102  int GetParametricCenter(double pcoords[3]) override;
103 
104  static void InterpolationFunctions(const double pcoords[3], double weights[12]);
105  static void InterpolationDerivs(const double pcoords[3], double derivs[36]);
107 
111  void InterpolateFunctions(const double pcoords[3], double weights[12]) override
112  {
114  }
115  void InterpolateDerivs(const double pcoords[3], double derivs[36]) override
116  {
118  }
121 
128  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
129  static const vtkIdType* GetFaceArray(vtkIdType faceId);
131 
137  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[36]);
138 
139 protected:
142 
148  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
149 
150 private:
152  void operator=(const vtkQuadraticLinearWedge&) = delete;
153 };
154 //----------------------------------------------------------------------------
155 // Return the center of the quadratic wedge in parametric coordinates.
156 inline int vtkQuadraticLinearWedge::GetParametricCenter(double pcoords[3])
157 {
158  pcoords[0] = pcoords[1] = 1. / 3;
159  pcoords[2] = 0.5;
160  return 0;
161 }
162 
163 VTK_ABI_NAMESPACE_END
164 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:176
represent and manipulate cell attribute data
Definition: vtkCellData.h:31
abstract class to specify cell behavior
Definition: vtkCell.h:50
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:45
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:23
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:29
cell represents a 1D line
Definition: vtkLine.h:23
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:30
represent and manipulate 3D points
Definition: vtkPoints.h:29
cell represents a parabolic, isoparametric edge
cell represents a quadratic-linear, 6-node isoparametric quad
cell represents a, 12-node isoparametric wedge
vtkQuadraticTriangle * TriangleFace
static vtkQuadraticLinearWedge * New()
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic linear wedge in parametric coordinates.
static void InterpolationDerivs(const double pcoords[3], double derivs[36])
~vtkQuadraticLinearWedge() override
int EvaluatePosition(const double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights) override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
double * GetParametricCoords() override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
void InterpolateDerivs(const double pcoords[3], double derivs[36]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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 Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
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
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
static void InterpolationFunctions(const double pcoords[3], double weights[12])
void InterpolateFunctions(const double pcoords[3], double weights[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetNumberOfFaces() override
Implement the vtkCell API.
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.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkQuadraticLinearQuad * Face
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 linear wedge using scalar value provided.
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
int GetCellType() override
Implement the vtkCell API.
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[36])
Given parametric coordinates compute inverse Jacobian transformation matrix.
int GetCellDimension() override
Implement the vtkCell API.
int GetNumberOfEdges() override
Implement the vtkCell API.
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:36
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:68
int vtkIdType
Definition: vtkType.h:315