VTK
vtkQuadraticLinearWedge.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkQuadraticLinearWedge.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 =========================================================================*/
41 #ifndef vtkQuadraticLinearWedge_h
42 #define vtkQuadraticLinearWedge_h
43 
44 #include "vtkCommonDataModelModule.h" // For export macro
45 #include "vtkNonLinearCell.h"
46 
47 class vtkQuadraticEdge;
48 class vtkLine;
51 class vtkWedge;
52 class vtkDoubleArray;
53 
55 {
56 public:
57  static vtkQuadraticLinearWedge *New ();
59  void PrintSelf (ostream & os, vtkIndent indent);
60 
62 
65  int GetCellDimension () { return 3; }
66  int GetNumberOfEdges () { return 9; }
67  int GetNumberOfFaces () { return 5; }
68  vtkCell *GetEdge (int edgeId);
69  vtkCell *GetFace (int faceId);
71 
72  int CellBoundary (int subId, double pcoords[3], vtkIdList * pts);
73 
75 
77  void Contour (double value, vtkDataArray * cellScalars,
78  vtkIncrementalPointLocator * locator, vtkCellArray * verts,
79  vtkCellArray * lines, vtkCellArray * polys,
80  vtkPointData * inPd, vtkPointData * outPd, vtkCellData * inCd,
81  vtkIdType cellId, vtkCellData * outCd);
82  int EvaluatePosition (double x[3], double *closestPoint,
83  int &subId, double pcoords[3], double &dist2, double *weights);
84  void EvaluateLocation (int &subId, double pcoords[3], double x[3], double *weights);
85  int Triangulate (int index, vtkIdList * ptIds, vtkPoints * pts);
86  void Derivatives (int subId, double pcoords[3], double *values, int dim, double *derivs);
87  virtual double *GetParametricCoords ();
89 
91 
94  void Clip (double value, vtkDataArray * cellScalars,
95  vtkIncrementalPointLocator * locator, vtkCellArray * tetras,
96  vtkPointData * inPd, vtkPointData * outPd,
97  vtkCellData * inCd, vtkIdType cellId, vtkCellData * outCd, int insideOut);
99 
101 
103  int IntersectWithLine (double p1[3], double p2[3], double tol, double &t,
104  double x[3], double pcoords[3], int &subId);
106 
109  int GetParametricCenter (double pcoords[3]);
110 
113  static void InterpolationFunctions (double pcoords[3], double weights[15]);
116  static void InterpolationDerivs (double pcoords[3], double derivs[45]);
118 
120  virtual void InterpolateFunctions (double pcoords[3], double weights[15])
121  {
123  }
124  virtual void InterpolateDerivs (double pcoords[3], double derivs[45])
125  {
127  }
129 
130 
132  static int *GetEdgeArray(int edgeId);
133  static int *GetFaceArray(int faceId);
135 
139  void JacobianInverse (double pcoords[3], double **inverse, double derivs[45]);
140 
141 protected:
144 
150  vtkDoubleArray *Scalars; //used to avoid New/Delete in contouring/clipping
151 
152 private:
153  vtkQuadraticLinearWedge (const vtkQuadraticLinearWedge &); // Not implemented.
154  void operator = (const vtkQuadraticLinearWedge &); // Not implemented.
155 };
156 //----------------------------------------------------------------------------
157 // Return the center of the quadratic wedge in parametric coordinates.
158 inline int vtkQuadraticLinearWedge::GetParametricCenter(double pcoords[3])
159 {
160  pcoords[0] = pcoords[1] = 1./3;
161  pcoords[2] = 0.5;
162  return 0;
163 }
164 
165 
166 #endif
int GetParametricCenter(double pcoords[3])
represent and manipulate point attribute data
Definition: vtkPointData.h:36
virtual double * GetParametricCoords()
represent and manipulate cell attribute data
Definition: vtkCellData.h:37
Abstract class in support of both point location and point insertion.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
virtual void EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights)=0
virtual int EvaluatePosition(double x[3], double *closestPoint, int &subId, double pcoords[3], double &dist2, double *weights)=0
abstract superclass for non-linear cells
int vtkIdType
Definition: vtkType.h:275
static void InterpolationFunctions(double pcoords[3], double weights[15])
virtual void InterpolateDerivs(double pcoords[3], double derivs[45])
static void InterpolationDerivs(double pcoords[3], double derivs[45])
dynamic, self-adjusting array of double
cell represents a 1D line
Definition: vtkLine.h:34
abstract class to specify cell behavior
Definition: vtkCell.h:61
a simple class to control print indentation
Definition: vtkIndent.h:38
list of point or cell ids
Definition: vtkIdList.h:35
virtual void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs)=0
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
virtual int IntersectWithLine(double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
virtual void Clip(double value, vtkDataArray *cellScalars, vtkIncrementalPointLocator *locator, vtkCellArray *connectivity, vtkPointData *inPd, vtkPointData *outPd, vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd, int insideOut)=0
void PrintSelf(ostream &os, vtkIndent indent)
virtual vtkCell * GetFace(int faceId)=0
vtkQuadraticLinearQuad * Face
object to represent cell connectivity
Definition: vtkCellArray.h:49
virtual vtkCell * GetEdge(int edgeId)=0
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
cell represents a parabolic, isoparametric triangle
vtkQuadraticTriangle * TriangleFace
virtual void InterpolateFunctions(double pcoords[3], double weights[15])
virtual int CellBoundary(int subId, double pcoords[3], vtkIdList *pts)=0
static vtkObject * New()
virtual int GetParametricCenter(double pcoords[3])
cell represents a, 12-node isoparametric wedge
cell represents a quadratic-linear, 6-node isoparametric quad
#define VTKCOMMONDATAMODEL_EXPORT
a 3D cell that represents a linear wedge
Definition: vtkWedge.h:48
represent and manipulate 3D points
Definition: vtkPoints.h:38