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 "vtkCommonDataModelModule.h" // For export macro 00045 #include "vtkNonLinearCell.h" 00046 00047 class vtkQuadraticEdge; 00048 class vtkQuadraticQuad; 00049 class vtkQuadraticTriangle; 00050 class vtkTetra; 00051 class vtkPyramid; 00052 class vtkDoubleArray; 00053 00054 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell 00055 { 00056 public: 00057 static vtkQuadraticPyramid *New(); 00058 vtkTypeMacro(vtkQuadraticPyramid,vtkNonLinearCell); 00059 void PrintSelf(ostream& os, vtkIndent indent); 00060 00062 00064 int GetCellType() {return VTK_QUADRATIC_PYRAMID;}; 00065 int GetCellDimension() {return 3;} 00066 int GetNumberOfEdges() {return 8;} 00067 int GetNumberOfFaces() {return 5;} 00068 vtkCell *GetEdge(int edgeId); 00069 vtkCell *GetFace(int faceId); 00071 00072 int CellBoundary(int subId, double pcoords[3], vtkIdList *pts); 00073 void Contour(double value, vtkDataArray *cellScalars, 00074 vtkIncrementalPointLocator *locator, vtkCellArray *verts, 00075 vtkCellArray *lines, vtkCellArray *polys, 00076 vtkPointData *inPd, vtkPointData *outPd, 00077 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd); 00078 int EvaluatePosition(double x[3], double* closestPoint, 00079 int& subId, double pcoords[3], 00080 double& dist2, double *weights); 00081 void EvaluateLocation(int& subId, double pcoords[3], double x[3], 00082 double *weights); 00083 int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts); 00084 void Derivatives(int subId, double pcoords[3], double *values, 00085 int dim, double *derivs); 00086 virtual double *GetParametricCoords(); 00087 00089 00092 void Clip(double value, vtkDataArray *cellScalars, 00093 vtkIncrementalPointLocator *locator, vtkCellArray *tets, 00094 vtkPointData *inPd, vtkPointData *outPd, 00095 vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, 00096 int insideOut); 00098 00100 00102 int IntersectWithLine(double p1[3], double p2[3], double tol, double& t, 00103 double x[3], double pcoords[3], int& subId); 00105 00106 00108 int GetParametricCenter(double pcoords[3]); 00109 00112 static void InterpolationFunctions(double pcoords[3], double weights[13]); 00115 static void InterpolationDerivs(double pcoords[3], double derivs[39]); 00117 00119 virtual void InterpolateFunctions(double pcoords[3], double weights[13]) 00120 { 00121 vtkQuadraticPyramid::InterpolationFunctions(pcoords,weights); 00122 } 00123 virtual void InterpolateDerivs(double pcoords[3], double derivs[39]) 00124 { 00125 vtkQuadraticPyramid::InterpolationDerivs(pcoords,derivs); 00126 } 00128 00129 00131 static int *GetEdgeArray(int edgeId); 00132 static int *GetFaceArray(int faceId); 00134 00138 void JacobianInverse(double pcoords[3], double **inverse, double derivs[39]); 00139 00140 protected: 00141 vtkQuadraticPyramid(); 00142 ~vtkQuadraticPyramid(); 00143 00144 vtkQuadraticEdge *Edge; 00145 vtkQuadraticTriangle *TriangleFace; 00146 vtkQuadraticQuad *Face; 00147 vtkTetra *Tetra; 00148 vtkPyramid *Pyramid; 00149 vtkPointData *PointData; 00150 vtkCellData *CellData; 00151 vtkDoubleArray *CellScalars; 00152 vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping 00153 00154 void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId, 00155 vtkDataArray *cellScalars); 00156 00157 private: 00158 vtkQuadraticPyramid(const vtkQuadraticPyramid&); // Not implemented. 00159 void operator=(const vtkQuadraticPyramid&); // Not implemented. 00160 }; 00161 //---------------------------------------------------------------------------- 00162 // Return the center of the quadratic pyramid in parametric coordinates. 00163 // 00164 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3]) 00165 { 00166 pcoords[0] = pcoords[1] = 6./13; 00167 pcoords[2] = 3./13; 00168 return 0; 00169 } 00170 00171 00172 #endif