VTK  9.2.5
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 =========================================================================*/
61 #ifndef vtkQuadraticLinearWedge_h
62 #define vtkQuadraticLinearWedge_h
63 
64 #include "vtkCommonDataModelModule.h" // For export macro
65 #include "vtkNonLinearCell.h"
66 
67 class vtkQuadraticEdge;
68 class vtkLine;
71 class vtkWedge;
72 class vtkDoubleArray;
73 
74 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearWedge : public vtkNonLinearCell
75 {
76 public:
79  void PrintSelf(ostream& os, vtkIndent indent) override;
80 
82 
86  int GetCellType() override { return VTK_QUADRATIC_LINEAR_WEDGE; }
87  int GetCellDimension() override { return 3; }
88  int GetNumberOfEdges() override { return 9; }
89  int GetNumberOfFaces() override { return 5; }
90  vtkCell* GetEdge(int edgeId) override;
91  vtkCell* GetFace(int faceId) override;
93 
94  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
95 
97 
101  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
102  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
103  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
104  int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
105  double& dist2, double* weights) override;
106  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
107  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
109  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
110  double* GetParametricCoords() override;
112 
118  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
119  vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
120  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
121 
126  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
127  double pcoords[3], int& subId) override;
128 
132  int GetParametricCenter(double pcoords[3]) override;
133 
134  static void InterpolationFunctions(const double pcoords[3], double weights[15]);
135  static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
137 
141  void InterpolateFunctions(const double pcoords[3], double weights[15]) override
142  {
144  }
145  void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
146  {
148  }
151 
158  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
159  static const vtkIdType* GetFaceArray(vtkIdType faceId);
161 
167  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
168 
169 protected:
172 
178  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
179 
180 private:
182  void operator=(const vtkQuadraticLinearWedge&) = delete;
183 };
184 //----------------------------------------------------------------------------
185 // Return the center of the quadratic wedge in parametric coordinates.
186 inline int vtkQuadraticLinearWedge::GetParametricCenter(double pcoords[3])
187 {
188  pcoords[0] = pcoords[1] = 1. / 3;
189  pcoords[2] = 0.5;
190  return 0;
191 }
192 
193 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:296
represent and manipulate cell attribute data
Definition: vtkCellData.h:151
abstract class to specify cell behavior
Definition: vtkCell.h:150
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:165
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:143
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:119
cell represents a 1D line
Definition: vtkLine.h:143
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:151
represent and manipulate 3D points
Definition: vtkPoints.h:149
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
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
double * GetParametricCoords() override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘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
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
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
Implement the vtkCell API.
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
The quadratic linear wedge is split into 4 linear wedges, each of them is contoured by a provided sca...
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
Implement the vtkCell API.
int GetCellType() override
Implement the vtkCell API.
int GetCellDimension() override
Implement the vtkCell API.
int GetNumberOfEdges() override
Implement the vtkCell API.
void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
cell represents a parabolic, isoparametric triangle
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:96
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:116
int vtkIdType
Definition: vtkType.h:332