VTK  9.3.1
vtkQuadraticWedge.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
27 #ifndef vtkQuadraticWedge_h
28 #define vtkQuadraticWedge_h
29 
30 #include "vtkCommonDataModelModule.h" // For export macro
31 #include "vtkNonLinearCell.h"
32 
33 VTK_ABI_NAMESPACE_BEGIN
34 class vtkQuadraticEdge;
35 class vtkQuadraticQuad;
37 class vtkWedge;
38 class vtkDoubleArray;
39 
40 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticWedge : public vtkNonLinearCell
41 {
42 public:
43  static vtkQuadraticWedge* New();
45  void PrintSelf(ostream& os, vtkIndent indent) override;
46 
48 
52  int GetCellType() override { return VTK_QUADRATIC_WEDGE; }
53  int GetCellDimension() override { return 3; }
54  int GetNumberOfEdges() override { return 9; }
55  int GetNumberOfFaces() override { return 5; }
56  vtkCell* GetEdge(int edgeId) override;
57  vtkCell* GetFace(int faceId) override;
59 
60  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
61  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
62  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
63  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
64  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
65  double& dist2, double weights[]) override;
66  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
67  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
68  void Derivatives(
69  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
70  double* GetParametricCoords() override;
71 
77  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
78  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
79  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
80 
85  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
86  double pcoords[3], int& subId) override;
87 
91  int GetParametricCenter(double pcoords[3]) override;
92 
93  static void InterpolationFunctions(const double pcoords[3], double weights[15]);
94  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
96 
100  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
101  {
103  }
104  void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
105  {
107  }
110 
117  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
118  static const vtkIdType* GetFaceArray(vtkIdType faceId);
120 
126  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
127 
128 protected:
130  ~vtkQuadraticWedge() override;
131 
139  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
140 
141  void Subdivide(
142  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
143 
144 private:
145  vtkQuadraticWedge(const vtkQuadraticWedge&) = delete;
146  void operator=(const vtkQuadraticWedge&) = delete;
147 };
148 //----------------------------------------------------------------------------
149 // Return the center of the quadratic wedge in parametric coordinates.
150 inline int vtkQuadraticWedge::GetParametricCenter(double pcoords[3])
151 {
152  pcoords[0] = pcoords[1] = 1. / 3;
153  pcoords[2] = 0.5;
154  return 0;
155 }
156 
157 VTK_ABI_NAMESPACE_END
158 #endif
static void InterpolationFunctions(const double pcoords[3], double weights[15])
int GetNumberOfFaces() override
Implement the vtkCell API.
vtkPointData * PointData
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic wedge in parametric coordinates.
vtkDoubleArray * Scalars
vtkQuadraticQuad * Face
represent and manipulate point attribute data
Definition: vtkPointData.h:29
vtkDoubleArray * CellScalars
void InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
represent and manipulate cell attribute data
Definition: vtkCellData.h:30
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:315
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[45]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
dynamic, self-adjusting array of double
abstract class to specify cell behavior
Definition: vtkCell.h:49
virtual void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
int GetNumberOfEdges() override
Implement the vtkCell API.
vtkQuadraticEdge * Edge
cell represents a parabolic, 8-node isoparametric quad
a simple class to control print indentation
Definition: vtkIndent.h:28
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
list of point or cell ids
Definition: vtkIdList.h:22
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:44
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
Cut (or clip) the cell based on the input cellScalars and the specified value.
vtkCellData * CellData
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
virtual int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[])=0
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
object to represent cell connectivity
Definition: vtkCellArray.h:175
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
cell represents a parabolic, isoparametric edge
virtual void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd)=0
Generate contouring primitives.
int GetCellDimension() override
Implement the vtkCell API.
virtual void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
cell represents a parabolic, isoparametric triangle
cell represents a parabolic, 15-node isoparametric wedge
vtkQuadraticTriangle * TriangleFace
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
int GetCellType() override
Implement the vtkCell API.
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:35
represent and manipulate 3D points
Definition: vtkPoints.h:28