VTK
vtkBiQuadraticQuad.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkBiQuadraticQuad.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 =========================================================================*/
40 #ifndef vtkBiQuadraticQuad_h
41 #define vtkBiQuadraticQuad_h
42 
43 #include "vtkCommonDataModelModule.h" // For export macro
44 #include "vtkNonLinearCell.h"
45 
46 class vtkQuadraticEdge;
47 class vtkQuad;
48 class vtkTriangle;
49 class vtkDoubleArray;
50 
51 class VTKCOMMONDATAMODEL_EXPORT vtkBiQuadraticQuad : public vtkNonLinearCell
52 {
53 public:
54  static vtkBiQuadraticQuad *New ();
56  void PrintSelf (ostream & os, vtkIndent indent) override;
57 
62  int GetCellType() override { return VTK_BIQUADRATIC_QUAD; }
63  int GetCellDimension() override { return 2; }
64  int GetNumberOfEdges() override { return 4; }
65  int GetNumberOfFaces() override { return 0; }
66  vtkCell *GetEdge (int) override;
67  vtkCell *GetFace (int) override { return nullptr; }
68 
69  int CellBoundary (int subId, double pcoords[3], vtkIdList * pts) override;
70  int EvaluatePosition (double x[3], double *closestPoint,
71  int &subId, double pcoords[3],
72  double &dist2, double *weights) override;
73  void EvaluateLocation (int &subId, double pcoords[3], double x[3],
74  double *weights) override;
75  int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts) override;
76  void Derivatives (int subId, double pcoords[3], double *values,
77  int dim, double *derivs) override;
78  double *GetParametricCoords() override;
79 
80  void Contour (double value, vtkDataArray * cellScalars,
81  vtkIncrementalPointLocator * locator, vtkCellArray * verts,
82  vtkCellArray * lines, vtkCellArray * polys,
83  vtkPointData * inPd, vtkPointData * outPd, vtkCellData * inCd,
84  vtkIdType cellId, vtkCellData * outCd) override;
85 
90  void Clip (double value, vtkDataArray * cellScalars,
91  vtkIncrementalPointLocator * locator, vtkCellArray * polys,
92  vtkPointData * inPd, vtkPointData * outPd,
93  vtkCellData * inCd, vtkIdType cellId, vtkCellData * outCd,
94  int insideOut) override;
95 
100  int IntersectWithLine (double p1[3], double p2[3], double tol, double &t,
101  double x[3], double pcoords[3],
102  int &subId) override;
103 
107  int GetParametricCenter(double pcoords[3]) override;
108 
109  void InterpolateFunctions (double pcoords[3], double weights[9]) override
110  {
111  vtkBiQuadraticQuad::InterpolationFunctionsPrivate(pcoords,weights);
112  }
113  void InterpolateDerivs (double pcoords[3], double derivs[18]) override
114  {
115  vtkBiQuadraticQuad::InterpolationDerivsPrivate(pcoords,derivs);
116  }
118 
119 protected:
121  ~vtkBiQuadraticQuad() override;
122 
127 
128 private:
129  vtkBiQuadraticQuad(const vtkBiQuadraticQuad&) = delete;
130  void operator=(const vtkBiQuadraticQuad&) = delete;
131 
132  static void InterpolationFunctionsPrivate (double pcoords[3], double weights[9]);
133  static void InterpolationDerivsPrivate (double pcoords[3], double derivs[18]);
134 };
135 //----------------------------------------------------------------------------
136 inline int vtkBiQuadraticQuad::GetParametricCenter(double pcoords[3])
137 {
138  pcoords[0] = pcoords[1] = 0.5;
139  pcoords[2] = 0.;
140  return 0;
141 }
142 
143 #endif
void InterpolateDerivs(double pcoords[3], double derivs[18]) override
represent and manipulate point attribute data
Definition: vtkPointData.h:31
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
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.
virtual void EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
virtual int EvaluatePosition(double x[3], double *closestPoint, 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...
abstract superclass for non-linear cells
a cell that represents a 2D quadrilateral
Definition: vtkQuad.h:35
int vtkIdType
Definition: vtkType.h:345
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
dynamic, self-adjusting array of double
abstract class to specify cell behavior
Definition: vtkCell.h:56
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
list of point or cell ids
Definition: vtkIdList.h:30
virtual void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
virtual int IntersectWithLine(double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
vtkQuadraticEdge * Edge
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 GetNumberOfEdges() override
Return the number of edges in the cell.
object to represent cell connectivity
Definition: vtkCellArray.h:44
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
cell represents a parabolic, isoparametric edge
a cell that represents a triangle
Definition: vtkTriangle.h:35
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.
void InterpolateFunctions(double pcoords[3], double weights[9]) override
virtual int CellBoundary(int subId, double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
int GetNumberOfFaces() override
Return the number of faces in the cell.
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.
int GetCellType() override
Implement the vtkCell API.
vtkDoubleArray * Scalars
represent and manipulate 3D points
Definition: vtkPoints.h:33
int GetParametricCenter(double pcoords[3]) override
Return the center of the pyramid in parametric coordinates.