VTK
vtkQuadraticHexahedron.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticHexahedron.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 vtkQuadraticHexahedron_h
38 #define vtkQuadraticHexahedron_h
39 
40 #include "vtkCommonDataModelModule.h" // For export macro
41 #include "vtkNonLinearCell.h"
42 
43 class vtkQuadraticEdge;
44 class vtkQuadraticQuad;
45 class vtkHexahedron;
46 class vtkDoubleArray;
47 
48 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticHexahedron : public vtkNonLinearCell
49 {
50 public:
51  static vtkQuadraticHexahedron *New();
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
56 
60  int GetCellType() override {return VTK_QUADRATIC_HEXAHEDRON;}
61  int GetCellDimension() override {return 3;}
62  int GetNumberOfEdges() override {return 12;}
63  int GetNumberOfFaces() override {return 6;}
64  vtkCell *GetEdge(int) override;
65  vtkCell *GetFace(int) override;
67 
68  int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts) override;
69  void Contour(double value, vtkDataArray *cellScalars,
71  vtkCellArray *lines, vtkCellArray *polys,
72  vtkPointData *inPd, vtkPointData *outPd,
73  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) 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(int& subId, const double pcoords[3], double x[3],
78  double *weights) override;
79  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
80  void Derivatives(int subId, const double pcoords[3], const double *values,
81  int dim, double *derivs) override;
82  double *GetParametricCoords() override;
83 
89  void Clip(double value, vtkDataArray *cellScalars,
90  vtkIncrementalPointLocator *locator, vtkCellArray *tetras,
91  vtkPointData *inPd, vtkPointData *outPd,
92  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
93  int insideOut) override;
94 
99  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
100  double x[3], double pcoords[3], int& subId) override;
101 
102 
106  static void InterpolationFunctions(const double pcoords[3], double weights[20]);
110  static void InterpolationDerivs(const double pcoords[3], double derivs[60]);
112 
116  void InterpolateFunctions(const double pcoords[3], double weights[20]) override
117  {
119  }
120  void InterpolateDerivs(const double pcoords[3], double derivs[60]) override
121  {
123  }
125 
126 
130  static int *GetEdgeArray(int edgeId);
131  static int *GetFaceArray(int faceId);
133 
139  void JacobianInverse(const double pcoords[3], double **inverse, double derivs[60]);
140 
141 protected:
143  ~vtkQuadraticHexahedron() override;
144 
152 
153  void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId,
154  vtkDataArray *cellScalars);
155 
156 private:
158  void operator=(const vtkQuadraticHexahedron&) = delete;
159 };
160 
161 #endif
162 
163 
void InterpolateDerivs(const double pcoords[3], double derivs[60]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
represent and manipulate point attribute data
Definition: vtkPointData.h:37
static void InterpolationDerivs(const double pcoords[3], double derivs[60])
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
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.
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 GetCellType() override
Implement the vtkCell API.
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...
dynamic, self-adjusting array of double
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.
cell represents a parabolic, 8-node isoparametric quad
a simple class to control print indentation
Definition: vtkIndent.h:39
list of point or cell ids
Definition: vtkIdList.h:36
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:47
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.
int GetNumberOfFaces() override
Implement the vtkCell API.
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
int GetNumberOfEdges() override
Implement the vtkCell API.
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.
void InterpolateFunctions(const double pcoords[3], double weights[20]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 20-node isoparametric hexahedron
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.
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
int GetCellDimension() override
Implement the vtkCell API.
virtual double * GetParametricCoords())
Return a contiguous array of parametric coordinates of the points defining this cell.
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.
represent and manipulate 3D points
Definition: vtkPoints.h:39
static void InterpolationFunctions(const double pcoords[3], double weights[20])