VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: vtkQuadraticPyramid.h 00005 00006 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen 00007 All rights reserved. 00008 See Copyright.txt or http://www.kitware.com/Copyright.htm for details. 00009 00010 This software is distributed WITHOUT ANY WARRANTY; without even 00011 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 00012 PURPOSE. See the above copyright notice for more information. 00013 00014 =========================================================================*/ 00041 #ifndef __vtkQuadraticPyramid_h 00042 #define __vtkQuadraticPyramid_h 00043 00044 #include "vtkNonLinearCell.h" 00045 00046 class vtkQuadraticEdge; 00047 class vtkQuadraticQuad; 00048 class vtkQuadraticTriangle; 00049 class vtkTetra; 00050 class vtkPyramid; 00051 class vtkDoubleArray; 00052 00053 class VTK_FILTERING_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell 00054 { 00055 public: 00056 static vtkQuadraticPyramid *New(); 00057 vtkTypeMacro(vtkQuadraticPyramid,vtkNonLinearCell); 00058 void PrintSelf(ostream& os, vtkIndent indent); 00059 00061 00063 int GetCellType() {return VTK_QUADRATIC_PYRAMID;}; 00064 int GetCellDimension() {return 3;} 00065 int GetNumberOfEdges() {return 8;} 00066 int GetNumberOfFaces() {return 5;} 00067 vtkCell *GetEdge(int edgeId); 00068 vtkCell *GetFace(int faceId); 00070 00071 int CellBoundary(int subId, double pcoords[3], vtkIdList *pts); 00072 void Contour(double value, vtkDataArray *cellScalars, 00073 vtkIncrementalPointLocator *locator, vtkCellArray *verts, 00074 vtkCellArray *lines, vtkCellArray *polys, 00075 vtkPointData *inPd, vtkPointData *outPd, 00076 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd); 00077 int EvaluatePosition(double x[3], double* closestPoint, 00078 int& subId, double pcoords[3], 00079 double& dist2, double *weights); 00080 void EvaluateLocation(int& subId, double pcoords[3], double x[3], 00081 double *weights); 00082 int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts); 00083 void Derivatives(int subId, double pcoords[3], double *values, 00084 int dim, double *derivs); 00085 virtual double *GetParametricCoords(); 00086 00088 00091 void Clip(double value, vtkDataArray *cellScalars, 00092 vtkIncrementalPointLocator *locator, vtkCellArray *tets, 00093 vtkPointData *inPd, vtkPointData *outPd, 00094 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, 00095 int insideOut); 00097 00099 00101 int IntersectWithLine(double p1[3], double p2[3], double tol, double& t, 00102 double x[3], double pcoords[3], int& subId); 00104 00105 00107 int GetParametricCenter(double pcoords[3]); 00108 00111 static void InterpolationFunctions(double pcoords[3], double weights[13]); 00114 static void InterpolationDerivs(double pcoords[3], double derivs[39]); 00116 00118 virtual void InterpolateFunctions(double pcoords[3], double weights[13]) 00119 { 00120 vtkQuadraticPyramid::InterpolationFunctions(pcoords,weights); 00121 } 00122 virtual void InterpolateDerivs(double pcoords[3], double derivs[39]) 00123 { 00124 vtkQuadraticPyramid::InterpolationDerivs(pcoords,derivs); 00125 } 00127 00128 00130 static int *GetEdgeArray(int edgeId); 00131 static int *GetFaceArray(int faceId); 00133 00137 void JacobianInverse(double pcoords[3], double **inverse, double derivs[39]); 00138 00139 protected: 00140 vtkQuadraticPyramid(); 00141 ~vtkQuadraticPyramid(); 00142 00143 vtkQuadraticEdge *Edge; 00144 vtkQuadraticTriangle *TriangleFace; 00145 vtkQuadraticQuad *Face; 00146 vtkTetra *Tetra; 00147 vtkPyramid *Pyramid; 00148 vtkPointData *PointData; 00149 vtkCellData *CellData; 00150 vtkDoubleArray *CellScalars; 00151 vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping 00152 00153 void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, 00154 vtkDataArray *cellScalars); 00155 00156 private: 00157 vtkQuadraticPyramid(const vtkQuadraticPyramid&); // Not implemented. 00158 void operator=(const vtkQuadraticPyramid&); // Not implemented. 00159 }; 00160 //---------------------------------------------------------------------------- 00161 // Return the center of the quadratic pyramid in parametric coordinates. 00162 // 00163 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3]) 00164 { 00165 pcoords[0] = pcoords[1] = 6./13; 00166 pcoords[2] = 3./13; 00167 return 0; 00168 } 00169 00170 00171 #endif