File: vtkQuadraticPyramid.h

package info (click to toggle)
vtk9 9.0.1%2Bdfsg1-8
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 133,688 kB
  • sloc: cpp: 1,568,287; ansic: 208,587; python: 87,847; xml: 8,022; java: 4,509; yacc: 4,027; sh: 2,515; perl: 2,183; lex: 1,766; objc: 143; makefile: 126; tcl: 59
file content (188 lines) | stat: -rw-r--r-- 7,195 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
/*=========================================================================

  Program:   Visualization Toolkit
  Module:    vtkQuadraticPyramid.h

  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
  All rights reserved.
  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.

     This software is distributed WITHOUT ANY WARRANTY; without even
     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
     PURPOSE.  See the above copyright notice for more information.

=========================================================================*/
/**
 * @class   vtkQuadraticPyramid
 * @brief   cell represents a parabolic, 13-node isoparametric pyramid
 *
 * vtkQuadraticPyramid is a concrete implementation of vtkNonLinearCell to
 * represent a three-dimensional, 13-node isoparametric parabolic
 * pyramid. The interpolation is the standard finite element, quadratic
 * isoparametric shape function. The cell includes a mid-edge node. The
 * ordering of the thirteen points defining the cell is point ids (0-4,5-12)
 * where point ids 0-4 are the five corner vertices of the pyramid; followed
 * by eight midedge nodes (5-12). Note that these midedge nodes lie
 * on the edges defined by (0,1), (1,2), (2,3), (3,0), (0,4), (1,4), (2,4),
 * (3,4), respectively. The parametric location of vertex #4 is [0, 0, 1].
 *
 * @sa
 * vtkQuadraticEdge vtkQuadraticTriangle vtkQuadraticTetra
 * vtkQuadraticHexahedron vtkQuadraticQuad vtkQuadraticWedge
 *
 * @par Thanks:
 * The shape functions and derivatives could be implemented thanks to
 * the report Pyramid Solid Elements Linear and Quadratic Iso-P Models
 * From Center For Aerospace Structures
 */

#ifndef vtkQuadraticPyramid_h
#define vtkQuadraticPyramid_h

#include "vtkCommonDataModelModule.h" // For export macro
#include "vtkNonLinearCell.h"

class vtkQuadraticEdge;
class vtkQuadraticQuad;
class vtkQuadraticTriangle;
class vtkTetra;
class vtkPyramid;
class vtkDoubleArray;

class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticPyramid : public vtkNonLinearCell
{
public:
  static vtkQuadraticPyramid* New();
  vtkTypeMacro(vtkQuadraticPyramid, vtkNonLinearCell);
  void PrintSelf(ostream& os, vtkIndent indent) override;

  //@{
  /**
   * Implement the vtkCell API. See the vtkCell API for descriptions
   * of these methods.
   */
  int GetCellType() override { return VTK_QUADRATIC_PYRAMID; }
  int GetCellDimension() override { return 3; }
  int GetNumberOfEdges() override { return 8; }
  int GetNumberOfFaces() override { return 5; }
  vtkCell* GetEdge(int edgeId) override;
  vtkCell* GetFace(int faceId) override;
  //@}

  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
  void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
    vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
    vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
  int EvaluatePosition(const double x[3], double closestPoint[3], int& subId, double pcoords[3],
    double& dist2, double weights[]) override;
  void EvaluateLocation(int& subId, const double pcoords[3], double x[3], double* weights) override;
  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
  void Derivatives(
    int subId, const double pcoords[3], const double* values, int dim, double* derivs) override;
  double* GetParametricCoords() override;

  /**
   * Clip this quadratic triangle using scalar value provided. Like
   * contouring, except that it cuts the triangle to produce linear
   * triangles.
   */
  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
    vtkCellArray* tets, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
    vtkIdType cellId, vtkCellData* outCd, int insideOut) override;

  /**
   * Line-edge intersection. Intersection has to occur within [0,1] parametric
   * coordinates and with specified tolerance.
   */
  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
    double pcoords[3], int& subId) override;

  /**
   * Return the center of the quadratic pyramid in parametric coordinates.
   */
  int GetParametricCenter(double pcoords[3]) override;

  static void InterpolationFunctions(const double pcoords[3], double weights[13]);
  static void InterpolationDerivs(const double pcoords[3], double derivs[39]);
  //@{
  /**
   * Compute the interpolation functions/derivatives
   * (aka shape functions/derivatives)
   */
  void InterpolateFunctions(const double pcoords[3], double weights[13]) override
  {
    vtkQuadraticPyramid::InterpolationFunctions(pcoords, weights);
  }
  void InterpolateDerivs(const double pcoords[3], double derivs[39]) override
  {
    vtkQuadraticPyramid::InterpolationDerivs(pcoords, derivs);
  }
  //@}
  //@{
  /**
   * Return the ids of the vertices defining edge/face (`edgeId`/`faceId').
   * Ids are related to the cell, not to the dataset.
   *
   * @note The return type changed. It used to be int*, it is now const vtkIdType*.
   * This is so ids are unified between vtkCell and vtkPoints.
   */
  static const vtkIdType* GetEdgeArray(vtkIdType edgeId);
  static const vtkIdType* GetFaceArray(vtkIdType faceId);
  //@}

  /**
   * Given parametric coordinates compute inverse Jacobian transformation
   * matrix. Returns 9 elements of 3x3 inverse Jacobian plus interpolation
   * function derivatives.
   */
  void JacobianInverse(const double pcoords[3], double** inverse, double derivs[39]);

protected:
  vtkQuadraticPyramid();
  ~vtkQuadraticPyramid() override;

  vtkQuadraticEdge* Edge;
  vtkQuadraticTriangle* TriangleFace;
  vtkQuadraticQuad* Face;
  vtkTetra* Tetra;
  vtkPyramid* Pyramid;
  vtkPointData* PointData;
  vtkCellData* CellData;
  vtkDoubleArray* CellScalars;
  vtkDoubleArray* Scalars; // used to avoid New/Delete in contouring/clipping

  //@{
  /**
   * This method adds in a point at the center of the quadrilateral face
   * and then interpolates values to that point. In order to do this it
   * also resizes certain member variable arrays. For safety should call
   * ResizeArrays() after the results of Subdivide() are not needed anymore.
   **/
  void Subdivide(
    vtkPointData* inPd, vtkCellData* inCd, vtkIdType cellId, vtkDataArray* cellScalars);
  //@}
  //@{
  /**
   * Resize the superclasses' member arrays to newSize where newSize should either be
   * 13 or 14. Call with 13 to reset the reallocation done in the Subdivide()
   * method or call with 14 to add one extra tuple for the generated point in
   * Subdivice. For efficiency it only resizes the superclasses' arrays.
   **/
  void ResizeArrays(vtkIdType newSize);
  //@}

private:
  vtkQuadraticPyramid(const vtkQuadraticPyramid&) = delete;
  void operator=(const vtkQuadraticPyramid&) = delete;
};
//----------------------------------------------------------------------------
// Return the center of the quadratic pyramid in parametric coordinates.
//
inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
{
  pcoords[0] = pcoords[1] = 6.0 / 13.0;
  pcoords[2] = 3.0 / 13.0;
  return 0;
}

#endif