VTK  9.1.0
vtkTriQuadraticPyramid.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkTriQuadraticPyramid.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=========================================================================*/
90#ifndef vtkTriQuadraticPyramid_h
91#define vtkTriQuadraticPyramid_h
92
93#include "vtkCommonDataModelModule.h" // For export macro
94#include "vtkNew.h" // initialize cells that are used for the implementation
95#include "vtkNonLinearCell.h"
96
100class vtkTetra;
101class vtkPyramid;
102class vtkDoubleArray;
103
104class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticPyramid : public vtkNonLinearCell
105{
106public:
109 void PrintSelf(ostream& os, vtkIndent indent) override;
110
112
116 int GetCellType() override { return VTK_TRIQUADRATIC_PYRAMID; }
117 int GetCellDimension() override { return 3; }
118 int GetNumberOfEdges() override { return 8; }
119 int GetNumberOfFaces() override { return 5; }
120 vtkCell* GetEdge(int edgeId) override;
121 vtkCell* GetFace(int faceId) override;
123
124 int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
125 void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
126 vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
127 vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
128 int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
129 double& dist2, double weights[]) override;
130 void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
131
136 int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
137 double pcoords[3], int& subId) override;
138
139 int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
141 int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
142 double* GetParametricCoords() override;
143
149 void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
150 vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
151 vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
152
156 int GetParametricCenter(double pcoords[3]) override;
157
162 double GetParametricDistance(const double pcoords[3]) override;
163
164 static void InterpolationFunctions(const double pcoords[3], double weights[19]);
165 static void InterpolationDerivs(const double pcoords[3], double derivs[57]);
167
171 void InterpolateFunctions(const double pcoords[3], double weights[19]) override
172 {
174 }
175 void InterpolateDerivs(const double pcoords[3], double derivs[57]) override
176 {
178 }
180
186 void JacobianInverse(const double pcoords[3], double** inverse, double derivs[57]);
187
189
196 static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
197 static const vtkIdType* GetFaceArray(vtkIdType faceId);
199
200protected:
203
210 vtkNew<vtkDoubleArray> Scalars; // used to avoid New/Delete in contouring/clipping
211
212private:
214 void operator=(const vtkTriQuadraticPyramid&) = delete;
215};
216//----------------------------------------------------------------------------
217// Return the center of the tri-quadratic pyramid in parametric coordinates.
218//
220{
221 pcoords[0] = pcoords[1] = 0.5;
222 // This is different compared to the last node, because the last node
223 // is the centroid of the nodes 0-4, and not the centroid of the nodes 0-17.
224 // So pcoords[2] is defined as followed to pass the requirement of TestGenericCell
225 pcoords[2] = 283.0 / 456.0;
226 return 0;
227}
228
229#endif
cell represents a parabolic, 9-node isoparametric quad
cell represents a parabolic, isoparametric triangle
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
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
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:103
cell represents a parabolic, isoparametric edge
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:111
cell represents a parabolic, 13-node isoparametric pyramid
static const vtkIdType * GetFaceArray(vtkIdType faceId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
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.
vtkNew< vtkPyramid > Pyramid
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.
@ value
Definition: vtkX3D.h:226
@ index
Definition: vtkX3D.h:252
@ VTK_TRIQUADRATIC_PYRAMID
Definition: vtkCellType.h:114
int vtkIdType
Definition: vtkType.h:332