VTK  9.0.1
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 =========================================================================*/
41 #ifndef vtkBiQuadraticQuadraticWedge_h
42 #define vtkBiQuadraticQuadraticWedge_h
43 
44 #include "vtkCommonDataModelModule.h" // For export macro
45 #include "vtkNonLinearCell.h"
46 
47 class vtkQuadraticEdge;
48 class vtkBiQuadraticQuad;
50 class vtkWedge;
51 class vtkDoubleArray;
52 
53 class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuadraticWedge : public vtkNonLinearCell
54 {
55 public:
58  void PrintSelf(ostream& os, vtkIndent indent) override;
59 
61 
65  int GetCellType() override { return VTK_BIQUADRATIC_QUADRATIC_WEDGE; }
66  int GetCellDimension() override { return 3; }
67  int GetNumberOfEdges() override { return 9; }
68  int GetNumberOfFaces() override { return 5; }
69  vtkCell* GetEdge(int edgeId) override;
70  vtkCell* GetFace(int faceId) override;
72 
73  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
74  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
75  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
76  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
77  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
78  double& dist2, double* weights) override;
79  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
80  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
81  void Derivatives(
82  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
83  double* GetParametricCoords() override;
84 
90  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
91  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
92  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
93 
98  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
99  double pcoords[3], int& subId) override;
100 
104  int GetParametricCenter(double pcoords[3]) override;
105 
106  static void InterpolationFunctions(const double pcoords[3], double weights[15]);
107  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
109 
113  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
114  {
116  }
117  void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
118  {
120  }
122 
123 
130  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
131  static const vtkIdType* GetFaceArray(vtkIdType faceId);
133 
139  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
140 
141 protected:
143  ~vtkBiQuadraticQuadraticWedge() override;
144 
149  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
150 
151 private:
153  void operator=(const vtkBiQuadraticQuadraticWedge&) = delete;
154 };
155 //----------------------------------------------------------------------------
156 // Return the center of the quadratic wedge in parametric coordinates.
158 {
159  pcoords[0] = pcoords[1] = 1. / 3;
160  pcoords[2] = 0.5;
161  return 0;
162 }
163 
164 #endif
represent and manipulate point attribute data
Definition: vtkPointData.h:31
int GetNumberOfEdges() override
Implement the vtkCell API.
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
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:338
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:56
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:33
list of point or cell ids
Definition: vtkIdList.h:30
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
void InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
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:179
static void InterpolationFunctions(const double pcoords[3], double weights[15])
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:43
represent and manipulate 3D points
Definition: vtkPoints.h:33