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 =========================================================================*/
41 #ifndef vtkQuadraticPyramid_h
42 #define vtkQuadraticPyramid_h
43 
44 #include "vtkCommonDataModelModule.h" // For export macro
45 #include "vtkNonLinearCell.h"
46 
47 class vtkQuadraticEdge;
48 class vtkQuadraticQuad;
50 class vtkTetra;
51 class vtkPyramid;
52 class vtkDoubleArray;
53 
55 {
56 public:
57  static vtkQuadraticPyramid *New();
59  void PrintSelf(ostream& os, vtkIndent indent);
60 
62 
65  int GetCellDimension() {return 3;}
66  int GetNumberOfEdges() {return 8;}
67  int GetNumberOfFaces() {return 5;}
68  vtkCell *GetEdge(int edgeId);
69  vtkCell *GetFace(int faceId);
71 
72  int CellBoundary(int subId, double pcoords[3], vtkIdList *pts);
73  void Contour(double value, vtkDataArray *cellScalars,
75  vtkCellArray *lines, vtkCellArray *polys,
76  vtkPointData *inPd, vtkPointData *outPd,
77  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd);
78  int EvaluatePosition(double x[3], double* closestPoint,
79  int& subId, double pcoords[3],
80  double& dist2, double *weights);
81  void EvaluateLocation(int& subId, double pcoords[3], double x[3],
82  double *weights);
83  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts);
84  void Derivatives(int subId, double pcoords[3], double *values,
85  int dim, double *derivs);
86  virtual double *GetParametricCoords();
87 
89 
92  void Clip(double value, vtkDataArray *cellScalars,
94  vtkPointData *inPd, vtkPointData *outPd,
95  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
96  int insideOut);
98 
100 
102  int IntersectWithLine(double p1[3], double p2[3], double tol, double& t,
103  double x[3], double pcoords[3], int& subId);
105 
106 
108  int GetParametricCenter(double pcoords[3]);
109 
112  static void InterpolationFunctions(double pcoords[3], double weights[13]);
115  static void InterpolationDerivs(double pcoords[3], double derivs[39]);
117 
119  virtual void InterpolateFunctions(double pcoords[3], double weights[13])
120  {
122  }
123  virtual void InterpolateDerivs(double pcoords[3], double derivs[39])
124  {
126  }
128 
129 
131  static int *GetEdgeArray(int edgeId);
132  static int *GetFaceArray(int faceId);
134 
138  void JacobianInverse(double pcoords[3], double **inverse, double derivs[39]);
139 
140 protected:
143 
152  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
153 
154  void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId,
155  vtkDataArray *cellScalars);
156 
157 private:
158  vtkQuadraticPyramid(const vtkQuadraticPyramid&); // Not implemented.
159  void operator=(const vtkQuadraticPyramid&); // Not implemented.
160 };
161 //----------------------------------------------------------------------------
162 // Return the center of the quadratic pyramid in parametric coordinates.
163 //
164 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
165 {
166  pcoords[0] = pcoords[1] = 6.0/13.0;
167  pcoords[2] = 3.0/13.0;
168  return 0;
169 }
170 
171 
172 #endif
static void InterpolationFunctions(double pcoords[3], double weights[13])
vtkDoubleArray * CellScalars
represent and manipulate point attribute data
Definition: vtkPointData.h:36
cell represents a parabolic, 13-node isoparametric pyramid
virtual double * GetParametricCoords()
a 3D cell that represents a linear pyramid
Definition: vtkPyramid.h:48
represent and manipulate cell attribute data
Definition: vtkCellData.h:37
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
virtual void EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights)=0
vtkQuadraticEdge * Edge
virtual int EvaluatePosition(double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights)=0
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:275
virtual void InterpolateDerivs(double pcoords[3], double derivs[39])
dynamic, self-adjusting array of double
a 3D cell that represents a tetrahedron
Definition: vtkTetra.h:46
abstract class to specify cell behavior
Definition: vtkCell.h:61
static void InterpolationDerivs(double pcoords[3], double derivs[39])
int GetParametricCenter(double pcoords[3])
cell represents a parabolic, 8-node isoparametric quad
a simple class to control print indentation
Definition: vtkIndent.h:38
vtkQuadraticTriangle * TriangleFace
list of point or cell ids
Definition: vtkIdList.h:35
virtual void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs)=0
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
virtual int IntersectWithLine(double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
void PrintSelf(ostream &os, vtkIndent indent)
virtual vtkCell * GetFace(int faceId)=0
object to represent cell connectivity
Definition: vtkCellArray.h:49
virtual vtkCell * GetEdge(int edgeId)=0
vtkQuadraticQuad * Face
cell represents a parabolic, isoparametric edge
virtual void InterpolateFunctions(double pcoords[3], double weights[13])
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
cell represents a parabolic, isoparametric triangle
virtual int CellBoundary(int subId, double pcoords[3], vtkIdList *pts)=0
static vtkObject * New()
vtkDoubleArray * Scalars
virtual int GetParametricCenter(double pcoords[3])
#define VTKCOMMONDATAMODEL_EXPORT
represent and manipulate 3D points
Definition: vtkPoints.h:38