VTK
vtkQuadraticPyramid.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticPyramid.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 =========================================================================*/
39 #ifndef vtkQuadraticPyramid_h
40 #define vtkQuadraticPyramid_h
41 
42 #include "vtkCommonDataModelModule.h" // For export macro
43 #include "vtkNonLinearCell.h"
44 
45 class vtkQuadraticEdge;
46 class vtkQuadraticQuad;
48 class vtkTetra;
49 class vtkPyramid;
50 class vtkDoubleArray;
51 
52 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
53 {
54 public:
55  static vtkQuadraticPyramid* New();
57  void PrintSelf(ostream& os, vtkIndent indent) override;
58 
60 
64  int GetCellType() override { return VTK_QUADRATIC_PYRAMID; }
65  int GetCellDimension() override { return 3; }
66  int GetNumberOfEdges() override { return 8; }
67  int GetNumberOfFaces() override { return 5; }
68  vtkCell* GetEdge(int edgeId) override;
69  vtkCell* GetFace(int faceId) override;
71 
72  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
73  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
74  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
75  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
76  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
77  double& dist2, double weights[]) override;
78  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
79  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
80  void Derivatives(
81  int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
82  double* GetParametricCoords() override;
83 
89  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
90  vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
91  vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
92 
97  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
98  double pcoords[3], int& subId) override;
99 
103  int GetParametricCenter(double pcoords[3]) override;
104 
108  static void InterpolationFunctions(const double pcoords[3], double weights[13]);
112  static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
114 
118  void InterpolateFunctions(const double pcoords[3], double weights[13]) override
119  {
121  }
122  void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
123  {
125  }
127 
128 
132  static int* GetEdgeArray(int edgeId);
133  static int* GetFaceArray(int faceId);
135 
141  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[39]);
142 
143 protected:
145  ~vtkQuadraticPyramid() override;
146 
155  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping
156 
158 
164  void Subdivide(
165  vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
167 
168 
174  void ResizeArrays(vtkIdType newSize);
176 
177 private:
178  vtkQuadraticPyramid(const vtkQuadraticPyramid&) = delete;
179  void operator=(const vtkQuadraticPyramid&) = delete;
180 };
181 //----------------------------------------------------------------------------
182 // Return the center of the quadratic pyramid in parametric coordinates.
183 //
184 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
185 {
186  pcoords[0] = pcoords[1] = 6.0 / 13.0;
187  pcoords[2] = 3.0 / 13.0;
188  return 0;
189 }
190 
191 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:33
vtkQuadraticPyramid::InterpolationDerivs
static void InterpolationDerivs(const double pcoords[3], double derivs[39])
vtkCell::IntersectWithLine
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.
vtkQuadraticPyramid::Edge
vtkQuadraticEdge * Edge
Definition: vtkQuadraticPyramid.h:147
vtkCell::Contour
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.
vtkQuadraticPyramid::InterpolationFunctions
static void InterpolationFunctions(const double pcoords[3], double weights[13])
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkX3D::value
Definition: vtkX3D.h:226
vtkIdType
int vtkIdType
Definition: vtkType.h:349
vtkQuadraticPyramid::CellData
vtkCellData * CellData
Definition: vtkQuadraticPyramid.h:153
vtkObject::New
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
vtkQuadraticPyramid::GetParametricCenter
int GetParametricCenter(double pcoords[3]) override
Return the center of the quadratic pyramid in parametric coordinates.
Definition: vtkQuadraticPyramid.h:184
vtkQuadraticTriangle
cell represents a parabolic, isoparametric triangle
Definition: vtkQuadraticTriangle.h:43
vtkPyramid
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:43
vtkQuadraticPyramid::Scalars
vtkDoubleArray * Scalars
Definition: vtkQuadraticPyramid.h:155
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
vtkCell::EvaluateLocation
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.
vtkQuadraticQuad
cell represents a parabolic, 8-node isoparametric quad
Definition: vtkQuadraticQuad.h:43
vtkQuadraticPyramid::InterpolateFunctions
void InterpolateFunctions(const double pcoords[3], double weights[13]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives)
Definition: vtkQuadraticPyramid.h:118
vtkQuadraticPyramid::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkQuadraticPyramid.h:66
VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:72
vtkQuadraticPyramid
cell represents a parabolic, 13-node isoparametric pyramid
Definition: vtkQuadraticPyramid.h:52
vtkQuadraticPyramid::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkQuadraticPyramid.h:67
vtkQuadraticPyramid::InterpolateDerivs
void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
Definition: vtkQuadraticPyramid.h:122
vtkCell::Triangulate
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:56
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
vtkQuadraticPyramid::Face
vtkQuadraticQuad * Face
Definition: vtkQuadraticPyramid.h:149
vtkCell::GetFace
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:33
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:163
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkQuadraticPyramid::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkQuadraticPyramid.h:65
vtkQuadraticPyramid::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkQuadraticPyramid.h:64
vtkQuadraticPyramid::PointData
vtkPointData * PointData
Definition: vtkQuadraticPyramid.h:152
vtkCell::CellBoundary
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...
vtkNonLinearCell.h
vtkCell::EvaluatePosition
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...
vtkCell::GetParametricCoords
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkNonLinearCell
abstract superclass for non-linear cells
Definition: vtkNonLinearCell.h:36
vtkCell::Clip
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.
vtkCell::GetEdge
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkQuadraticPyramid::Pyramid
vtkPyramid * Pyramid
Definition: vtkQuadraticPyramid.h:151
vtkNonLinearCell::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkCell::Derivatives
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.
vtkCell::GetParametricCenter
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
vtkQuadraticPyramid::Tetra
vtkTetra * Tetra
Definition: vtkQuadraticPyramid.h:150
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:35
vtkQuadraticPyramid::CellScalars
vtkDoubleArray * CellScalars
Definition: vtkQuadraticPyramid.h:154
vtkTetra
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:41
vtkX3D::index
Definition: vtkX3D.h:252
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:40
vtkQuadraticPyramid::TriangleFace
vtkQuadraticTriangle * TriangleFace
Definition: vtkQuadraticPyramid.h:148