VTK  9.1.0
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=========================================================================*/
58#ifndef vtkQuadraticLinearWedge_h
59#define vtkQuadraticLinearWedge_h
60
61#include "vtkCommonDataModelModule.h" // For export macro
62#include "vtkNonLinearCell.h"
63
65class vtkLine;
68class vtkWedge;
69class vtkDoubleArray;
70
71class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticLinearWedge : public vtkNonLinearCell
72{
73public:
76 void PrintSelf(ostream& os, vtkIndent indent) override;
77
79
83 int GetCellType() override { return VTK_QUADRATIC_LINEAR_WEDGE; }
84 int GetCellDimension() override { return 3; }
85 int GetNumberOfEdges() override { return 9; }
86 int GetNumberOfFaces() override { return 5; }
87 vtkCell* GetEdge(int edgeId) override;
88 vtkCell* GetFace(int faceId) override;
90
91 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
92
94
98 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
99 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
100 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
101 int EvaluatePosition(const double x[3], double* closestPoint, int& subId, double pcoords[3],
102 double& dist2, double* weights) override;
103 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
104 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
106 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
107 double* GetParametricCoords() override;
109
115 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
116 vtkCellArray* tetras, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
117 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
118
123 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
124 double pcoords[3], int& subId) override;
125
129 int GetParametricCenter(double pcoords[3]) override;
130
131 static void InterpolationFunctions(const double pcoords[3], double weights[15]);
132 static void InterpolationDerivs(const double pcoords[3], double derivs[45]);
134
138 void InterpolateFunctions(const double pcoords[3], double weights[15]) override
139 {
141 }
142 void InterpolateDerivs(const double pcoords[3], double derivs[45]) override
143 {
145 }
148
155 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
156 static const vtkIdType* GetFaceArray(vtkIdType faceId);
158
164 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[45]);
165
166protected:
169
175 vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
176
177private:
179 void operator=(const vtkQuadraticLinearWedge&) = delete;
180};
181//----------------------------------------------------------------------------
182// Return the center of the quadratic wedge in parametric coordinates.
184{
185 pcoords[0] = pcoords[1] = 1. / 3;
186 pcoords[2] = 0.5;
187 return 0;
188}
189
190#endif
object to represent cell connectivity
Definition: vtkCellArray.h:290
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
abstract class to specify cell behavior
Definition: vtkCell.h:147
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:159
dynamic, self-adjusting array of double
list of point or cell ids
Definition: vtkIdList.h:140
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:113
cell represents a 1D line
Definition: vtkLine.h:140
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:142
represent and manipulate 3D points
Definition: vtkPoints.h:143
cell represents a parabolic, isoparametric edge
cell represents a quadratic-linear, 6-node isoparametric quad
cell represents a, 12-node isoparametric wedge
vtkQuadraticTriangle * TriangleFace
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...
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
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 * GetEdgeArray(vtkIdType edgeId)
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...
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])
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static vtkQuadraticLinearWedge * New()
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:93
@ VTK_QUADRATIC_LINEAR_WEDGE
Definition: vtkCellType.h:116
int vtkIdType
Definition: vtkType.h:332