VTK  9.3.20240420
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
49#ifndef vtkQuadraticLinearWedge_h
50#define vtkQuadraticLinearWedge_h
51
52#include "vtkCommonDataModelModule.h" // For export macro
53#include "vtkNonLinearCell.h"
54
55VTK_ABI_NAMESPACE_BEGIN
57class vtkLine;
60class vtkWedge;
61class vtkDoubleArray;
62
63class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearWedge : public vtkNonLinearCell
64{
65public:
68 void PrintSelf(ostream& os, vtkIndent indent) override;
69
71
75 int GetCellType() override { return VTK_QUADRATIC_LINEAR_WEDGE; }
76 int GetCellDimension() override { return 3; }
77 int GetNumberOfEdges() override { return 9; }
78 int GetNumberOfFaces() override { return 5; }
79 vtkCell* GetEdge(int edgeId) override;
80 vtkCell* GetFace(int faceId) override;
82
83 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
84
86
90 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
91 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
92 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
93 int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
94 double& dist2, double* weights) override;
95 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
96 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
98 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
99 double* GetParametricCoords() override;
101
107 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
108 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
109 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
110
115 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
116 double pcoords[3], int& subId) override;
117
121 int GetParametricCenter(double pcoords[3]) override;
122
123 static void InterpolationFunctions(const double pcoords[3], double weights[12]);
124 static void InterpolationDerivs(const double pcoords[3], double derivs[36]);
126
130 void InterpolateFunctions(const double pcoords[3], double weights[12]) override
131 {
133 }
134 void InterpolateDerivs(const double pcoords[3], double derivs[36]) override
135 {
137 }
140
147 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
148 static const vtkIdType* GetFaceArray(vtkIdType faceId);
150
156 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[36]);
157
158protected:
161
167 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
168
169private:
171 void operator=(const vtkQuadraticLinearWedge&) = delete;
172};
173//----------------------------------------------------------------------------
174// Return the center of the quadratic wedge in parametric coordinates.
176{
177 pcoords[0] = pcoords[1] = 1. / 3;
178 pcoords[2] = 0.5;
179 return 0;
180}
181
182VTK_ABI_NAMESPACE_END
183#endif
object to represent cell connectivity
represent and manipulate cell attribute data
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
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
cell represents a 1D line
Definition vtkLine.h:132
abstract superclass for non-linear cells
represent and manipulate point attribute data
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.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
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:85
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition vtkCellType.h:87
int vtkIdType
Definition vtkType.h:315