VTK
vtkLagrangeWedge.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkLagrangeWedge.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 =========================================================================*/
34 #ifndef vtkLagrangeWedge_h
35 #define vtkLagrangeWedge_h
36 
37 #include "vtkCommonDataModelModule.h" // For export macro
38 #include "vtkNonLinearCell.h"
39 #include "vtkSmartPointer.h" // For member variable.
40 #include "vtkCellType.h" // For GetCellType.
41 #include "vtkNew.h" // For member variable.
42 
43 class vtkCellData;
44 class vtkDoubleArray;
45 class vtkWedge;
46 class vtkIdList;
47 class vtkPointData;
48 class vtkPoints;
49 class vtkVector3d;
50 class vtkVector3i;
51 class vtkLagrangeCurve;
55 
56 class VTKCOMMONDATAMODEL_EXPORT vtkLagrangeWedge : public vtkNonLinearCell
57 {
58 public:
59  static vtkLagrangeWedge* New();
61  void PrintSelf(ostream& os, vtkIndent indent) override;
62 
63  int GetCellType() override { return VTK_LAGRANGE_WEDGE; }
64  int GetCellDimension() override { return 3; }
65  int RequiresInitialization() override { return 1; }
66  int GetNumberOfEdges() override { return 9; }
67  int GetNumberOfFaces() override { return 5; }
68  vtkCell* GetEdge(int edgeId) override;
69  vtkCell* GetFace(int faceId) override;
70 
71  void Initialize() override;
72 
73  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
74  int EvaluatePosition(const double x[3], double closestPoint[3],
75  int& subId, double pcoords[3],
76  double& dist2, double weights[]) override;
77  void EvaluateLocation(
78  int& subId, const double pcoords[3], double x[3],
79  double* weights) override;
80  void Contour(
81  double value, vtkDataArray* cellScalars,
83  vtkCellArray* lines, vtkCellArray* polys,
84  vtkPointData* inPd, vtkPointData* outPd,
85  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
86  void Clip(
87  double value, vtkDataArray* cellScalars,
89  vtkPointData* inPd, vtkPointData* outPd,
90  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd,
91  int insideOut) override;
92  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
93  double x[3], double pcoords[3], int& subId) override;
94  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
95  void Derivatives(
96  int subId, const double pcoords[3], const double* values,
97  int dim, double* derivs) override;
98  double* GetParametricCoords() override;
99  int GetParametricCenter(double center[3]) override;
100 
101  double GetParametricDistance(const double pcoords[3]) override;
102 
103  const int* GetOrder();
104  int GetOrder(int i) { return this->GetOrder()[i]; }
105 
106  void InterpolateFunctions(const double pcoords[3], double* weights) override;
107  void InterpolateDerivs(const double pcoords[3], double* derivs) override;
108 
109  bool SubCellCoordinatesFromId(vtkVector3i& ijk, int subId);
110  bool SubCellCoordinatesFromId(int& i, int& j, int& k, int subId);
111  static int PointIndexFromIJK(int i, int j, int k, const int* order);
112  int PointIndexFromIJK(int i, int j, int k);
113  bool TransformApproxToCellParams(int subCell, double* pcoords);
114  bool TransformFaceToCellParams(int bdyFace, double* pcoords);
115 
116  static int GetNumberOfApproximatingWedges(const int* order);
118  { return vtkLagrangeWedge::GetNumberOfApproximatingWedges(this->GetOrder()); }
119 
120 protected:
122  ~vtkLagrangeWedge() override;
123 
124  vtkWedge* GetApprox();
125  void PrepareApproxData(vtkPointData* pd, vtkCellData* cd, vtkIdType cellId, vtkDataArray* cellScalars);
126  vtkWedge* GetApproximateWedge(
127  int subId, vtkDataArray* scalarsIn = nullptr, vtkDataArray* scalarsOut = nullptr);
128 
129  vtkLagrangeTriangle* GetTriangularFace(int iAxis, int k);
130  vtkLagrangeQuadrilateral* GetQuadrilateralFace(int di, int dj);
131 
132  int Order[4];
145 
146 private:
147  vtkLagrangeWedge(const vtkLagrangeWedge&) = delete;
148  void operator=(const vtkLagrangeWedge&) = delete;
149 };
150 
152 {
153  center[0] = center[1] = 1./3.;
154  center[2] = 0.5;
155  return 0;
156 }
157 
158 #endif // vtkLagrangeWedge_h
int GetParametricCenter(double center[3]) override
Return center of the cell in parametric coordinates.
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkNew< vtkLagrangeCurve > BdyEdge
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
vtkSmartPointer< vtkWedge > Approx
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this leve...
Definition: vtkCell.h:354
int GetCellType() override
Return the type of cell.
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:347
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetNumberOfApproximatingWedges()
virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:357
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
dynamic, self-adjusting array of double
vtkNew< vtkDoubleArray > Scalars
virtual double GetParametricDistance(const double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
abstract class to specify cell behavior
Definition: vtkCell.h:56
virtual void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkNew< vtkDoubleArray > CellScalars
A 3D cell that represents an arbitrary order Lagrange wedge.
A 2D cell that represents an arbitrary order Lagrange triangle.
a simple class to control print indentation
Definition: vtkIndent.h:33
list of point or cell ids
Definition: vtkIdList.h:30
vtkNew< vtkLagrangeTriangle > BdyTri
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
Cut (or clip) the cell based on the input cellScalars and the specified value.
vtkNew< vtkLagrangeInterpolation > Interp
vtkSmartPointer< vtkCellData > ApproxCD
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
virtual int EvaluatePosition(const double x[3], double closestPoint[3], int &subId, double pcoords[3], double &dist2, double weights[])=0
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
object to represent cell connectivity
Definition: vtkCellArray.h:44
int RequiresInitialization() override
Some cells require initialization prior to access.
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
virtual void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd)=0
Generate contouring primitives.
int GetNumberOfFaces() override
Return the number of faces in the cell.
virtual void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
vtkNew< vtkIdList > TmpIds
vtkSmartPointer< vtkPointData > ApproxPD
virtual void Initialize()
Definition: vtkCell.h:111
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
vtkSmartPointer< vtkPoints > PointParametricCoordinates
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
vtkNew< vtkPoints > TmpPts
int GetNumberOfEdges() override
Return the number of edges in the cell.
vtkNew< vtkLagrangeQuadrilateral > BdyQuad
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:43
represent and manipulate 3D points
Definition: vtkPoints.h:33