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 =========================================================================*/
29 #ifndef vtkLagrangeHexahedron_h
30 #define vtkLagrangeHexahedron_h
31 
32 #include "vtkCommonDataModelModule.h" // For export macro
33 #include "vtkNonLinearCell.h"
34 #include "vtkSmartPointer.h" // For member variable.
35 #include "vtkCellType.h" // For GetCellType.
36 #include "vtkNew.h" // For member variable.
37 
38 class vtkCellData;
39 class vtkDoubleArray;
40 class vtkHexahedron;
41 class vtkIdList;
42 class vtkLagrangeCurve;
45 class vtkPointData;
46 class vtkPoints;
47 class vtkVector3d;
48 class vtkVector3i;
49 
50 class VTKCOMMONDATAMODEL_EXPORT vtkLagrangeHexahedron : public vtkNonLinearCell
51 {
52 public:
53  static vtkLagrangeHexahedron* New();
55  void PrintSelf(ostream& os, vtkIndent indent) override;
56 
57  int GetCellType() override { return VTK_LAGRANGE_HEXAHEDRON; }
58  int GetCellDimension() override { return 3; }
59  int RequiresInitialization() override { return 1; }
60  int GetNumberOfEdges() override { return 12; }
61  int GetNumberOfFaces() override { return 6; }
62  vtkCell* GetEdge(int edgeId) override;
63  vtkCell* GetFace(int faceId) override;
64 
65  void Initialize() override;
66 
67  int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) override;
68  int EvaluatePosition(const double x[3], double closestPoint[3],
69  int& subId, double pcoords[3],
70  double& dist2, double weights[]) override;
71  void EvaluateLocation(
72  int& subId, const double pcoords[3], double x[3],
73  double* weights) override;
74  void Contour(
75  double value, vtkDataArray* cellScalars,
77  vtkCellArray* lines, vtkCellArray* polys,
78  vtkPointData* inPd, vtkPointData* outPd,
79  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) override;
80  void Clip(
81  double value, vtkDataArray* cellScalars,
83  vtkPointData* inPd, vtkPointData* outPd,
84  vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd,
85  int insideOut) override;
86  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
87  double x[3], double pcoords[3], int& subId) override;
88  int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) override;
89  void Derivatives(
90  int subId, const double pcoords[3], const double* values,
91  int dim, double* derivs) override;
92  double* GetParametricCoords() override;
93  int GetParametricCenter(double center[3]) override;
94 
95  double GetParametricDistance(const double pcoords[3]) override;
96 
97  const int* GetOrder();
98  int GetOrder(int i) { return this->GetOrder()[i]; }
99 
100  void InterpolateFunctions(const double pcoords[3], double* weights) override;
101  void InterpolateDerivs(const double pcoords[3], double* derivs) override;
102 
103  bool SubCellCoordinatesFromId(vtkVector3i& ijk, int subId);
104  bool SubCellCoordinatesFromId(int& i, int& j, int& k, int subId);
105  static int PointIndexFromIJK(int i, int j, int k, const int* order);
106  int PointIndexFromIJK(int i, int j, int k);
107  bool TransformApproxToCellParams(int subCell, double* pcoords);
108  bool TransformFaceToCellParams(int bdyFace, double* pcoords);
109 
110 protected:
112  ~vtkLagrangeHexahedron() override;
113 
114  vtkHexahedron* GetApprox();
115  void PrepareApproxData(vtkPointData* pd, vtkCellData* cd, vtkIdType cellId, vtkDataArray* cellScalars);
116  vtkHexahedron* GetApproximateHex(
117  int subId, vtkDataArray* scalarsIn = nullptr, vtkDataArray* scalarsOut = nullptr);
118 
119  int Order[4];
131 
132 private:
134  void operator=(const vtkLagrangeHexahedron&) = delete;
135 };
136 
138 {
139  center[0] = center[1] = center[2] = 0.5;
140  return 0;
141 }
142 
143 #endif // vtkLagrangeHexahedron_h
represent and manipulate point attribute data
Definition: vtkPointData.h:37
vtkNew< vtkLagrangeInterpolation > Interp
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
vtkSmartPointer< vtkPoints > PointParametricCoordinates
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:357
abstract superclass for non-linear cells
vtkNew< vtkPoints > TmpPts
vtkNew< vtkLagrangeQuadrilateral > FaceCell
int vtkIdType
Definition: vtkType.h:347
vtkSmartPointer< vtkPointData > ApproxPD
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
int GetCellType() override
Return the type of cell.
vtkNew< vtkDoubleArray > CellScalars
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:360
dynamic, self-adjusting array of double
vtkSmartPointer< vtkCellData > ApproxCD
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:59
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.
int RequiresInitialization() override
Some cells require initialization prior to access.
a simple class to control print indentation
Definition: vtkIndent.h:39
list of point or cell ids
Definition: vtkIdList.h:36
vtkNew< vtkIdList > TmpIds
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:47
int GetNumberOfEdges() override
Return the number of edges in the cell.
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.
vtkSmartPointer< vtkHexahedron > Approx
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:50
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.
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.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
vtkNew< vtkDoubleArray > Scalars
virtual void Initialize()
Definition: vtkCell.h:114
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
A 3D cell that represents an arbitrary order Lagrange hex.
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.
int GetNumberOfFaces() override
Return the number of faces in the cell.
int GetParametricCenter(double center[3]) override
Return center of the cell in parametric coordinates.
vtkNew< vtkLagrangeCurve > EdgeCell
represent and manipulate 3D points
Definition: vtkPoints.h:39