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) override;
61 
63 
67  int GetCellType() override { return VTK_QUADRATIC_LINEAR_WEDGE; }
68  int GetCellDimension() override { return 3; }
69  int GetNumberOfEdges() override { return 9; }
70  int GetNumberOfFaces() override { return 5; }
71  vtkCell *GetEdge (int edgeId) override;
72  vtkCell *GetFace (int faceId) override;
74 
75  int CellBoundary(int subId, const double pcoords[3], vtkIdList * pts) 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) override;
87  int EvaluatePosition(const double x[3], double *closestPoint,
88  int &subId, double pcoords[3], double &dist2, double *weights) override;
89  void EvaluateLocation(int &subId, const double pcoords[3], double x[3],
90  double *weights) override;
91  int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts) override;
92  void Derivatives(int subId, const double pcoords[3], const double *values,
93  int dim, double *derivs) override;
94  double *GetParametricCoords () 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) override;
107 
112  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t,
113  double x[3], double pcoords[3], int &subId) override;
114 
118  int GetParametricCenter (double pcoords[3]) override;
119 
123  static void InterpolationFunctions(const double pcoords[3], double weights[15]);
127  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
129 
133  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
134  {
136  }
137  void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
138  {
140  }
142 
143 
147  static int *GetEdgeArray(int edgeId);
148  static int *GetFaceArray(int faceId);
150 
156  void JacobianInverse(const double pcoords[3], double **inverse, double derivs[45]);
157 
158 protected:
160  ~vtkQuadraticLinearWedge () override;
161 
167  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
168 
169 private:
171  void operator = (const vtkQuadraticLinearWedge &) = delete;
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
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:39
vtkCell::IntersectWithLine
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.
vtkCell::Contour
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.
VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:75
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:37
vtkX3D::value
Definition: vtkX3D.h:220
vtkIdType
int vtkIdType
Definition: vtkType.h:347
vtkObject::New
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
vtkQuadraticLinearWedge::InterpolateDerivs
void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
Definition: vtkQuadraticLinearWedge.h:137
vtkQuadraticTriangle
cell represents a parabolic, isoparametric triangle
Definition: vtkQuadraticTriangle.h:46
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
vtkCell::EvaluateLocation
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.
vtkLine
cell represents a 1D line
Definition: vtkLine.h:35
vtkCell::Triangulate
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
vtkWedge
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:49
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:59
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
vtkQuadraticLinearWedge::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkQuadraticLinearWedge.h:70
vtkCell::GetFace
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:39
vtkQuadraticLinearWedge::Face
vtkQuadraticLinearQuad * Face
Definition: vtkQuadraticLinearWedge.h:165
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:50
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:36
vtkQuadraticLinearWedge
cell represents a, 12-node isoparametric wedge
Definition: vtkQuadraticLinearWedge.h:55
vtkQuadraticLinearWedge::Scalars
vtkDoubleArray * Scalars
Definition: vtkQuadraticLinearWedge.h:167
vtkQuadraticLinearWedge::QuadEdge
vtkQuadraticEdge * QuadEdge
Definition: vtkQuadraticLinearWedge.h:162
vtkCell::CellBoundary
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...
vtkNonLinearCell.h
vtkCell::EvaluatePosition
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...
vtkCell::GetParametricCoords
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkQuadraticLinearWedge::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkQuadraticLinearWedge.h:69
vtkNonLinearCell
abstract superclass for non-linear cells
Definition: vtkNonLinearCell.h:36
vtkCell::Clip
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.
vtkQuadraticLinearQuad
cell represents a quadratic-linear, 6-node isoparametric quad
Definition: vtkQuadraticLinearQuad.h:50
vtkCell::GetEdge
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkQuadraticLinearWedge::TriangleFace
vtkQuadraticTriangle * TriangleFace
Definition: vtkQuadraticLinearWedge.h:164
vtkNonLinearCell::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkQuadraticLinearWedge::InterpolationFunctions
static void InterpolationFunctions(const double pcoords[3], double weights[15])
vtkQuadraticLinearWedge::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkQuadraticLinearWedge.h:68
vtkQuadraticLinearWedge::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkQuadraticLinearWedge.h:67
vtkQuadraticLinearWedge::Edge
vtkLine * Edge
Definition: vtkQuadraticLinearWedge.h:163
vtkCell::Derivatives
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.
vtkCell::GetParametricCenter
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
vtkQuadraticLinearWedge::GetParametricCenter
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic linear wedge in parametric coordinates.
Definition: vtkQuadraticLinearWedge.h:175
vtkQuadraticLinearWedge::InterpolationDerivs
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:41
vtkX3D::index
Definition: vtkX3D.h:246
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:43
vtkQuadraticLinearWedge::InterpolateFunctions
void InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkQuadraticLinearWedge.h:133
vtkQuadraticLinearWedge::Wedge
vtkWedge * Wedge
Definition: vtkQuadraticLinearWedge.h:166