VTK  9.2.6
vtkQuadraticWedge.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticWedge.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 =========================================================================*/
58 #ifndef vtkQuadraticWedge_h
59 #define vtkQuadraticWedge_h
60 
61 #include "vtkCommonDataModelModule.h" // For export macro
62 #include "vtkNonLinearCell.h"
63 
64 class vtkQuadraticEdge;
65 class vtkQuadraticQuad;
67 class vtkWedge;
68 class vtkDoubleArray;
69 
70 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticWedge : public vtkNonLinearCell
71 {
72 public:
75  void PrintSelf(ostream& os, vtkIndent indent) override;
76 
78 
82  int GetCellType() override { return VTK_QUADRATIC_WEDGE; }
83  int GetCellDimension() override { return 3; }
84  int GetNumberOfEdges() override { return 9; }
85  int GetNumberOfFaces() override { return 5; }
86  vtkCell* GetEdge(int edgeId) override;
87  vtkCell* GetFace(int faceId) override;
89 
90  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
91  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
92  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
93  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
94  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
95  double& dist2, double weights[]) override;
96  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
97  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
99  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
100  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[15]);
124  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
126 
130  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
131  {
133  }
134  void InterpolateDerivs(const double pcoords[3], double derivs[45]) 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[45]);
157 
158 protected:
160  ~vtkQuadraticWedge() override;
161 
169  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
170 
171  void Subdivide(
172  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
173 
174 private:
175  vtkQuadraticWedge(const vtkQuadraticWedge&) = delete;
176  void operator=(const vtkQuadraticWedge&) = delete;
177 };
178 //----------------------------------------------------------------------------
179 // Return the center of the quadratic wedge in parametric coordinates.
180 inline int vtkQuadraticWedge::GetParametricCenter(double pcoords[3])
181 {
182  pcoords[0] = pcoords[1] = 1. / 3;
183  pcoords[2] = 0.5;
184  return 0;
185 }
186 
187 #endif
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, 8-node isoparametric quad
cell represents a parabolic, isoparametric triangle
cell represents a parabolic, 15-node isoparametric wedge
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
vtkDoubleArray * Scalars
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
vtkQuadraticTriangle * TriangleFace
static vtkQuadraticWedge * New()
~vtkQuadraticWedge() override
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[45])
Given parametric coordinates compute inverse Jacobian transformation matrix.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, vtkDataArray *cellScalars)
vtkQuadraticEdge * Edge
vtkPointData * PointData
int GetNumberOfEdges() override
Implement the vtkCell API.
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
int GetCellType() override
Implement the vtkCell API.
vtkDoubleArray * CellScalars
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
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.
vtkCellData * CellData
vtkQuadraticQuad * Face
static void InterpolationFunctions(const double pcoords[3], double weights[15])
void InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
int GetCellDimension() override
Implement the vtkCell API.
int GetNumberOfFaces() override
Implement the vtkCell API.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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 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 hexahedron using scalar value provided.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
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.
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:96
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_WEDGE
Definition: vtkCellType.h:110
int vtkIdType
Definition: vtkType.h:332