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 =========================================================================*/
42 #ifndef vtkQuadraticPyramid_h
43 #define vtkQuadraticPyramid_h
44 
45 #include "vtkCommonDataModelModule.h" // For export macro
46 #include "vtkNonLinearCell.h"
47 
48 class vtkQuadraticEdge;
49 class vtkQuadraticQuad;
51 class vtkTetra;
52 class vtkPyramid;
53 class vtkDoubleArray;
54 
55 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
56 {
57 public:
60  void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
61 
63 
67  int GetCellType() VTK_OVERRIDE {return VTK_QUADRATIC_PYRAMID;};
68  int GetCellDimension() VTK_OVERRIDE {return 3;}
69  int GetNumberOfEdges() VTK_OVERRIDE {return 8;}
70  int GetNumberOfFaces() VTK_OVERRIDE {return 5;}
71  vtkCell *GetEdge(int edgeId) VTK_OVERRIDE;
72  vtkCell *GetFace(int faceId) VTK_OVERRIDE;
74 
75  int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) VTK_OVERRIDE;
76  void Contour(double value, vtkDataArray *cellScalars,
78  vtkCellArray *lines, vtkCellArray *polys,
79  vtkPointData *inPd, vtkPointData *outPd,
80  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) VTK_OVERRIDE;
81  int EvaluatePosition(double x[3], double* closestPoint,
82  int& subId, double pcoords[3],
83  double& dist2, double *weights) VTK_OVERRIDE;
84  void EvaluateLocation(int& subId, double pcoords[3], double x[3],
85  double *weights) VTK_OVERRIDE;
86  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) VTK_OVERRIDE;
87  void Derivatives(int subId, double pcoords[3], double *values,
88  int dim, double *derivs) VTK_OVERRIDE;
89  double *GetParametricCoords() VTK_OVERRIDE;
90 
96  void Clip(double value, vtkDataArray *cellScalars,
98  vtkPointData *inPd, vtkPointData *outPd,
99  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
100  int insideOut) VTK_OVERRIDE;
101 
106  int IntersectWithLine(double p1[3], double p2[3], double tol, double& t,
107  double x[3], double pcoords[3], int& subId) VTK_OVERRIDE;
108 
109 
113  int GetParametricCenter(double pcoords[3]) VTK_OVERRIDE;
114 
118  static void InterpolationFunctions(double pcoords[3], double weights[13]);
122  static void InterpolationDerivs(double pcoords[3], double derivs[39]);
124 
128  void InterpolateFunctions(double pcoords[3], double weights[13]) VTK_OVERRIDE
129  {
131  }
132  void InterpolateDerivs(double pcoords[3], double derivs[39]) VTK_OVERRIDE
133  {
135  }
137 
138 
142  static int *GetEdgeArray(int edgeId);
143  static int *GetFaceArray(int faceId);
145 
151  void JacobianInverse(double pcoords[3], double **inverse, double derivs[39]);
152 
153 protected:
155  ~vtkQuadraticPyramid() VTK_OVERRIDE;
156 
158  vtkQuadraticTriangle *TriangleFace;
160  vtkTetra *Tetra;
161  vtkPyramid *Pyramid;
162  vtkPointData *PointData;
163  vtkCellData *CellData;
164  vtkDoubleArray *CellScalars;
165  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
166 
167  void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId,
168  vtkDataArray *cellScalars);
169 
170 private:
171  vtkQuadraticPyramid(const vtkQuadraticPyramid&) VTK_DELETE_FUNCTION;
172  void operator=(const vtkQuadraticPyramid&) VTK_DELETE_FUNCTION;
173 };
174 //----------------------------------------------------------------------------
175 // Return the center of the quadratic pyramid in parametric coordinates.
176 //
177 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
178 {
179  pcoords[0] = pcoords[1] = 6.0/13.0;
180  pcoords[2] = 3.0/13.0;
181  return 0;
182 }
183 
184 
185 #endif
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:40
vtkQuadraticPyramid::JacobianInverse
void JacobianInverse(double pcoords[3], double **inverse, double derivs[39])
Given parametric coordinates compute inverse Jacobian transformation matrix.
vtkQuadraticPyramid::Triangulate
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:38
vtkX3D::value
@ value
Definition: vtkX3D.h:220
vtkIdType
int vtkIdType
Definition: vtkType.h:287
vtkQuadraticPyramid::GetFaceArray
static int * GetFaceArray(int faceId)
vtkQuadraticPyramid::EvaluatePosition
int EvaluatePosition(double x[3], double *closestPoint, 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...
vtkQuadraticTriangle
cell represents a parabolic, isoparametric triangle
Definition: vtkQuadraticTriangle.h:47
vtkPyramid
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:50
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
vtkQuadraticQuad
cell represents a parabolic, 8-node isoparametric quad
Definition: vtkQuadraticQuad.h:47
vtkQuadraticPyramid::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkQuadraticPyramid.h:69
VTK_QUADRATIC_PYRAMID
@ VTK_QUADRATIC_PYRAMID
Definition: vtkCellType.h:71
vtkQuadraticPyramid::InterpolationDerivs
static void InterpolationDerivs(double pcoords[3], double derivs[39])
vtkQuadraticPyramid
cell represents a parabolic, 13-node isoparametric pyramid
Definition: vtkQuadraticPyramid.h:56
vtkQuadraticPyramid::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkQuadraticPyramid.h:70
vtkQuadraticPyramid::CellBoundary
int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:60
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:39
vtkQuadraticPyramid::GetParametricCoords
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkQuadraticPyramid::New
static vtkQuadraticPyramid * New()
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:40
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:51
vtkQuadraticPyramid::EvaluateLocation
void EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:37
vtkQuadraticPyramid::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkQuadraticPyramid::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkQuadraticPyramid.h:68
vtkQuadraticPyramid::GetCellType
int GetCellType() override
Implement the vtkCell API.
Definition: vtkQuadraticPyramid.h:67
vtkNonLinearCell.h
vtkQuadraticPyramid::Derivatives
void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
vtkQuadraticPyramid::InterpolateDerivs
void InterpolateDerivs(double pcoords[3], double derivs[39]) override
Definition: vtkQuadraticPyramid.h:132
vtkNonLinearCell
abstract superclass for non-linear cells
Definition: vtkNonLinearCell.h:37
vtkQuadraticPyramid::vtkQuadraticPyramid
vtkQuadraticPyramid()
vtkQuadraticPyramid::InterpolationFunctions
static void InterpolationFunctions(double pcoords[3], double weights[13])
vtkQuadraticPyramid::~vtkQuadraticPyramid
~vtkQuadraticPyramid() override
vtkDoubleArray
dynamic, self-adjusting array of double
Definition: vtkDoubleArray.h:42
vtkQuadraticPyramid::Contour
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.
vtkQuadraticPyramid::GetFace
vtkCell * GetFace(int faceId) override
Return the face cell from the faceId of the cell.
vtkTetra
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:48
vtkQuadraticPyramid::GetEdgeArray
static int * GetEdgeArray(int edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
vtkX3D::index
@ index
Definition: vtkX3D.h:246
vtkQuadraticEdge
cell represents a parabolic, isoparametric edge
Definition: vtkQuadraticEdge.h:44
vtkQuadraticPyramid::GetEdge
vtkCell * GetEdge(int edgeId) override
Return the edge cell from the edgeId of the cell.