VTK  9.0.1
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, const double pcoords[3], vtkIdList* pts) override;
73 
75 
79  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
80  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
81  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
82  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
83  double& dist2, double* weights) override;
84  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
85  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
86  void Derivatives(
87  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
88  double* GetParametricCoords() override;
90 
96  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
97  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
98  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
99 
104  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
105  double pcoords[3], int& subId) override;
106 
110  int GetParametricCenter(double pcoords[3]) override;
111 
112  static void InterpolationFunctions(const double pcoords[3], double weights[15]);
113  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
115 
119  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
120  {
122  }
123  void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
124  {
126  }
128 
129 
136  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
137  static const vtkIdType* GetFaceArray(vtkIdType faceId);
139 
145  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
146 
147 protected:
149  ~vtkQuadraticLinearWedge() override;
150 
156  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
157 
158 private:
160  void operator=(const vtkQuadraticLinearWedge&) = delete;
161 };
162 //----------------------------------------------------------------------------
163 // Return the center of the quadratic wedge in parametric coordinates.
164 inline int vtkQuadraticLinearWedge::GetParametricCenter(double pcoords[3])
165 {
166  pcoords[0] = pcoords[1] = 1. / 3;
167  pcoords[2] = 0.5;
168  return 0;
169 }
170 
171 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:33
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
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:76
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkX3D::value
@ value
Definition: vtkX3D.h:226
vtkIdType
int vtkIdType
Definition: vtkType.h:338
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:123
vtkQuadraticTriangle
cell represents a parabolic, isoparametric triangle
Definition: vtkQuadraticTriangle.h:43
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
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:29
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:43
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:56
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
vtkQuadraticLinearWedge::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkQuadraticLinearWedge.h:67
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:33
vtkQuadraticLinearWedge::Face
vtkQuadraticLinearQuad * Face
Definition: vtkQuadraticLinearWedge.h:154
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:179
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:30
vtkQuadraticLinearWedge
cell represents a, 12-node isoparametric wedge
Definition: vtkQuadraticLinearWedge.h:52
vtkQuadraticLinearWedge::Scalars
vtkDoubleArray * Scalars
Definition: vtkQuadraticLinearWedge.h:156
vtkQuadraticLinearWedge::QuadEdge
vtkQuadraticEdge * QuadEdge
Definition: vtkQuadraticLinearWedge.h:151
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:66
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:47
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:153
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:65
vtkQuadraticLinearWedge::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkQuadraticLinearWedge.h:64
vtkQuadraticLinearWedge::Edge
vtkLine * Edge
Definition: vtkQuadraticLinearWedge.h:152
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:164
vtkQuadraticLinearWedge::InterpolationDerivs
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:35
vtkX3D::index
@ index
Definition: vtkX3D.h:252
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:40
vtkQuadraticLinearWedge::InterpolateFunctions
void InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkQuadraticLinearWedge.h:119
vtkQuadraticLinearWedge::Wedge
vtkWedge * Wedge
Definition: vtkQuadraticLinearWedge.h:155