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 =========================================================================*/
34 #ifndef vtkQuadraticHexahedron_h
35 #define vtkQuadraticHexahedron_h
36 
37 #include "vtkCommonDataModelModule.h" // For export macro
38 #include "vtkNonLinearCell.h"
39 
40 class vtkQuadraticEdge;
41 class vtkQuadraticQuad;
42 class vtkHexahedron;
43 class vtkDoubleArray;
44 
45 class VTKCOMMONDATAMODEL_EXPORT vtkQuadraticHexahedron : public vtkNonLinearCell
46 {
47 public:
48  static vtkQuadraticHexahedron *New();
50  void PrintSelf(ostream& os, vtkIndent indent) override;
51 
53 
57  int GetCellType() override {return VTK_QUADRATIC_HEXAHEDRON;}
58  int GetCellDimension() override {return 3;}
59  int GetNumberOfEdges() override {return 12;}
60  int GetNumberOfFaces() override {return 6;}
61  vtkCell *GetEdge(int) override;
62  vtkCell *GetFace(int) override;
64 
65  int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) override;
66  void Contour(double value, vtkDataArray *cellScalars,
68  vtkCellArray *lines, vtkCellArray *polys,
69  vtkPointData *inPd, vtkPointData *outPd,
70  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override;
71  int EvaluatePosition(double x[3], double* closestPoint,
72  int& subId, double pcoords[3],
73  double& dist2, double *weights) override;
74  void EvaluateLocation(int& subId, double pcoords[3], double x[3],
75  double *weights) override;
76  int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override;
77  void Derivatives(int subId, double pcoords[3], double *values,
78  int dim, double *derivs) override;
79  double *GetParametricCoords() override;
80 
86  void Clip(double value, vtkDataArray *cellScalars,
87  vtkIncrementalPointLocator *locator, vtkCellArray *tetras,
88  vtkPointData *inPd, vtkPointData *outPd,
89  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
90  int insideOut) override;
91 
96  int IntersectWithLine(double p1[3], double p2[3], double tol, double& t,
97  double x[3], double pcoords[3], int& subId) override;
98 
99 
103  static void InterpolationFunctions(double pcoords[3], double weights[20]);
107  static void InterpolationDerivs(double pcoords[3], double derivs[60]);
109 
113  void InterpolateFunctions(double pcoords[3], double weights[20]) override
114  {
116  }
117  void InterpolateDerivs(double pcoords[3], double derivs[60]) override
118  {
120  }
122 
123 
127  static int *GetEdgeArray(int edgeId);
128  static int *GetFaceArray(int faceId);
130 
136  void JacobianInverse(double pcoords[3], double **inverse, double derivs[60]);
137 
138 protected:
140  ~vtkQuadraticHexahedron() override;
141 
149 
150  void Subdivide(vtkPointData *inPd, vtkCellData *inCd, vtkIdType cellId,
151  vtkDataArray *cellScalars);
152 
153 private:
155  void operator=(const vtkQuadraticHexahedron&) = delete;
156 };
157 
158 #endif
159 
160 
represent and manipulate point attribute data
Definition: vtkPointData.h:31
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
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 EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
virtual int EvaluatePosition(double x[3], double *closestPoint, 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...
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:345
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.
void InterpolateDerivs(double pcoords[3], double derivs[60]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
dynamic, self-adjusting array of double
abstract class to specify cell behavior
Definition: vtkCell.h:56
static void InterpolationDerivs(double pcoords[3], double derivs[60])
cell represents a parabolic, 8-node isoparametric quad
a simple class to control print indentation
Definition: vtkIndent.h:33
list of point or cell ids
Definition: vtkIdList.h:30
virtual void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:48
virtual int IntersectWithLine(double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:41
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.
object to represent cell connectivity
Definition: vtkCellArray.h:44
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
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.
static void InterpolationFunctions(double pcoords[3], double weights[20])
virtual int CellBoundary(int subId, double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
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.
void InterpolateFunctions(double pcoords[3], double weights[20]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
represent and manipulate 3D points
Definition: vtkPoints.h:33