VTK  9.2.5
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 =========================================================================*/
63 #ifndef vtkBiQuadraticQuadraticWedge_h
64 #define vtkBiQuadraticQuadraticWedge_h
65 
66 #include "vtkCommonDataModelModule.h" // For export macro
67 #include "vtkNonLinearCell.h"
68 
69 class vtkQuadraticEdge;
70 class vtkBiQuadraticQuad;
72 class vtkWedge;
73 class vtkDoubleArray;
74 
75 class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticWedge : public vtkNonLinearCell
76 {
77 public:
80  void PrintSelf(ostream& os, vtkIndent indent) override;
81 
83 
87  int GetCellType() override { return VTK_BIQUADRATIC_QUADRATIC_WEDGE; }
88  int GetCellDimension() override { return 3; }
89  int GetNumberOfEdges() override { return 9; }
90  int GetNumberOfFaces() override { return 5; }
91  vtkCell* GetEdge(int edgeId) override;
92  vtkCell* GetFace(int faceId) override;
94 
95  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
96  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
97  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
98  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
99  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
100  double& dist2, double* weights) override;
101  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
102  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
104  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
105  double* GetParametricCoords() override;
106 
112  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
113  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
114  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
115 
120  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
121  double pcoords[3], int& subId) override;
122 
126  int GetParametricCenter(double pcoords[3]) override;
127 
128  static void InterpolationFunctions(const double pcoords[3], double weights[15]);
129  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
131 
135  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
136  {
138  }
139  void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
140  {
142  }
145 
152  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
153  static const vtkIdType* GetFaceArray(vtkIdType faceId);
155 
161  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
162 
163 protected:
166 
171  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
172 
173 private:
175  void operator=(const vtkBiQuadraticQuadraticWedge&) = delete;
176 };
177 //----------------------------------------------------------------------------
178 // Return the center of the quadratic wedge in parametric coordinates.
180 {
181  pcoords[0] = pcoords[1] = 1. / 3;
182  pcoords[2] = 0.5;
183  return 0;
184 }
185 
186 #endif
cell represents a parabolic, 9-node isoparametric quad
cell represents a parabolic, 18-node isoparametric wedge
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)
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
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’).
static vtkBiQuadraticQuadraticWedge * New()
static void InterpolationFunctions(const double pcoords[3], double weights[15])
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
int GetNumberOfEdges() override
Implement the vtkCell API.
int GetCellDimension() override
Implement the vtkCell API.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
int GetNumberOfFaces() override
Implement the vtkCell API.
int GetCellType() override
Implement the vtkCell API.
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
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
~vtkBiQuadraticQuadraticWedge() override
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:296
represent and manipulate cell attribute data
Definition: vtkCellData.h:151
abstract class to specify cell behavior
Definition: vtkCell.h:150
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:143
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:119
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:151
represent and manipulate 3D points
Definition: vtkPoints.h:149
cell represents a parabolic, isoparametric edge
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:96
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_BIQUADRATIC_QUADRATIC_WEDGE
Definition: vtkCellType.h:117
int vtkIdType
Definition: vtkType.h:332