VTK  9.1.0
vtkBiQuadraticQuadraticWedge.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkBiQuadraticQuadraticWedge.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
60#ifndef vtkBiQuadraticQuadraticWedge_h
61#define vtkBiQuadraticQuadraticWedge_h
62
63#include "vtkCommonDataModelModule.h" // For export macro
64#include "vtkNonLinearCell.h"
65
69class vtkWedge;
70class vtkDoubleArray;
71
72class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticWedge : public vtkNonLinearCell
73{
74public:
77 void PrintSelf(ostream& os, vtkIndent indent) override;
78
80
85 int GetCellDimension() override { return 3; }
86 int GetNumberOfEdges() override { return 9; }
87 int GetNumberOfFaces() override { return 5; }
88 vtkCell* GetEdge(int edgeId) override;
89 vtkCell* GetFace(int faceId) override;
91
92 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
93 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
94 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
95 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
96 int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
97 double& dist2, double* weights) override;
98 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
99 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
101 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
102 double* GetParametricCoords() override;
103
109 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
110 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
111 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
112
117 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
118 double pcoords[3], int& subId) override;
119
123 int GetParametricCenter(double pcoords[3]) override;
124
125 static void InterpolationFunctions(const double pcoords[3], double weights[15]);
126 static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
128
132 void InterpolateFunctions(const double pcoords[3], double weights[15]) override
133 {
135 }
136 void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
137 {
139 }
142
149 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
150 static const vtkIdType* GetFaceArray(vtkIdType faceId);
152
158 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
159
160protected:
163
168 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
169
170private:
172 void operator=(const vtkBiQuadraticQuadraticWedge&) = delete;
173};
174//----------------------------------------------------------------------------
175// Return the center of the quadratic wedge in parametric coordinates.
177{
178 pcoords[0] = pcoords[1] = 1. / 3;
179 pcoords[2] = 0.5;
180 return 0;
181}
182
183#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 InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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.
void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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 void InterpolationFunctions(const double pcoords[3], double weights[15])
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.
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
int GetNumberOfEdges() override
Implement the vtkCell API.
int GetCellDimension() override
Implement the vtkCell API.
static vtkBiQuadraticQuadraticWedge * New()
int GetNumberOfFaces() override
Implement the vtkCell API.
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 JacobianInverse(const double pcoords[3], double **inverse, double derivs[45])
Given parametric coordinates compute inverse Jacobian transformation matrix.
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
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.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:140
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:113
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
cell represents a parabolic, isoparametric edge
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:93
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
Definition: vtkCellType.h:117
int vtkIdType
Definition: vtkType.h:332