VTK
vtkQuadraticLinearQuad.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticLinearQuad.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 =========================================================================*/
39 #ifndef vtkQuadraticLinearQuad_h
40 #define vtkQuadraticLinearQuad_h
41 
42 #include "vtkCommonDataModelModule.h" // For export macro
43 #include "vtkNonLinearCell.h"
44 
45 class vtkQuadraticEdge;
46 class vtkLine;
47 class vtkQuad;
48 class vtkDoubleArray;
49 
50 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearQuad : public vtkNonLinearCell
51 {
52 public:
53  static vtkQuadraticLinearQuad *New();
55  void PrintSelf(ostream & os, vtkIndent indent) VTK_OVERRIDE;
56 
58 
62  int GetCellType() VTK_OVERRIDE { return VTK_QUADRATIC_LINEAR_QUAD; };
63  int GetCellDimension() VTK_OVERRIDE { return 2; }
64  int GetNumberOfEdges() VTK_OVERRIDE { return 4; }
65  int GetNumberOfFaces() VTK_OVERRIDE { return 0; }
66  vtkCell *GetEdge (int) VTK_OVERRIDE;
67  vtkCell *GetFace (int) VTK_OVERRIDE { return 0; }
69 
70  int CellBoundary (int subId, double pcoords[3], vtkIdList * pts) VTK_OVERRIDE;
71  void Contour (double value, vtkDataArray * cellScalars,
72  vtkIncrementalPointLocator * locator, vtkCellArray * verts,
73  vtkCellArray * lines, vtkCellArray * polys,
74  vtkPointData * inPd, vtkPointData * outPd, vtkCellData * inCd,
75  vtkIdType cellId, vtkCellData * outCd) VTK_OVERRIDE;
76  int EvaluatePosition (double x[3], double *closestPoint,
77  int &subId, double pcoords[3], double &dist2, double *weights) VTK_OVERRIDE;
78  void EvaluateLocation (int &subId, double pcoords[3], double x[3],
79  double *weights) VTK_OVERRIDE;
80  int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts) VTK_OVERRIDE;
81  void Derivatives (int subId, double pcoords[3], double *values, int dim,
82  double *derivs) VTK_OVERRIDE;
83  double *GetParametricCoords () VTK_OVERRIDE;
84 
89  void Clip (double value, vtkDataArray * cellScalars,
90  vtkIncrementalPointLocator * locator, vtkCellArray * polys,
91  vtkPointData * inPd, vtkPointData * outPd,
92  vtkCellData * inCd, vtkIdType cellId, vtkCellData * outCd,
93  int insideOut) VTK_OVERRIDE;
94 
99  int IntersectWithLine (double p1[3], double p2[3], double tol, double &t,
100  double x[3], double pcoords[3], int &subId) VTK_OVERRIDE;
101 
105  int GetParametricCenter(double pcoords[3]) VTK_OVERRIDE;
106 
110  static void InterpolationFunctions (double pcoords[3], double weights[6]);
114  static void InterpolationDerivs (double pcoords[3], double derivs[12]);
116 
120  void InterpolateFunctions (double pcoords[3], double weights[6]) VTK_OVERRIDE
121  {
123  }
124  void InterpolateDerivs (double pcoords[3], double derivs[12]) VTK_OVERRIDE
125  {
127  }
129 
133  static int *GetEdgeArray(int edgeId);
134 
135 protected:
137  ~vtkQuadraticLinearQuad () VTK_OVERRIDE;
138 
140  vtkLine *LinEdge;
141  vtkQuad *Quad;
142  vtkDoubleArray *Scalars;
143 
144 private:
145  vtkQuadraticLinearQuad (const vtkQuadraticLinearQuad &) VTK_DELETE_FUNCTION;
146  void operator = (const vtkQuadraticLinearQuad &) VTK_DELETE_FUNCTION;
147 };
148 //----------------------------------------------------------------------------
149 inline int vtkQuadraticLinearQuad::GetParametricCenter(double pcoords[3])
150 {
151  pcoords[0] = pcoords[1] = 0.5;
152  pcoords[2] = 0.;
153  return 0;
154 }
155 
156 #endif
int GetCellType() override
Implement the vtkCell API.
represent and manipulate point attribute data
Definition: vtkPointData.h:37
virtual double * GetParametricCoords()
Return a contiguous array of parametric coordinates of the points defining this cell.
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
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:41
int GetCellDimension() override
Implement the vtkCell API.
int vtkIdType
Definition: vtkType.h:287
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
void InterpolateDerivs(double pcoords[3], double derivs[12]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
cell represents a 1D line
Definition: vtkLine.h:35
abstract class to specify cell behavior
Definition: vtkCell.h:59
a simple class to control print indentation
Definition: vtkIndent.h:39
int GetNumberOfFaces() override
Implement the vtkCell API.
list of point or cell ids
Definition: vtkIdList.h:36
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:54
static void InterpolationDerivs(double pcoords[3], double derivs[12])
object to represent cell connectivity
Definition: vtkCellArray.h:50
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 GetNumberOfEdges() override
Implement the vtkCell API.
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...
static void InterpolationFunctions(double pcoords[3], double weights[6])
cell represents a quadratic-linear, 6-node isoparametric quad
represent and manipulate 3D points
Definition: vtkPoints.h:39