1 /*=========================================================================
2 
3   Program:   Visualization Toolkit
4   Module:    vtkQuadraticPyramid.h
5 
6   Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7   All rights reserved.
8   See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10      This software is distributed WITHOUT ANY WARRANTY; without even
11      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12      PURPOSE.  See the above copyright notice for more information.
13 
14 =========================================================================*/
15 /**
16  * @class   vtkQuadraticPyramid
17  * @brief   cell represents a parabolic, 13-node isoparametric pyramid
18  *
19  * vtkQuadraticPyramid is a concrete implementation of vtkNonLinearCell to
20  * represent a three-dimensional, 13-node isoparametric parabolic
21  * pyramid. The interpolation is the standard finite element, quadratic
22  * isoparametric shape function. The cell includes a mid-edge node. The
23  * ordering of the thirteen points defining the cell is point ids (0-4,5-12)
24  * where point ids 0-4 are the five corner vertices of the pyramid; followed
25  * by eight midedge nodes (5-12). Note that these midedge nodes lie
26  * on the edges defined by (0,1), (1,2), (2,3), (3,0), (0,4), (1,4), (2,4),
27  * (3,4), respectively. The parametric location of vertex #4 is [0, 0, 1].
28  *
29  * @sa
30  * vtkQuadraticEdge vtkQuadraticTriangle vtkQuadraticTetra
31  * vtkQuadraticHexahedron vtkQuadraticQuad vtkQuadraticWedge
32  *
33  * @par Thanks:
34  * The shape functions and derivatives could be implemented thanks to
35  * the report Pyramid Solid Elements Linear and Quadratic Iso-P Models
36  * From Center For Aerospace Structures
37 */
38 
39 #ifndef vtkQuadraticPyramid_h
40 #define vtkQuadraticPyramid_h
41 
42 #include "vtkCommonDataModelModule.h" // For export macro
43 #include "vtkNonLinearCell.h"
44 
45 class vtkQuadraticEdge;
46 class vtkQuadraticQuad;
47 class vtkQuadraticTriangle;
48 class vtkTetra;
49 class vtkPyramid;
50 class vtkDoubleArray;
51 
52 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
53 {
54 public:
55   static vtkQuadraticPyramid *New();
56   vtkTypeMacro(vtkQuadraticPyramid,vtkNonLinearCell);
57   void PrintSelf(ostream& os, vtkIndent indent) override;
58 
59   //@{
60   /**
61    * Implement the vtkCell API. See the vtkCell API for descriptions
62    * of these methods.
63    */
GetCellType()64   int GetCellType() override {return VTK_QUADRATIC_PYRAMID;};
GetCellDimension()65   int GetCellDimension() override {return 3;}
GetNumberOfEdges()66   int GetNumberOfEdges() override {return 8;}
GetNumberOfFaces()67   int GetNumberOfFaces() override {return 5;}
68   vtkCell *GetEdge(int edgeId) override;
69   vtkCell *GetFace(int faceId) override;
70   //@}
71 
72   int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
73   void Contour(double value, vtkDataArray *cellScalars,
74                vtkIncrementalPointLocator *locator, vtkCellArray *verts,
75                vtkCellArray *lines, vtkCellArray *polys,
76                vtkPointData *inPd, vtkPointData *outPd,
77                vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override;
78   int EvaluatePosition(const double x[3], double closestPoint[3],
79                        int& subId, double pcoords[3],
80                        double& dist2, double weights[]) override;
81   void EvaluateLocation(int& subId, const double pcoords[3], double x[3],
82                         double *weights) override;
83   int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
84   void Derivatives(int subId, const double pcoords[3], const double *values,
85                    int dim, double *derivs) override;
86   double *GetParametricCoords() override;
87 
88   /**
89    * Clip this quadratic triangle using scalar value provided. Like
90    * contouring, except that it cuts the triangle to produce linear
91    * triangles.
92    */
93   void Clip(double value, vtkDataArray *cellScalars,
94             vtkIncrementalPointLocator *locator, vtkCellArray *tets,
95             vtkPointData *inPd, vtkPointData *outPd,
96             vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
97             int insideOut) override;
98 
99   /**
100    * Line-edge intersection. Intersection has to occur within [0,1] parametric
101    * coordinates and with specified tolerance.
102    */
103   int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
104                         double x[3], double pcoords[3], int& subId) override;
105 
106 
107   /**
108    * Return the center of the quadratic pyramid in parametric coordinates.
109    */
110   int GetParametricCenter(double pcoords[3]) override;
111 
112   /**
113    * @deprecated Replaced by vtkQuadraticPyramid::InterpolateFunctions as of VTK 5.2
114    */
115   static void InterpolationFunctions(const double pcoords[3], double weights[13]);
116   /**
117    * @deprecated Replaced by vtkQuadraticPyramid::InterpolateDerivs as of VTK 5.2
118    */
119   static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
120   //@{
121   /**
122    * Compute the interpolation functions/derivatives
123    * (aka shape functions/derivatives)
124    */
InterpolateFunctions(const double pcoords[3],double weights[13])125   void InterpolateFunctions(const double pcoords[3], double weights[13]) override
126   {
127     vtkQuadraticPyramid::InterpolationFunctions(pcoords,weights);
128   }
InterpolateDerivs(const double pcoords[3],double derivs[39])129   void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
130   {
131     vtkQuadraticPyramid::InterpolationDerivs(pcoords,derivs);
132   }
133   //@}
134   //@{
135   /**
136    * Return the ids of the vertices defining edge/face (`edgeId`/`faceId').
137    * Ids are related to the cell, not to the dataset.
138    */
139   static int *GetEdgeArray(int edgeId);
140   static int *GetFaceArray(int faceId);
141   //@}
142 
143   /**
144    * Given parametric coordinates compute inverse Jacobian transformation
145    * matrix. Returns 9 elements of 3x3 inverse Jacobian plus interpolation
146    * function derivatives.
147    */
148   void JacobianInverse(const double pcoords[3], double **inverse, double derivs[39]);
149 
150 protected:
151   vtkQuadraticPyramid();
152   ~vtkQuadraticPyramid() override;
153 
154   vtkQuadraticEdge *Edge;
155   vtkQuadraticTriangle *TriangleFace;
156   vtkQuadraticQuad *Face;
157   vtkTetra         *Tetra;
158   vtkPyramid       *Pyramid;
159   vtkPointData     *PointData;
160   vtkCellData      *CellData;
161   vtkDoubleArray   *CellScalars;
162   vtkDoubleArray   *Scalars; //used to avoid New/Delete in contouring/clipping
163 
164   //@{
165   /**
166    * This method adds in a point at the center of the quadrilateral face
167    * and then interpolates values to that point. In order to do this it
168    * also resizes certain member variable arrays. For safety should call
169    * ResizeArrays() after the results of Subdivide() are not needed anymore.
170    **/
171   void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId,
172     vtkDataArray *cellScalars);
173   //@}
174   //@{
175   /**
176    * Resize the superclasses' member arrays to newSize where newSize should either be
177    * 13 or 14. Call with 13 to reset the reallocation done in the Subdivide()
178    * method or call with 14 to add one extra tuple for the generated point in
179    * Subdivice. For efficiency it only resizes the superclasses' arrays.
180    **/
181   void ResizeArrays(vtkIdType newSize);
182   //@}
183 
184 private:
185   vtkQuadraticPyramid(const vtkQuadraticPyramid&) = delete;
186   void operator=(const vtkQuadraticPyramid&) = delete;
187 };
188 //----------------------------------------------------------------------------
189 // Return the center of the quadratic pyramid in parametric coordinates.
190 //
GetParametricCenter(double pcoords[3])191 inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
192 {
193   pcoords[0] = pcoords[1] = 6.0/13.0;
194   pcoords[2] = 3.0/13.0;
195   return 0;
196 }
197 
198 
199 #endif
200