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