VTK
vtkLagrangeHexahedron.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkLagrangeHexahedron.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 =========================================================================*/
26 #ifndef vtkLagrangeHexahedron_h
27 #define vtkLagrangeHexahedron_h
28 
29 #include "vtkCommonDataModelModule.h" // For export macro
30 #include "vtkNonLinearCell.h"
31 #include "vtkSmartPointer.h" // For member variable.
32 #include "vtkCellType.h" // For GetCellType.
33 #include "vtkNew.h" // For member variable.
34 
35 class vtkCellData;
36 class vtkDoubleArray;
37 class vtkHexahedron;
38 class vtkIdList;
39 class vtkLagrangeCurve;
42 class vtkPointData;
43 class vtkPoints;
44 class vtkVector3d;
45 class vtkVector3i;
46 
47 class VTKCOMMONDATAMODEL_EXPORT vtkLagrangeHexahedron : public vtkNonLinearCell
48 {
49 public:
50  static vtkLagrangeHexahedron* New();
52  void PrintSelf(ostream& os, vtkIndent indent) override;
53 
54  int GetCellType() override { return VTK_LAGRANGE_HEXAHEDRON; }
55  int GetCellDimension() override { return 3; }
56  int RequiresInitialization() override { return 1; }
57  int GetNumberOfEdges() override { return 12; }
58  int GetNumberOfFaces() override { return 6; }
59  vtkCell* GetEdge(int edgeId) override;
60  vtkCell* GetFace(int faceId) override;
61 
62  void Initialize() override;
63 
64  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
65  int EvaluatePosition(const double x[3], double closestPoint[3],
66  int& subId, double pcoords[3],
67  double& dist2, double weights[]) override;
68  void EvaluateLocation(
69  int& subId, const double pcoords[3], double x[3],
70  double* weights) override;
71  void Contour(
72  double value, vtkDataArray* cellScalars,
74  vtkCellArray* lines, vtkCellArray* polys,
75  vtkPointData* inPd, vtkPointData* outPd,
76  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
77  void Clip(
78  double value, vtkDataArray* cellScalars,
80  vtkPointData* inPd, vtkPointData* outPd,
81  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd,
82  int insideOut) override;
83  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
84  double x[3], double pcoords[3], int& subId) override;
85  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
86  void Derivatives(
87  int subId, const double pcoords[3], const double* values,
88  int dim, double* derivs) override;
89  double* GetParametricCoords() override;
90  int GetParametricCenter(double center[3]) override;
91 
92  double GetParametricDistance(const double pcoords[3]) override;
93 
94  const int* GetOrder();
95  int GetOrder(int i) { return this->GetOrder()[i]; }
96 
97  void InterpolateFunctions(const double pcoords[3], double* weights) override;
98  void InterpolateDerivs(const double pcoords[3], double* derivs) override;
99 
100  bool SubCellCoordinatesFromId(vtkVector3i& ijk, int subId);
101  bool SubCellCoordinatesFromId(int& i, int& j, int& k, int subId);
102  static int PointIndexFromIJK(int i, int j, int k, const int* order);
103  int PointIndexFromIJK(int i, int j, int k);
104  bool TransformApproxToCellParams(int subCell, double* pcoords);
105  bool TransformFaceToCellParams(int bdyFace, double* pcoords);
106 
107 protected:
109  ~vtkLagrangeHexahedron() override;
110 
111  vtkHexahedron* GetApprox();
112  void PrepareApproxData(vtkPointData* pd, vtkCellData* cd, vtkIdType cellId, vtkDataArray* cellScalars);
113  vtkHexahedron* GetApproximateHex(
114  int subId, vtkDataArray* scalarsIn = nullptr, vtkDataArray* scalarsOut = nullptr);
115 
116  int Order[4];
128 
129 private:
131  void operator=(const vtkLagrangeHexahedron&) = delete;
132 };
133 
135 {
136  center[0] = center[1] = center[2] = 0.5;
137  return 0;
138 }
139 
140 #endif // vtkLagrangeHexahedron_h
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:33
vtkLagrangeHexahedron::GetNumberOfEdges
int GetNumberOfEdges() override
Return the number of edges in the cell.
Definition: vtkLagrangeHexahedron.h:57
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.
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.
vtkHexahedron
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:41
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:31
vtkX3D::value
Definition: vtkX3D.h:220
vtkIdType
int vtkIdType
Definition: vtkType.h:347
vtkLagrangeHexahedron::TmpIds
vtkNew< vtkIdList > TmpIds
Definition: vtkLagrangeHexahedron.h:124
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:111
vtkSmartPointer< vtkPoints >
VTK_LAGRANGE_HEXAHEDRON
Definition: vtkCellType.h:112
vtkLagrangeHexahedron::ApproxCD
vtkSmartPointer< vtkCellData > ApproxCD
Definition: vtkLagrangeHexahedron.h:120
vtkLagrangeHexahedron::FaceCell
vtkNew< vtkLagrangeQuadrilateral > FaceCell
Definition: vtkLagrangeHexahedron.h:125
vtkLagrangeHexahedron::Approx
vtkSmartPointer< vtkHexahedron > Approx
Definition: vtkLagrangeHexahedron.h:118
vtkLagrangeHexahedron::ApproxPD
vtkSmartPointer< vtkPointData > ApproxPD
Definition: vtkLagrangeHexahedron.h:119
vtkX3D::center
Definition: vtkX3D.h:230
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
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.
vtkLagrangeHexahedron::EdgeCell
vtkNew< vtkLagrangeCurve > EdgeCell
Definition: vtkLagrangeHexahedron.h:126
vtkLagrangeHexahedron::RequiresInitialization
int RequiresInitialization() override
Some cells require initialization prior to access.
Definition: vtkLagrangeHexahedron.h:56
vtkVector3i
Definition: vtkVector.h:439
vtkLagrangeHexahedron::GetOrder
int GetOrder(int i)
Definition: vtkLagrangeHexahedron.h:95
vtkLagrangeHexahedron::TmpPts
vtkNew< vtkPoints > TmpPts
Definition: vtkLagrangeHexahedron.h:123
vtkCell::Triangulate
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:56
vtkLagrangeHexahedron::CellScalars
vtkNew< vtkDoubleArray > CellScalars
Definition: vtkLagrangeHexahedron.h:121
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
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:33
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:44
vtkCellType.h
vtkLagrangeQuadrilateral
Definition: vtkLagrangeQuadrilateral.h:38
vtkSmartPointer.h
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:51
vtkLagrangeHexahedron::Interp
vtkNew< vtkLagrangeInterpolation > Interp
Definition: vtkLagrangeHexahedron.h:127
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:30
vtkNew< vtkDoubleArray >
vtkLagrangeInterpolation
Definition: vtkLagrangeInterpolation.h:34
vtkLagrangeHexahedron::PointParametricCoordinates
vtkSmartPointer< vtkPoints > PointParametricCoordinates
Definition: vtkLagrangeHexahedron.h:117
vtkLagrangeHexahedron::Scalars
vtkNew< vtkDoubleArray > Scalars
Definition: vtkLagrangeHexahedron.h:122
vtkX3D::order
Definition: vtkX3D.h:440
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...
vtkNonLinearCell.h
vtkCell::InterpolateDerivs
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:357
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.
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:354
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:35
vtkLagrangeHexahedron::GetParametricCenter
int GetParametricCenter(double center[3]) override
Return center of the cell in parametric coordinates.
Definition: vtkLagrangeHexahedron.h:134
vtkLagrangeHexahedron
A 3D cell that represents an arbitrary order Lagrange hex.
Definition: vtkLagrangeHexahedron.h:47
vtkX3D::index
Definition: vtkX3D.h:246
vtkLagrangeHexahedron::GetNumberOfFaces
int GetNumberOfFaces() override
Return the number of faces in the cell.
Definition: vtkLagrangeHexahedron.h:58
vtkVector3d
Definition: vtkVector.h:459
vtkLagrangeHexahedron::GetCellType
int GetCellType() override
Return the type of cell.
Definition: vtkLagrangeHexahedron.h:54
vtkLagrangeHexahedron::GetCellDimension
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
Definition: vtkLagrangeHexahedron.h:55