VTK  9.3.1
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
32 #ifndef vtkBiQuadraticQuadraticWedge_h
33 #define vtkBiQuadraticQuadraticWedge_h
34 
35 #include "vtkCommonDataModelModule.h" // For export macro
36 #include "vtkNonLinearCell.h"
37 
38 VTK_ABI_NAMESPACE_BEGIN
39 class vtkQuadraticEdge;
40 class vtkBiQuadraticQuad;
42 class vtkWedge;
43 class vtkDoubleArray;
44 
45 class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticWedge : public vtkNonLinearCell
46 {
47 public:
50  void PrintSelf(ostream& os, vtkIndent indent) override;
51 
53 
57  int GetCellType() override { return VTK_BIQUADRATIC_QUADRATIC_WEDGE; }
58  int GetCellDimension() override { return 3; }
59  int GetNumberOfEdges() override { return 9; }
60  int GetNumberOfFaces() override { return 5; }
61  vtkCell* GetEdge(int edgeId) override;
62  vtkCell* GetFace(int faceId) override;
64 
65  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
66  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
67  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
68  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
69  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
70  double& dist2, double* weights) override;
71  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
72  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
73  void Derivatives(
74  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
75  double* GetParametricCoords() override;
76 
82  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
83  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
84  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
85 
90  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
91  double pcoords[3], int& subId) override;
92 
96  int GetParametricCenter(double pcoords[3]) override;
97 
98  static void InterpolationFunctions(const double pcoords[3], double weights[18]);
99  static void InterpolationDerivs(const double pcoords[3], double derivs[54]);
101 
105  void InterpolateFunctions(const double pcoords[3], double weights[18]) override
106  {
108  }
109  void InterpolateDerivs(const double pcoords[3], double derivs[54]) override
110  {
112  }
115 
122  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
123  static const vtkIdType* GetFaceArray(vtkIdType faceId);
125 
131  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[54]);
132 
133 protected:
135  ~vtkBiQuadraticQuadraticWedge() override;
136 
141  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
142 
143 private:
145  void operator=(const vtkBiQuadraticQuadraticWedge&) = delete;
146 };
147 //----------------------------------------------------------------------------
148 // Return the center of the quadratic wedge in parametric coordinates.
150 {
151  pcoords[0] = pcoords[1] = 1. / 3;
152  pcoords[2] = 0.5;
153  return 0;
154 }
155 
156 VTK_ABI_NAMESPACE_END
157 #endif
void InterpolateFunctions(const double pcoords[3], double weights[18]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
void InterpolateDerivs(const double pcoords[3], double derivs[54]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
represent and manipulate point attribute data
Definition: vtkPointData.h:29
int GetNumberOfEdges() override
Implement the vtkCell API.
represent and manipulate cell attribute data
Definition: vtkCellData.h:30
cell represents a parabolic, 9-node isoparametric quad
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.
cell represents a parabolic, 18-node isoparametric wedge
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.
a simple class to control print indentation
Definition: vtkIndent.h:28
list of point or cell ids
Definition: vtkIdList.h:22
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:44
static void InterpolationDerivs(const double pcoords[3], double derivs[54])
static void InterpolationFunctions(const double pcoords[3], double weights[18])
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.
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic wedge in parametric coordinates.
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.
int GetCellType() override
Implement the vtkCell API.
cell represents a parabolic, isoparametric edge
int GetCellDimension() override
Implement the vtkCell API.
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.
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
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 GetNumberOfFaces() 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