VTK  9.2.6
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 =========================================================================*/
93 #ifndef vtkTriQuadraticPyramid_h
94 #define vtkTriQuadraticPyramid_h
95 
96 #include "vtkCommonDataModelModule.h" // For export macro
97 #include "vtkNew.h" // initialize cells that are used for the implementation
98 #include "vtkNonLinearCell.h"
99 
100 class vtkQuadraticEdge;
101 class vtkBiQuadraticQuad;
103 class vtkTetra;
104 class vtkPyramid;
105 class vtkDoubleArray;
106 
107 class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticPyramid : public vtkNonLinearCell
108 {
109 public:
110  static vtkTriQuadraticPyramid* New();
112  void PrintSelf(ostream& os, vtkIndent indent) override;
113 
115 
119  int GetCellType() override { return VTK_TRIQUADRATIC_PYRAMID; }
120  int GetCellDimension() override { return 3; }
121  int GetNumberOfEdges() override { return 8; }
122  int GetNumberOfFaces() override { return 5; }
123  vtkCell* GetEdge(int edgeId) override;
124  vtkCell* GetFace(int faceId) override;
126 
127  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
128  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
129  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
130  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
131  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
132  double& dist2, double weights[]) override;
133  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
134 
139  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
140  double pcoords[3], int& subId) override;
141 
142  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
143  void Derivatives(
144  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
145  double* GetParametricCoords() override;
146 
152  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
153  vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
154  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
155 
159  int GetParametricCenter(double pcoords[3]) override;
160 
165  double GetParametricDistance(const double pcoords[3]) override;
166 
167  static void InterpolationFunctions(const double pcoords[3], double weights[19]);
168  static void InterpolationDerivs(const double pcoords[3], double derivs[57]);
170 
174  void InterpolateFunctions(const double pcoords[3], double weights[19]) override
175  {
177  }
178  void InterpolateDerivs(const double pcoords[3], double derivs[57]) override
179  {
181  }
183 
189  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[57]);
190 
192 
199  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
200  static const vtkIdType* GetFaceArray(vtkIdType faceId);
202 
203 protected:
205  ~vtkTriQuadraticPyramid() override;
206 
213  vtkNew<vtkDoubleArray> Scalars; // used to avoid New/Delete in contouring/clipping
214 
215 private:
217  void operator=(const vtkTriQuadraticPyramid&) = delete;
218 };
219 //----------------------------------------------------------------------------
220 // Return the center of the tri-quadratic pyramid in parametric coordinates.
221 //
222 inline int vtkTriQuadraticPyramid::GetParametricCenter(double pcoords[3])
223 {
224  pcoords[0] = pcoords[1] = 0.5;
225  // This is different compared to the last node, because the last node
226  // is the centroid of the nodes 0-4, and not the centroid of the nodes 0-17.
227  // So pcoords[2] is defined as followed to pass the requirement of TestGenericCell
228  pcoords[2] = 283.0 / 456.0;
229  return 0;
230 }
231 
232 #endif
represent and manipulate point attribute data
Definition: vtkPointData.h:41
vtkNew< vtkDoubleArray > Scalars
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:46
represent and manipulate cell attribute data
Definition: vtkCellData.h:41
cell represents a parabolic, 9-node isoparametric quad
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
void InterpolateDerivs(const double pcoords[3], double derivs[57]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
abstract superclass for non-linear cells
int GetCellDimension() override
Implement the vtkCell API.
cell represents a parabolic, 19-node isoparametric pyramid
int vtkIdType
Definition: vtkType.h:332
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static void InterpolationDerivs(const double pcoords[3], double derivs[57])
virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
dynamic, self-adjusting array of double
int GetCellType() override
Implement the vtkCell API.
vtkNew< vtkBiQuadraticTriangle > TriangleFace2
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:44
virtual double GetParametricDistance(const double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
abstract class to specify cell behavior
Definition: vtkCell.h:60
virtual void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
a simple class to control print indentation
Definition: vtkIndent.h:39
int GetNumberOfFaces() override
Implement the vtkCell API.
list of point or cell ids
Definition: vtkIdList.h:33
int GetNumberOfEdges() override
Implement the vtkCell API.
vtkNew< vtkBiQuadraticTriangle > TriangleFace
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
vtkNew< vtkBiQuadraticQuad > QuadFace
static void InterpolationFunctions(const double pcoords[3], double weights[19])
cell represents a parabolic, isoparametric triangle
vtkNew< vtkPyramid > Pyramid
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
Cut (or clip) the cell based on the input cellScalars and the specified value.
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
virtual int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[])=0
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
object to represent cell connectivity
Definition: vtkCellArray.h:186
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkNew< vtkQuadraticEdge > Edge
cell represents a parabolic, isoparametric edge
virtual void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd)=0
Generate contouring primitives.
virtual void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
int GetParametricCenter(double pcoords[3]) override
Return the center of the tri-quadratic pyramid in parametric coordinates.
represent and manipulate 3D points
Definition: vtkPoints.h:39
void InterpolateFunctions(const double pcoords[3], double weights[19]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...