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 =========================================================================*/
39 #ifndef vtkQuadraticLinearWedge_h
40 #define vtkQuadraticLinearWedge_h
41 
42 #include "vtkCommonDataModelModule.h" // For export macro
43 #include "vtkNonLinearCell.h"
44 
45 class vtkQuadraticEdge;
46 class vtkLine;
49 class vtkWedge;
50 class vtkDoubleArray;
51 
52 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearWedge : public vtkNonLinearCell
53 {
54 public:
55  static vtkQuadraticLinearWedge *New ();
57  void PrintSelf (ostream & os, vtkIndent indent) override;
58 
60 
64  int GetCellType() override { return VTK_QUADRATIC_LINEAR_WEDGE; }
65  int GetCellDimension() override { return 3; }
66  int GetNumberOfEdges() override { return 9; }
67  int GetNumberOfFaces() override { return 5; }
68  vtkCell *GetEdge (int edgeId) override;
69  vtkCell *GetFace (int faceId) override;
71 
72  int CellBoundary (int subId, double pcoords[3], vtkIdList * pts) override;
73 
75 
79  void Contour (double value, vtkDataArray * cellScalars,
80  vtkIncrementalPointLocator * locator, vtkCellArray * verts,
81  vtkCellArray * lines, vtkCellArray * polys,
82  vtkPointData * inPd, vtkPointData * outPd, vtkCellData * inCd,
83  vtkIdType cellId, vtkCellData * outCd) override;
84  int EvaluatePosition (double x[3], double *closestPoint,
85  int &subId, double pcoords[3], double &dist2, double *weights) override;
86  void EvaluateLocation (int &subId, double pcoords[3], double x[3],
87  double *weights) override;
88  int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts) override;
89  void Derivatives (int subId, double pcoords[3], double *values,
90  int dim, double *derivs) override;
91  double *GetParametricCoords () override;
93 
99  void Clip (double value, vtkDataArray * cellScalars,
100  vtkIncrementalPointLocator * locator, vtkCellArray * tetras,
101  vtkPointData * inPd, vtkPointData * outPd,
102  vtkCellData * inCd, vtkIdType cellId, vtkCellData * outCd,
103  int insideOut) override;
104 
109  int IntersectWithLine (double p1[3], double p2[3], double tol, double &t,
110  double x[3], double pcoords[3], int &subId) override;
111 
115  int GetParametricCenter (double pcoords[3]) override;
116 
120  static void InterpolationFunctions (double pcoords[3], double weights[15]);
124  static void InterpolationDerivs (double pcoords[3], double derivs[45]);
126 
130  void InterpolateFunctions (double pcoords[3], double weights[15]) override
131  {
133  }
134  void InterpolateDerivs (double pcoords[3], double derivs[45]) override
135  {
137  }
139 
140 
144  static int *GetEdgeArray(int edgeId);
145  static int *GetFaceArray(int faceId);
147 
153  void JacobianInverse (double pcoords[3], double **inverse, double derivs[45]);
154 
155 protected:
157  ~vtkQuadraticLinearWedge () override;
158 
164  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
165 
166 private:
168  void operator = (const vtkQuadraticLinearWedge &) = delete;
169 };
170 //----------------------------------------------------------------------------
171 // Return the center of the quadratic wedge in parametric coordinates.
172 inline int vtkQuadraticLinearWedge::GetParametricCenter(double pcoords[3])
173 {
174  pcoords[0] = pcoords[1] = 1./3;
175  pcoords[2] = 0.5;
176  return 0;
177 }
178 
179 
180 #endif
represent and manipulate point attribute data
Definition: vtkPointData.h:31
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
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:345
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:29
abstract class to specify cell behavior
Definition: vtkCell.h:56
int GetNumberOfFaces() override
Implement the vtkCell API.
a simple class to control print indentation
Definition: vtkIndent.h:33
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.
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.
void InterpolateFunctions(double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
vtkQuadraticLinearQuad * Face
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.
int GetCellDimension() override
Implement the vtkCell API.
cell represents a parabolic, isoparametric edge
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic linear wedge in parametric coordinates.
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
vtkQuadraticTriangle * TriangleFace
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...
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 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:43
represent and manipulate 3D points
Definition: vtkPoints.h:33