VTK
vtkTriQuadraticHexahedron.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkTriQuadraticHexahedron.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 =========================================================================*/
74 #ifndef vtkTriQuadraticHexahedron_h
75 #define vtkTriQuadraticHexahedron_h
76 
77 #include "vtkCommonDataModelModule.h" // For export macro
78 #include "vtkNonLinearCell.h"
79 
80 class vtkQuadraticEdge;
81 class vtkBiQuadraticQuad;
82 class vtkHexahedron;
83 class vtkDoubleArray;
84 
85 class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticHexahedron : public vtkNonLinearCell
86 {
87 public:
88  static vtkTriQuadraticHexahedron *New ();
90  void PrintSelf (ostream & os, vtkIndent indent) override;
91 
93 
97  int GetCellType() override { return VTK_TRIQUADRATIC_HEXAHEDRON; }
98  int GetCellDimension() override { return 3; }
99  int GetNumberOfEdges() override { return 12; }
100  int GetNumberOfFaces() override { return 6; }
101  vtkCell *GetEdge (int) override;
102  vtkCell *GetFace (int) override;
104 
105  int CellBoundary(int subId, const double pcoords[3], vtkIdList * pts) override;
106  void Contour (double value, vtkDataArray * cellScalars,
107  vtkIncrementalPointLocator * locator, vtkCellArray * verts,
108  vtkCellArray * lines, vtkCellArray * polys,
109  vtkPointData * inPd, vtkPointData * outPd, vtkCellData * inCd,
110  vtkIdType cellId, vtkCellData * outCd) override;
111  int EvaluatePosition(const double x[3], double *closestPoint,
112  int &subId, double pcoords[3], double &dist2, double *weights) override;
113  void EvaluateLocation(int &subId, const double pcoords[3],
114  double x[3], double *weights) override;
115  int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts) override;
116  void Derivatives(int subId, const double pcoords[3], const double *values,
117  int dim, double *derivs) override;
118  double *GetParametricCoords () override;
119 
125  void Clip (double value, vtkDataArray * cellScalars,
126  vtkIncrementalPointLocator * locator, vtkCellArray * tetras,
127  vtkPointData * inPd, vtkPointData * outPd,
128  vtkCellData * inCd, vtkIdType cellId, vtkCellData * outCd,
129  int insideOut) override;
130 
135  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t,
136  double x[3], double pcoords[3], int &subId) override;
137 
141  static void InterpolationFunctions(const double pcoords[3], double weights[27]);
145  static void InterpolationDerivs(const double pcoords[3], double derivs[81]);
147 
151  void InterpolateFunctions(const double pcoords[3], double weights[27]) override
152  {
154  }
155  void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
156  {
158  }
160 
161 
165  static int *GetEdgeArray(int edgeId);
166  static int *GetFaceArray(int faceId);
168 
174  void JacobianInverse(const double pcoords[3], double **inverse, double derivs[81]);
175 
176 protected:
178  ~vtkTriQuadraticHexahedron () override;
179 
184 
185 private:
187  void operator = (const vtkTriQuadraticHexahedron &) = delete;
188 };
189 
190 #endif
represent and manipulate point attribute data
Definition: vtkPointData.h:37
int GetNumberOfEdges() override
Implement the vtkCell API.
represent and manipulate cell attribute data
Definition: vtkCellData.h:38
cell represents a parabolic, 9-node isoparametric quad
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
static void InterpolationFunctions(const double pcoords[3], double weights[27])
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.
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.
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
void InterpolateFunctions(const double pcoords[3], double weights[27]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
int GetCellDimension() override
Implement the vtkCell API.
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:47
void InterpolateDerivs(const double pcoords[3], double derivs[81]) override
Compute the interpolation functions/derivatives (aka shape functions/derivatives) ...
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.
static void InterpolationDerivs(const double pcoords[3], double derivs[81])
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.
cell represents a parabolic, 27-node isoparametric hexahedron
cell represents a parabolic, isoparametric edge
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 GetCellType() override
Implement the vtkCell API.
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
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
int GetNumberOfFaces() override
Implement the vtkCell API.