VTK
|
00001 /*========================================================================= 00002 00003 Program: Visualization Toolkit 00004 Module: vtkTetra.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 =========================================================================*/ 00035 #ifndef __vtkTetra_h 00036 #define __vtkTetra_h 00037 00038 #include "vtkCell3D.h" 00039 00040 class vtkLine; 00041 class vtkTriangle; 00042 class vtkUnstructuredGrid; 00043 class vtkIncrementalPointLocator; 00044 00045 class VTK_FILTERING_EXPORT vtkTetra : public vtkCell3D 00046 { 00047 public: 00048 static vtkTetra *New(); 00049 vtkTypeMacro(vtkTetra,vtkCell3D); 00050 void PrintSelf(ostream& os, vtkIndent indent); 00051 00053 00054 virtual void GetEdgePoints(int edgeId, int* &pts); 00055 virtual void GetFacePoints(int faceId, int* &pts); 00057 00059 00060 int GetCellType() {return VTK_TETRA;} 00061 int GetNumberOfEdges() {return 6;} 00062 int GetNumberOfFaces() {return 4;} 00063 vtkCell *GetEdge(int edgeId); 00064 vtkCell *GetFace(int faceId); 00065 void Contour(double value, vtkDataArray *cellScalars, 00066 vtkIncrementalPointLocator *locator, vtkCellArray *verts, 00067 vtkCellArray *lines, vtkCellArray *polys, 00068 vtkPointData *inPd, vtkPointData *outPd, 00069 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd); 00070 void Clip(double value, vtkDataArray *cellScalars, 00071 vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, 00072 vtkPointData *inPd, vtkPointData *outPd, 00073 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, 00074 int insideOut); 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 00091 int CellBoundary(int subId, double pcoords[3], vtkIdList *pts); 00092 00094 int GetParametricCenter(double pcoords[3]); 00095 00098 double GetParametricDistance(double pcoords[3]); 00099 00101 00102 static void TetraCenter(double p1[3], double p2[3], double p3[3], double p4[3], 00103 double center[3]); 00105 00107 00110 static double Circumsphere(double p1[3], double p2[3], double p3[3], 00111 double p4[3], double center[3]); 00113 00115 00118 static double Insphere(double p1[3], double p2[3], double p3[3], 00119 double p4[3], double center[3]); 00121 00123 00134 static int BarycentricCoords(double x[3], double x1[3], double x2[3], 00135 double x3[3], double x4[3], double bcoords[4]); 00137 00139 00141 static double ComputeVolume(double p1[3], double p2[3], double p3[3], 00142 double p4[3]); 00144 00148 int JacobianInverse(double **inverse, double derivs[12]); 00149 00151 static void InterpolationFunctions(double pcoords[3], double weights[4]); 00153 static void InterpolationDerivs(double pcoords[3], double derivs[12]); 00155 00157 virtual void InterpolateFunctions(double pcoords[3], double weights[4]) 00158 { 00159 vtkTetra::InterpolationFunctions(pcoords,weights); 00160 } 00161 virtual void InterpolateDerivs(double pcoords[3], double derivs[12]) 00162 { 00163 vtkTetra::InterpolationDerivs(pcoords,derivs); 00164 } 00166 00168 00170 static int *GetEdgeArray(int edgeId); 00171 static int *GetFaceArray(int faceId); 00173 00174 protected: 00175 vtkTetra(); 00176 ~vtkTetra(); 00177 00178 vtkLine *Line; 00179 vtkTriangle *Triangle; 00180 00181 private: 00182 vtkTetra(const vtkTetra&); // Not implemented. 00183 void operator=(const vtkTetra&); // Not implemented. 00184 }; 00185 00186 inline int vtkTetra::GetParametricCenter(double pcoords[3]) 00187 { 00188 pcoords[0] = pcoords[1] = pcoords[2] = 0.25; 00189 return 0; 00190 } 00191 00192 #endif 00193 00194 00195