VTK
vtkQuadraticLinearWedge.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticLinearWedge.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 =========================================================================*/
42 #ifndef vtkQuadraticLinearWedge_h
43 #define vtkQuadraticLinearWedge_h
44 
45 #include "vtkCommonDataModelModule.h" // For export macro
46 #include "vtkNonLinearCell.h"
47 
48 class vtkQuadraticEdge;
49 class vtkLine;
52 class vtkWedge;
53 class vtkDoubleArray;
54 
55 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearWedge : public vtkNonLinearCell
56 {
57 public:
58  static vtkQuadraticLinearWedge *New ();
60  void PrintSelf (ostream & os, vtkIndent indent) VTK_OVERRIDE;
61 
63 
67  int GetCellType() VTK_OVERRIDE { return VTK_QUADRATIC_LINEAR_WEDGE; }
68  int GetCellDimension() VTK_OVERRIDE { return 3; }
69  int GetNumberOfEdges() VTK_OVERRIDE { return 9; }
70  int GetNumberOfFaces() VTK_OVERRIDE { return 5; }
71  vtkCell *GetEdge (int edgeId) VTK_OVERRIDE;
72  vtkCell *GetFace (int faceId) VTK_OVERRIDE;
74 
75  int CellBoundary (int subId, double pcoords[3], vtkIdList * pts) VTK_OVERRIDE;
76 
78 
82  void Contour (double value, vtkDataArray * cellScalars,
83  vtkIncrementalPointLocator * locator, vtkCellArray * verts,
84  vtkCellArray * lines, vtkCellArray * polys,
85  vtkPointData * inPd, vtkPointData * outPd, vtkCellData * inCd,
86  vtkIdType cellId, vtkCellData * outCd) VTK_OVERRIDE;
87  int EvaluatePosition (double x[3], double *closestPoint,
88  int &subId, double pcoords[3], double &dist2, double *weights) VTK_OVERRIDE;
89  void EvaluateLocation (int &subId, double pcoords[3], double x[3],
90  double *weights) VTK_OVERRIDE;
91  int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts) VTK_OVERRIDE;
92  void Derivatives (int subId, double pcoords[3], double *values,
93  int dim, double *derivs) VTK_OVERRIDE;
94  double *GetParametricCoords () VTK_OVERRIDE;
96 
102  void Clip (double value, vtkDataArray * cellScalars,
103  vtkIncrementalPointLocator * locator, vtkCellArray * tetras,
104  vtkPointData * inPd, vtkPointData * outPd,
105  vtkCellData * inCd, vtkIdType cellId, vtkCellData * outCd,
106  int insideOut) VTK_OVERRIDE;
107 
112  int IntersectWithLine (double p1[3], double p2[3], double tol, double &t,
113  double x[3], double pcoords[3], int &subId) VTK_OVERRIDE;
114 
118  int GetParametricCenter (double pcoords[3]) VTK_OVERRIDE;
119 
123  static void InterpolationFunctions (double pcoords[3], double weights[15]);
127  static void InterpolationDerivs (double pcoords[3], double derivs[45]);
129 
133  void InterpolateFunctions (double pcoords[3], double weights[15]) VTK_OVERRIDE
134  {
136  }
137  void InterpolateDerivs (double pcoords[3], double derivs[45]) VTK_OVERRIDE
138  {
140  }
142 
143 
147  static int *GetEdgeArray(int edgeId);
148  static int *GetFaceArray(int faceId);
150 
156  void JacobianInverse (double pcoords[3], double **inverse, double derivs[45]);
157 
158 protected:
160  ~vtkQuadraticLinearWedge () VTK_OVERRIDE;
161 
162  vtkQuadraticEdge *QuadEdge;
163  vtkLine *Edge;
164  vtkQuadraticTriangle *TriangleFace;
166  vtkWedge *Wedge;
167  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
168 
169 private:
170  vtkQuadraticLinearWedge (const vtkQuadraticLinearWedge &) VTK_DELETE_FUNCTION;
171  void operator = (const vtkQuadraticLinearWedge &) VTK_DELETE_FUNCTION;
172 };
173 //----------------------------------------------------------------------------
174 // Return the center of the quadratic wedge in parametric coordinates.
175 inline int vtkQuadraticLinearWedge::GetParametricCenter(double pcoords[3])
176 {
177  pcoords[0] = pcoords[1] = 1./3;
178  pcoords[2] = 0.5;
179  return 0;
180 }
181 
182 
183 #endif
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
int vtkIdType
Definition: vtkType.h:287
static void InterpolationFunctions(double pcoords[3], double weights[15])
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static void InterpolationDerivs(double pcoords[3], double derivs[45])
dynamic, self-adjusting array of double
cell represents a 1D line
Definition: vtkLine.h:35
abstract class to specify cell behavior
Definition: vtkCell.h:59
int GetNumberOfFaces() override
Implement the vtkCell API.
a simple class to control print indentation
Definition: vtkIndent.h:39
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
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
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.
int GetCellDimension() override
Implement the vtkCell API.
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.
cell represents a parabolic, isoparametric triangle
int GetCellType() 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...
int GetNumberOfEdges() override
Implement the vtkCell API.
cell represents a, 12-node isoparametric wedge
cell represents a quadratic-linear, 6-node isoparametric quad
void InterpolateDerivs(double pcoords[3], double derivs[45]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:49
represent and manipulate 3D points
Definition: vtkPoints.h:39