VTK  9.0.2
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:
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;
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:
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
object to represent cell connectivity
Definition: vtkCellArray.h:180
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
abstract class to specify cell behavior
Definition: vtkCell.h:57
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:50
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:31
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:34
cell represents a 1D line
Definition: vtkLine.h:30
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:32
represent and manipulate 3D points
Definition: vtkPoints.h:34
cell represents a parabolic, isoparametric edge
cell represents a quadratic-linear, 6-node isoparametric quad
cell represents a, 12-node isoparametric wedge
vtkQuadraticTriangle * TriangleFace
static vtkQuadraticLinearWedge * New()
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic linear wedge in parametric coordinates.
~vtkQuadraticLinearWedge() override
static void InterpolationFunctions(const double pcoords[3], double weights[15])
int EvaluatePosition(const double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights) override
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCell * GetEdge(int edgeId) override
Return the edge cell from the edgeId of the cell.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
static const vtkIdType * GetFaceArray(vtkIdType faceId)
int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[45])
Given parametric coordinates compute inverse Jacobian transformation matrix.
int GetNumberOfFaces() override
Return the number of faces in the cell.
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId) override
Line-edge intersection.
void InterpolateFunctions(const double pcoords[3], double weights[15]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkQuadraticLinearQuad * Face
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tetras, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this quadratic linear wedge using scalar value provided.
static void InterpolationDerivs(const double pcoords[3], double derivs[45])
vtkCell * GetFace(int faceId) override
Return the face cell from the faceId of the cell.
int GetCellType() override
Implement the vtkCell API.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
int GetNumberOfEdges() override
Return the number of edges in the cell.
void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:44
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:76
int vtkIdType
Definition: vtkType.h:338