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
39#ifndef vtkQuadraticLinearWedge_h
40#define vtkQuadraticLinearWedge_h
41
42#include "vtkCommonDataModelModule.h" // For export macro
43#include "vtkNonLinearCell.h"
44
45VTK_ABI_NAMESPACE_BEGIN
47class vtkLine;
50class vtkWedge;
51class vtkDoubleArray;
52
53class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearWedge : public vtkNonLinearCell
54{
55public:
58 void PrintSelf(ostream& os, vtkIndent indent) override;
59
61
65 int GetCellType() override { return VTK_QUADRATIC_LINEAR_WEDGE; }
66 int GetCellDimension() override { return 3; }
67 int GetNumberOfEdges() override { return 9; }
68 int GetNumberOfFaces() override { return 5; }
69 vtkCell* GetEdge(int edgeId) override;
70 vtkCell* GetFace(int faceId) override;
72
73 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
74
76
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, 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 Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
88 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
89 double* GetParametricCoords() override;
91
97 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
98 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
99 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
100
105 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
106 double pcoords[3], int& subId) override;
107
111 int GetParametricCenter(double pcoords[3]) override;
112
113 static void InterpolationFunctions(const double pcoords[3], double weights[12]);
114 static void InterpolationDerivs(const double pcoords[3], double derivs[36]);
116
120 void InterpolateFunctions(const double pcoords[3], double weights[12]) override
121 {
123 }
124 void InterpolateDerivs(const double pcoords[3], double derivs[36]) override
125 {
127 }
130
137 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
138 static const vtkIdType* GetFaceArray(vtkIdType faceId);
140
146 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[36]);
147
148protected:
151
157 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
158
159private:
161 void operator=(const vtkQuadraticLinearWedge&) = delete;
162};
163//----------------------------------------------------------------------------
164// Return the center of the quadratic wedge in parametric coordinates.
166{
167 pcoords[0] = pcoords[1] = 1. / 3;
168 pcoords[2] = 0.5;
169 return 0;
170}
171
172VTK_ABI_NAMESPACE_END
173#endif
object to represent cell connectivity
represent and manipulate cell attribute data
Definition vtkCellData.h:40
abstract class to specify cell behavior
Definition vtkCell.h:59
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:32
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:38
cell represents a 1D line
Definition vtkLine.h:32
abstract superclass for non-linear cells
represent and manipulate point attribute data
represent and manipulate 3D points
Definition vtkPoints.h:38
cell represents a parabolic, isoparametric edge
cell represents a quadratic-linear, 6-node isoparametric quad
cell represents a, 12-node isoparametric wedge
vtkQuadraticTriangle * TriangleFace
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...
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
void InterpolateDerivs(const double pcoords[3], double derivs[36]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
double * GetParametricCoords() override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
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 * 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
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...
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.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static vtkQuadraticLinearWedge * New()
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:45
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition vtkCellType.h:68
int vtkIdType
Definition vtkType.h:315