VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: vtkWedge.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 __vtkWedge_h 00037 #define __vtkWedge_h 00038 00039 #include "vtkCommonDataModelModule.h" // For export macro 00040 #include "vtkCell3D.h" 00041 00042 class vtkLine; 00043 class vtkTriangle; 00044 class vtkQuad; 00045 class vtkUnstructuredGrid; 00046 class vtkIncrementalPointLocator; 00047 00048 class VTKCOMMONDATAMODEL_EXPORT vtkWedge : public vtkCell3D 00049 { 00050 public: 00051 static vtkWedge *New(); 00052 vtkTypeMacro(vtkWedge,vtkCell3D); 00053 void PrintSelf(ostream& os, vtkIndent indent); 00054 00056 00057 virtual void GetEdgePoints(int edgeId, int* &pts); 00058 virtual void GetFacePoints(int faceId, int* &pts); 00060 00062 00063 int GetCellType() {return VTK_WEDGE;} 00064 int GetCellDimension() {return 3;} 00065 int GetNumberOfEdges() {return 9;} 00066 int GetNumberOfFaces() {return 5;} 00067 vtkCell *GetEdge(int edgeId); 00068 vtkCell *GetFace(int faceId); 00069 int CellBoundary(int subId, double pcoords[3], vtkIdList *pts); 00070 void Contour(double value, vtkDataArray *cellScalars, 00071 vtkIncrementalPointLocator *locator, vtkCellArray *verts, 00072 vtkCellArray *lines, vtkCellArray *polys, 00073 vtkPointData *inPd, vtkPointData *outPd, 00074 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd); 00075 int EvaluatePosition(double x[3], double* closestPoint, 00076 int& subId, double pcoords[3], 00077 double& dist2, double *weights); 00078 void EvaluateLocation(int& subId, double pcoords[3], double x[3], 00079 double *weights); 00080 int IntersectWithLine(double p1[3], double p2[3], double tol, double& t, 00081 double x[3], double pcoords[3], int& subId); 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(); 00087 00089 int GetParametricCenter(double pcoords[3]); 00090 00092 static void InterpolationFunctions(double pcoords[3], double weights[6]); 00094 static void InterpolationDerivs(double pcoords[3], double derivs[18]); 00096 00098 virtual void InterpolateFunctions(double pcoords[3], double weights[6]) 00099 { 00100 vtkWedge::InterpolationFunctions(pcoords,weights); 00101 } 00102 virtual void InterpolateDerivs(double pcoords[3], double derivs[18]) 00103 { 00104 vtkWedge::InterpolationDerivs(pcoords,derivs); 00105 } 00106 int JacobianInverse(double pcoords[3], double **inverse, double derivs[18]); 00108 00110 00112 static int *GetEdgeArray(int edgeId); 00113 static int *GetFaceArray(int faceId); 00115 00116 protected: 00117 vtkWedge(); 00118 ~vtkWedge(); 00119 00120 vtkLine *Line; 00121 vtkTriangle *Triangle; 00122 vtkQuad *Quad; 00123 00124 private: 00125 vtkWedge(const vtkWedge&); // Not implemented. 00126 void operator=(const vtkWedge&); // Not implemented. 00127 }; 00128 00129 inline int vtkWedge::GetParametricCenter(double pcoords[3]) 00130 { 00131 pcoords[0] = pcoords[1] = 0.333333; 00132 pcoords[2] = 0.5; 00133 return 0; 00134 } 00135 00136 #endif 00137 00138 00139