VTK  9.4.20250201
vtkTriQuadraticPyramid.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
81#ifndef vtkTriQuadraticPyramid_h
82#define vtkTriQuadraticPyramid_h
83
84#include "vtkCommonDataModelModule.h" // For export macro
85#include "vtkNew.h" // initialize cells that are used for the implementation
86#include "vtkNonLinearCell.h"
87
88VTK_ABI_NAMESPACE_BEGIN
92class vtkTetra;
93class vtkPyramid;
94class vtkDoubleArray;
95
96class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticPyramid : public vtkNonLinearCell
97{
98public:
101 void PrintSelf(ostream& os, vtkIndent indent) override;
102
104
108 int GetCellType() override { return VTK_TRIQUADRATIC_PYRAMID; }
109 int GetCellDimension() override { return 3; }
110 int GetNumberOfEdges() override { return 8; }
111 int GetNumberOfFaces() override { return 5; }
112 vtkCell* GetEdge(int edgeId) override;
113 vtkCell* GetFace(int faceId) override;
115
116 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
117 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
118 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
119 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
120 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
121 double& dist2, double weights[]) override;
122 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
123
128 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
129 double pcoords[3], int& subId) override;
130
131 int TriangulateLocalIds(int index, vtkIdList* ptIds) override;
133 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
134 double* GetParametricCoords() override;
135
141 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
142 vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
143 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
144
148 int GetParametricCenter(double pcoords[3]) override;
149
154 double GetParametricDistance(const double pcoords[3]) override;
155
156 static void InterpolationFunctions(const double pcoords[3], double weights[19]);
157 static void InterpolationDerivs(const double pcoords[3], double derivs[57]);
159
163 void InterpolateFunctions(const double pcoords[3], double weights[19]) override
164 {
166 }
167 void InterpolateDerivs(const double pcoords[3], double derivs[57]) override
168 {
170 }
172
178 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[57]);
179
181
188 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
189 static const vtkIdType* GetFaceArray(vtkIdType faceId);
191
192protected:
195
202 vtkNew<vtkDoubleArray> Scalars; // used to avoid New/Delete in contouring/clipping
203
204private:
206 void operator=(const vtkTriQuadraticPyramid&) = delete;
207};
208//----------------------------------------------------------------------------
209// Return the center of the tri-quadratic pyramid in parametric coordinates.
210//
212{
213 pcoords[0] = pcoords[1] = 0.5;
214 // This is different compared to the last node, because the last node
215 // is the centroid of the nodes 0-4, and not the centroid of the nodes 0-17.
216 // So pcoords[2] is defined as followed to pass the requirement of TestGenericCell
217 pcoords[2] = 283.0 / 456.0;
218 return 0;
219}
220
221VTK_ABI_NAMESPACE_END
222#endif
cell represents a parabolic, 9-node isoparametric quad
cell represents a parabolic, isoparametric triangle
object to represent cell connectivity
represent and manipulate cell attribute data
abstract class to specify cell behavior
Definition vtkCell.h:130
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
abstract superclass for arrays of numeric data
dynamic, self-adjusting array of double
list of point or cell ids
Definition vtkIdList.h:133
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
Allocate and hold a VTK object.
Definition vtkNew.h:167
abstract superclass for non-linear cells
represent and manipulate point attribute data
a 3D cell that represents a linear pyramid
Definition vtkPyramid.h:95
cell represents a parabolic, isoparametric edge
a 3D cell that represents a tetrahedron
Definition vtkTetra.h:113
cell represents a parabolic, 19-node isoparametric pyramid
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static void InterpolationDerivs(const double pcoords[3], double derivs[57])
int GetNumberOfEdges() override
Implement the vtkCell API.
vtkCell * GetFace(int faceId) override
Implement the vtkCell API.
int GetCellDimension() override
Implement the vtkCell API.
vtkNew< vtkDoubleArray > Scalars
int GetCellType() override
Implement the vtkCell API.
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
Generate contouring primitives.
vtkCell * GetEdge(int edgeId) override
Implement the vtkCell API.
double GetParametricDistance(const double pcoords[3]) override
Return the distance of the parametric coordinate provided to the cell.
void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkNew< vtkBiQuadraticTriangle > TriangleFace
void InterpolateDerivs(const double pcoords[3], double derivs[57]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
vtkNew< vtkBiQuadraticQuad > QuadFace
vtkNew< vtkQuadraticEdge > Edge
~vtkTriQuadraticPyramid() override
void InterpolateFunctions(const double pcoords[3], double weights[19]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
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...
static void InterpolationFunctions(const double pcoords[3], double weights[19])
static const vtkIdType * GetEdgeArray(vtkIdType edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
static vtkTriQuadraticPyramid * New()
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.
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
void JacobianInverse(const double pcoords[3], double **inverse, double derivs[57])
Given parametric coordinates compute inverse Jacobian transformation matrix.
int GetNumberOfFaces() override
Implement the vtkCell API.
int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[]) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
vtkNew< vtkBiQuadraticTriangle > TriangleFace2
int GetParametricCenter(double pcoords[3]) override
Return the center of the tri-quadratic pyramid in parametric coordinates.
int TriangulateLocalIds(int index, vtkIdList *ptIds) override
Generate simplices of proper dimension.
void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *tets, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut) override
Clip this quadratic triangle using scalar value provided.
@ VTK_TRIQUADRATIC_PYRAMID
Definition vtkCellType.h:66
int vtkIdType
Definition vtkType.h:315