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
82class vtkHexahedron;
83class vtkDoubleArray;
84
85class VTKCOMMONDATAMODEL_EXPORT vtkTriQuadraticHexahedron : public vtkNonLinearCell
86{
87public:
90 void PrintSelf (ostream & os, vtkIndent indent) VTK_OVERRIDE;
91
93
97 int GetCellType() VTK_OVERRIDE { return VTK_TRIQUADRATIC_HEXAHEDRON; }
98 int GetCellDimension() VTK_OVERRIDE { return 3; }
99 int GetNumberOfEdges() VTK_OVERRIDE { return 12; }
100 int GetNumberOfFaces() VTK_OVERRIDE { return 6; }
101 vtkCell *GetEdge (int) VTK_OVERRIDE;
102 vtkCell *GetFace (int) VTK_OVERRIDE;
104
105 int CellBoundary (int subId, double pcoords[3], vtkIdList * pts) VTK_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) VTK_OVERRIDE;
111 int EvaluatePosition (double x[3], double *closestPoint,
112 int &subId, double pcoords[3], double &dist2, double *weights) VTK_OVERRIDE;
113 void EvaluateLocation (int &subId, double pcoords[3],
114 double x[3], double *weights) VTK_OVERRIDE;
115 int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts) VTK_OVERRIDE;
116 void Derivatives (int subId, double pcoords[3], double *values,
117 int dim, double *derivs) VTK_OVERRIDE;
118 double *GetParametricCoords () VTK_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) VTK_OVERRIDE;
130
135 int IntersectWithLine (double p1[3], double p2[3], double tol, double &t,
136 double x[3], double pcoords[3], int &subId) VTK_OVERRIDE;
137
141 static void InterpolationFunctions (double pcoords[3], double weights[27]);
145 static void InterpolationDerivs (double pcoords[3], double derivs[81]);
147
151 void InterpolateFunctions (double pcoords[3], double weights[27]) VTK_OVERRIDE
152 {
154 }
155 void InterpolateDerivs (double pcoords[3], double derivs[81]) VTK_OVERRIDE
156 {
158 }
160
161
165 static int *GetEdgeArray(int edgeId);
166 static int *GetFaceArray(int faceId);
168
174 void JacobianInverse (double pcoords[3], double **inverse, double derivs[81]);
175
176protected:
179
184
185private:
186 vtkTriQuadraticHexahedron (const vtkTriQuadraticHexahedron &) VTK_DELETE_FUNCTION;
187 void operator = (const vtkTriQuadraticHexahedron &) VTK_DELETE_FUNCTION;
188};
189
190#endif
cell represents a parabolic, 9-node isoparametric quad
object to represent cell connectivity
Definition: vtkCellArray.h:51
represent and manipulate cell attribute data
Definition: vtkCellData.h:39
abstract class to specify cell behavior
Definition: vtkCell.h:60
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
dynamic, self-adjusting array of double
a cell that represents a linear 3D hexahedron
Definition: vtkHexahedron.h:48
list of point or cell ids
Definition: vtkIdList.h:37
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:40
abstract superclass for non-linear cells
represent and manipulate point attribute data
Definition: vtkPointData.h:38
represent and manipulate 3D points
Definition: vtkPoints.h:40
cell represents a parabolic, isoparametric edge
cell represents a parabolic, 27-node isoparametric hexahedron
static void InterpolationDerivs(double pcoords[3], double derivs[81])
int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) override
Generate simplices of proper dimension.
void Contour(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *verts, vtkCellArray *lines, vtkCellArray *polys, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd) override
Generate contouring primitives.
static void InterpolationFunctions(double pcoords[3], double weights[27])
int GetCellType() override
Implement the vtkCell API.
int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) override
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
double * GetParametricCoords() override
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCell * GetEdge(int) override
Return the edge cell from the edgeId of the cell.
vtkCell * GetFace(int) override
Return the face cell from the faceId of the cell.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static int * GetFaceArray(int faceId)
int GetNumberOfEdges() override
Return the number of edges in the cell.
int GetCellDimension() override
Return the topological dimensional of the cell (0,1,2, or 3).
int EvaluatePosition(double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights) override
Given a point x[3] return inside(=1), outside(=0) cell, or (-1) computational problem encountered; ev...
static vtkTriQuadraticHexahedron * New()
void EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights) override
Determine global coordinate (x[3]) from subId and parametric coordinates.
static int * GetEdgeArray(int edgeId)
Return the ids of the vertices defining edge/face (edgeId/‘faceId’).
void InterpolateDerivs(double pcoords[3], double derivs[81]) override
void JacobianInverse(double pcoords[3], double **inverse, double derivs[81])
Given parametric coordinates compute inverse Jacobian transformation matrix.
int GetNumberOfFaces() override
Return the number of faces in the cell.
void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs) override
Compute derivatives given cell subId and parametric coordinates.
~vtkTriQuadraticHexahedron() override
@ value
Definition: vtkX3D.h:220
@ index
Definition: vtkX3D.h:246
@ VTK_TRIQUADRATIC_HEXAHEDRON
Definition: vtkCellType.h:73
int vtkIdType
Definition: vtkType.h:287