VTK  9.3.20240424
vtkBiQuadraticQuadraticWedge.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
51#ifndef vtkBiQuadraticQuadraticWedge_h
52#define vtkBiQuadraticQuadraticWedge_h
53
54#include "vtkCommonDataModelModule.h" // For export macro
55#include "vtkNonLinearCell.h"
56
57VTK_ABI_NAMESPACE_BEGIN
61class vtkWedge;
62class vtkDoubleArray;
63
64class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticWedge : public vtkNonLinearCell
65{
66public:
69 void PrintSelf(ostream& os, vtkIndent indent) override;
70
72
77 int GetCellDimension() override { return 3; }
78 int GetNumberOfEdges() override { return 9; }
79 int GetNumberOfFaces() override { return 5; }
80 vtkCell* GetEdge(int edgeId) override;
81 vtkCell* GetFace(int faceId) override;
83
84 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
85 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
86 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
87 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
88 int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
89 double& dist2, double* weights) override;
90 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
91 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
93 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
94 double* GetParametricCoords() override;
95
101 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
102 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
103 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
104
109 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
110 double pcoords[3], int& subId) override;
111
115 int GetParametricCenter(double pcoords[3]) override;
116
117 static void InterpolationFunctions(const double pcoords[3], double weights[18]);
118 static void InterpolationDerivs(const double pcoords[3], double derivs[54]);
120
124 void InterpolateFunctions(const double pcoords[3], double weights[18]) override
125 {
127 }
128 void InterpolateDerivs(const double pcoords[3], double derivs[54]) override
129 {
131 }
134
141 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
142 static const vtkIdType* GetFaceArray(vtkIdType faceId);
144
150 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[54]);
151
152protected:
155
160 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
161
162private:
164 void operator=(const vtkBiQuadraticQuadraticWedge&) = delete;
165};
166//----------------------------------------------------------------------------
167// Return the center of the quadratic wedge in parametric coordinates.
169{
170 pcoords[0] = pcoords[1] = 1. / 3;
171 pcoords[2] = 0.5;
172 return 0;
173}
174
175VTK_ABI_NAMESPACE_END
176#endif
cell represents a parabolic, 9-node isoparametric quad
cell represents a parabolic, 18-node isoparametric wedge
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
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 Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
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.
static void InterpolationFunctions(const double pcoords[3], double weights[18])
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void InterpolateDerivs(const double pcoords[3], double derivs[54]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetNumberOfEdges() override
Implement the vtkCell API.
int GetCellDimension() override
Implement the vtkCell API.
static vtkBiQuadraticQuadraticWedge * New()
int GetNumberOfFaces() override
Implement the vtkCell API.
static void InterpolationDerivs(const double pcoords[3], double derivs[54])
int GetCellType() override
Implement the vtkCell API.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
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.
void InterpolateFunctions(const double pcoords[3], double weights[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int EvaluatePosition(const double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights) override
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
~vtkBiQuadraticQuadraticWedge() override
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[54])
Given parametric coordinates compute inverse Jacobian transformation matrix.
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
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 Wedge using scalar value provided.
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
abstract superclass for non-linear cells
represent and manipulate point attribute data
cell represents a parabolic, isoparametric edge
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a linear wedge
Definition vtkWedge.h:85
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
Definition vtkCellType.h:88
int vtkIdType
Definition vtkType.h:315