VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: vtkPyramid.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 =========================================================================*/ 00036 #ifndef __vtkPyramid_h 00037 #define __vtkPyramid_h 00038 00039 #include "vtkCell3D.h" 00040 00041 class vtkLine; 00042 class vtkQuad; 00043 class vtkTriangle; 00044 class vtkUnstructuredGrid; 00045 class vtkIncrementalPointLocator; 00046 00047 class VTK_FILTERING_EXPORT vtkPyramid : public vtkCell3D 00048 { 00049 public: 00050 static vtkPyramid *New(); 00051 vtkTypeMacro(vtkPyramid,vtkCell3D); 00052 void PrintSelf(ostream& os, vtkIndent indent); 00053 00055 00056 virtual void GetEdgePoints(int edgeId, int* &pts); 00057 virtual void GetFacePoints(int faceId, int* &pts); 00059 00061 00062 int GetCellType() {return VTK_PYRAMID;} 00063 int GetCellDimension() {return 3;} 00064 int GetNumberOfEdges() {return 8;} 00065 int GetNumberOfFaces() {return 5;} 00066 vtkCell *GetEdge(int edgeId); 00067 vtkCell *GetFace(int faceId); 00068 int CellBoundary(int subId, double pcoords[3], vtkIdList *pts); 00069 void Contour(double value, vtkDataArray *cellScalars, 00070 vtkIncrementalPointLocator *locator, vtkCellArray *verts, 00071 vtkCellArray *lines, vtkCellArray *polys, 00072 vtkPointData *inPd, vtkPointData *outPd, 00073 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd); 00074 int EvaluatePosition(double x[3], double* closestPoint, 00075 int& subId, double pcoords[3], 00076 double& dist2, double *weights); 00077 void EvaluateLocation(int& subId, double pcoords[3], double x[3], 00078 double *weights); 00079 int IntersectWithLine(double p1[3], double p2[3], double tol, double& t, 00080 double x[3], double pcoords[3], int& subId); 00081 int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts); 00082 void Derivatives(int subId, double pcoords[3], double *values, 00083 int dim, double *derivs); 00084 virtual double *GetParametricCoords(); 00086 00088 int GetParametricCenter(double pcoords[3]); 00089 00092 static void InterpolationFunctions(double pcoords[3], double weights[5]); 00094 static void InterpolationDerivs(double pcoords[3], double derivs[15]); 00096 00098 virtual void InterpolateFunctions(double pcoords[3], double weights[5]) 00099 { 00100 vtkPyramid::InterpolationFunctions(pcoords,weights); 00101 } 00102 virtual void InterpolateDerivs(double pcoords[3], double derivs[15]) 00103 { 00104 vtkPyramid::InterpolationDerivs(pcoords,derivs); 00105 } 00107 00108 int JacobianInverse(double pcoords[3], double **inverse, double derivs[15]); 00109 00111 00113 static int *GetEdgeArray(int edgeId); 00114 static int *GetFaceArray(int faceId); 00116 00117 protected: 00118 vtkPyramid(); 00119 ~vtkPyramid(); 00120 00121 vtkLine *Line; 00122 vtkTriangle *Triangle; 00123 vtkQuad *Quad; 00124 00125 private: 00126 vtkPyramid(const vtkPyramid&); // Not implemented. 00127 void operator=(const vtkPyramid&); // Not implemented. 00128 }; 00129 00130 //---------------------------------------------------------------------------- 00131 inline int vtkPyramid::GetParametricCenter(double pcoords[3]) 00132 { 00133 pcoords[0] = pcoords[1] = 0.4; 00134 pcoords[2] = 0.2; 00135 return 0; 00136 } 00137 00138 #endif 00139 00140 00141