File: vtkHigherOrderWedge.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 (164 lines) | stat: -rw-r--r-- 7,131 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
/*=========================================================================

  Program:   Visualization Toolkit
  Module:    vtkHigherOrderWedge.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   vtkHigherOrderWedge
 * @brief   A 3D cell that represents an arbitrary order HigherOrder wedge
 *
 * vtkHigherOrderWedge is a concrete implementation of vtkCell to represent a
 * 3D wedge using HigherOrder shape functions of user specified order.
 * A wedge consists of two triangular and three quadrilateral faces.
 * The first six points of the wedge (0-5) are the "corner" points
 * where the first three points are the base of the wedge. This wedge
 * point ordering is opposite the vtkWedge ordering though in that
 * the base of the wedge defined by the first three points (0,1,2) form
 * a triangle whose normal points inward (toward the triangular face (3,4,5)).
 * While this is opposite the vtkWedge convention it is consistent with
 * every other cell type in VTK. The first 2 parametric coordinates of the
 * HigherOrder wedge or for the triangular base and vary between 0 and 1. The
 * third parametric coordinate is between the two triangular faces and goes
 * from 0 to 1 as well.
 */

#ifndef vtkHigherOrderWedge_h
#define vtkHigherOrderWedge_h

#include <functional> //For std::function

#include "vtkCellType.h"              // For GetCellType.
#include "vtkCommonDataModelModule.h" // For export macro
#include "vtkNew.h"                   // For member variable.
#include "vtkNonLinearCell.h"
#include "vtkSmartPointer.h" // For member variable.

class vtkCellData;
class vtkDoubleArray;
class vtkWedge;
class vtkIdList;
class vtkPointData;
class vtkPoints;
class vtkVector3d;
class vtkVector3i;
class vtkHigherOrderCurve;
class vtkHigherOrderInterpolation;
class vtkHigherOrderQuadrilateral;
class vtkHigherOrderTriangle;

class VTKCOMMONDATAMODEL_EXPORT vtkHigherOrderWedge : public vtkNonLinearCell
{
public:
  vtkTypeMacro(vtkHigherOrderWedge, vtkNonLinearCell);
  void PrintSelf(ostream& os, vtkIndent indent) override;

  int GetCellType() override = 0;
  int GetCellDimension() override { return 3; }
  int RequiresInitialization() override { return 1; }
  int GetNumberOfEdges() override { return 9; }
  int GetNumberOfFaces() override { return 5; }
  vtkCell* GetEdge(int edgeId) override = 0;
  void SetEdgeIdsAndPoints(int edgeId,
    const std::function<void(const vtkIdType&)>& set_number_of_ids_and_points,
    const std::function<void(const vtkIdType&, const vtkIdType&)>& set_ids_and_points);
  vtkCell* GetFace(int faceId) override = 0;

  void Initialize() override;

  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) 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;
  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;
  void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
    vtkCellArray* polys, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
    vtkIdType cellId, vtkCellData* outCd, int insideOut) override;
  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
    double pcoords[3], int& subId) 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;
  void SetParametricCoords();
  double* GetParametricCoords() override;
  int GetParametricCenter(double center[3]) override;

  double GetParametricDistance(const double pcoords[3]) override;

  virtual void SetOrderFromCellData(
    vtkCellData* cell_data, const vtkIdType numPts, const vtkIdType cell_id);
  virtual void SetUniformOrderFromNumPoints(const vtkIdType numPts);
  virtual void SetOrder(const int s, const int t, const int u, const vtkIdType numPts);
  virtual const int* GetOrder();
  virtual int GetOrder(int i) { return this->GetOrder()[i]; }

  void InterpolateFunctions(const double pcoords[3], double* weights) override = 0;
  void InterpolateDerivs(const double pcoords[3], double* derivs) override = 0;

  bool SubCellCoordinatesFromId(vtkVector3i& ijk, int subId);
  bool SubCellCoordinatesFromId(int& i, int& j, int& k, int subId);
  static int PointIndexFromIJK(int i, int j, int k, const int* order);
  int PointIndexFromIJK(int i, int j, int k);
  bool TransformApproxToCellParams(int subCell, double* pcoords);
  bool TransformFaceToCellParams(int bdyFace, double* pcoords);

  static int GetNumberOfApproximatingWedges(const int* order);
  int GetNumberOfApproximatingWedges()
  {
    return vtkHigherOrderWedge::GetNumberOfApproximatingWedges(this->GetOrder());
  }
  virtual vtkHigherOrderQuadrilateral* getBdyQuad() = 0;
  virtual vtkHigherOrderTriangle* getBdyTri() = 0;
  virtual vtkHigherOrderCurve* getEdgeCell() = 0;
  virtual vtkHigherOrderInterpolation* getInterp() = 0;

protected:
  vtkHigherOrderWedge();
  ~vtkHigherOrderWedge() override;

  vtkWedge* GetApprox();
  void PrepareApproxData(
    vtkPointData* pd, vtkCellData* cd, vtkIdType cellId, vtkDataArray* cellScalars);
  vtkWedge* GetApproximateWedge(
    int subId, vtkDataArray* scalarsIn = nullptr, vtkDataArray* scalarsOut = nullptr);

  void GetTriangularFace(vtkHigherOrderTriangle* result, int faceId,
    const std::function<void(const vtkIdType&)>& set_number_of_ids_and_points,
    const std::function<void(const vtkIdType&, const vtkIdType&)>& set_ids_and_points);
  void GetQuadrilateralFace(vtkHigherOrderQuadrilateral* result, int faceId,
    const std::function<void(const vtkIdType&)>& set_number_of_ids_and_points,
    const std::function<void(const vtkIdType&, const vtkIdType&)>& set_ids_and_points);

  int Order[4];
  vtkSmartPointer<vtkPoints> PointParametricCoordinates;
  vtkSmartPointer<vtkWedge> Approx;
  vtkSmartPointer<vtkPointData> ApproxPD;
  vtkSmartPointer<vtkCellData> ApproxCD;
  vtkNew<vtkDoubleArray> CellScalars;
  vtkNew<vtkDoubleArray> Scalars;
  vtkNew<vtkPoints> TmpPts;
  vtkNew<vtkIdList> TmpIds;

private:
  vtkHigherOrderWedge(const vtkHigherOrderWedge&) = delete;
  void operator=(const vtkHigherOrderWedge&) = delete;
};

inline int vtkHigherOrderWedge::GetParametricCenter(double center[3])
{
  center[0] = center[1] = 1. / 3.;
  center[2] = 0.5;
  return 0;
}

#endif // vtkHigherOrderWedge_h