VTK
vtkCell.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkCell.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 =========================================================================*/
40 #ifndef vtkCell_h
41 #define vtkCell_h
42 
43 #define VTK_CELL_SIZE 512
44 #define VTK_TOL 1.e-05 // Tolerance for geometric calculation
45 
46 #include "vtkCommonDataModelModule.h" // For export macro
47 #include "vtkObject.h"
48 
49 #include "vtkIdList.h" // Needed for inline methods
50 #include "vtkCellType.h" // Needed to define cell types
51 
52 class vtkCellArray;
53 class vtkCellData;
54 class vtkDataArray;
55 class vtkPointData;
57 class vtkPoints;
58 
59 class VTKCOMMONDATAMODEL_EXPORT vtkCell : public vtkObject
60 {
61 public:
62  vtkTypeMacro(vtkCell,vtkObject);
63  void PrintSelf(ostream& os, vtkIndent indent) VTK_OVERRIDE;
64 
69  void Initialize(int npts, vtkIdType *pts, vtkPoints *p);
70 
76  virtual void ShallowCopy(vtkCell *c);
77 
82  virtual void DeepCopy(vtkCell *c);
83 
87  virtual int GetCellType() = 0;
88 
92  virtual int GetCellDimension() = 0;
93 
99  virtual int IsLinear() {return 1;}
100 
105  virtual int RequiresInitialization() {return 0;}
106  virtual void Initialize() {}
107 
113  virtual int IsExplicitCell() {return 0;}
114 
120  virtual int RequiresExplicitFaceRepresentation() {return 0;}
121  virtual void SetFaces(vtkIdType *vtkNotUsed(faces)) {}
122  virtual vtkIdType *GetFaces() {return NULL;}
123 
127  vtkPoints *GetPoints() {return this->Points;}
128 
132  vtkIdType GetNumberOfPoints() {return this->PointIds->GetNumberOfIds();}
133 
137  virtual int GetNumberOfEdges() = 0;
138 
142  virtual int GetNumberOfFaces() = 0;
143 
147  vtkIdList *GetPointIds() {return this->PointIds;}
148 
152  vtkIdType GetPointId(int ptId) {return this->PointIds->GetId(ptId);}
153 
157  virtual vtkCell *GetEdge(int edgeId) = 0;
158 
162  virtual vtkCell *GetFace(int faceId) = 0;
163 
171  virtual int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) = 0;
172 
190  virtual int EvaluatePosition(double x[3], double* closestPoint,
191  int& subId, double pcoords[3],
192  double& dist2, double *weights) = 0;
193 
199  virtual void EvaluateLocation(int& subId, double pcoords[3],
200  double x[3], double *weights) = 0;
201 
215  virtual void Contour(double value, vtkDataArray *cellScalars,
216  vtkIncrementalPointLocator *locator, vtkCellArray *verts,
217  vtkCellArray *lines, vtkCellArray *polys,
218  vtkPointData *inPd, vtkPointData *outPd,
219  vtkCellData *inCd, vtkIdType cellId,
220  vtkCellData *outCd) = 0;
221 
234  virtual void Clip(double value, vtkDataArray *cellScalars,
235  vtkIncrementalPointLocator *locator, vtkCellArray *connectivity,
236  vtkPointData *inPd, vtkPointData *outPd,
237  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
238  int insideOut) = 0;
239 
248  virtual int IntersectWithLine(double p1[3], double p2[3],
249  double tol, double& t, double x[3],
250  double pcoords[3], int& subId) = 0;
251 
262  virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) = 0;
263 
278  virtual void Derivatives(int subId, double pcoords[3], double *values,
279  int dim, double *derivs) = 0;
280 
281 
286  void GetBounds(double bounds[6]);
287 
288 
293  double *GetBounds();
294 
295 
299  double GetLength2();
300 
301 
308  virtual int GetParametricCenter(double pcoords[3]);
309 
310 
318  virtual double GetParametricDistance(double pcoords[3]);
319 
320 
328  virtual int IsPrimaryCell() {return 1;}
329 
330 
340  virtual double *GetParametricCoords();
341 
347  virtual void InterpolateFunctions(double vtkNotUsed(pcoords)[3], double* vtkNotUsed(weight))
348  {
349  }
350  virtual void InterpolateDerivs(double vtkNotUsed(pcoords)[3], double* vtkNotUsed(derivs))
351  {
352  }
353 
354  // left public for quick computational access
357 
358 protected:
360  ~vtkCell() VTK_OVERRIDE;
361 
362  double Bounds[6];
363 
364 private:
365  vtkCell(const vtkCell&) VTK_DELETE_FUNCTION;
366  void operator=(const vtkCell&) VTK_DELETE_FUNCTION;
367 };
368 
369 #endif
370 
371 
vtkPoints
represent and manipulate 3D points
Definition: vtkPoints.h:40
vtkCell::Initialize
void Initialize(int npts, vtkIdType *pts, vtkPoints *p)
Initialize cell from outside with point ids and point coordinates specified.
vtkCell::CellBoundary
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...
vtkCell::GetBounds
double * GetBounds()
Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax).
vtkCell::DeepCopy
virtual void DeepCopy(vtkCell *c)
Copy this cell by completely copying internal data structures.
vtkCell::Contour
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.
vtkCell::ShallowCopy
virtual void ShallowCopy(vtkCell *c)
Copy this cell by reference counting the internal data structures.
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:38
vtkX3D::value
@ value
Definition: vtkX3D.h:220
vtkCell::RequiresExplicitFaceRepresentation
virtual int RequiresExplicitFaceRepresentation()
Determine whether the cell requires explicit face representation, and methods for setting and getting...
Definition: vtkCell.h:120
vtkIdType
int vtkIdType
Definition: vtkType.h:287
vtkCell::Points
vtkPoints * Points
Definition: vtkCell.h:355
vtkCell::GetCellDimension
virtual int GetCellDimension()=0
Return the topological dimensional of the cell (0,1,2, or 3).
vtkCell::Initialize
virtual void Initialize()
Definition: vtkCell.h:106
vtkX3D::weight
@ weight
Definition: vtkX3D.h:532
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:60
vtkCell::GetParametricDistance
virtual double GetParametricDistance(double pcoords[3])
Return the distance of the parametric coordinate provided to the cell.
vtkDataArray
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
vtkCell::GetCellType
virtual int GetCellType()=0
Return the type of cell.
vtkCell::InterpolateFunctions
virtual void InterpolateFunctions(double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
Compute the interpolation functions/derivatives (aka shape functions/derivatives) No-ops at this leve...
Definition: vtkCell.h:347
vtkCell::EvaluateLocation
virtual void EvaluateLocation(int &subId, double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
vtkCell::IsLinear
virtual int IsLinear()
Non-linear cells require special treatment beyond the usual cell type and connectivity list informati...
Definition: vtkCell.h:99
vtkCell::~vtkCell
~vtkCell() override
vtkCell::Triangulate
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
vtkCell::GetPointIds
vtkIdList * GetPointIds()
Return the list of point ids defining the cell.
Definition: vtkCell.h:147
vtkCell::GetFaces
virtual vtkIdType * GetFaces()
Definition: vtkCell.h:122
vtkCell::InterpolateDerivs
virtual void InterpolateDerivs(double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:350
vtkCell::GetBounds
void GetBounds(double bounds[6])
Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax).
vtkCell
abstract class to specify cell behavior
Definition: vtkCell.h:60
vtkCell::Derivatives
virtual void Derivatives(int subId, double pcoords[3], double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:39
vtkCell::GetFace
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:40
vtkCell::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:51
vtkCellType.h
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:37
vtkCell::GetParametricCoords
virtual double * GetParametricCoords()
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkCell::IsPrimaryCell
virtual int IsPrimaryCell()
Return whether this cell type has a fixed topology or whether the topology varies depending on the da...
Definition: vtkCell.h:328
vtkCell::GetNumberOfEdges
virtual int GetNumberOfEdges()=0
Return the number of edges in the cell.
vtkCell::GetNumberOfFaces
virtual int GetNumberOfFaces()=0
Return the number of faces in the cell.
vtkCell::IsExplicitCell
virtual int IsExplicitCell()
Explicit cells require additional representational information beyond the usual cell type and connect...
Definition: vtkCell.h:113
vtkCell::GetNumberOfPoints
vtkIdType GetNumberOfPoints()
Return the number of points in the cell.
Definition: vtkCell.h:132
vtkObject.h
vtkCell::RequiresInitialization
virtual int RequiresInitialization()
Some cells require initialization prior to access.
Definition: vtkCell.h:105
vtkCell::SetFaces
virtual void SetFaces(vtkIdType *vtkNotUsed(faces))
Definition: vtkCell.h:121
vtkCell::PointIds
vtkIdList * PointIds
Definition: vtkCell.h:356
vtkCell::Clip
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.
vtkCell::GetEdge
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
vtkCell::GetPointId
vtkIdType GetPointId(int ptId)
For cell point i, return the actual point id.
Definition: vtkCell.h:152
vtkCell::GetParametricCenter
virtual int GetParametricCenter(double pcoords[3])
Return center of the cell in parametric coordinates.
vtkCell::GetLength2
double GetLength2()
Compute Length squared of cell (i.e., bounding box diagonal squared).
vtkCell::GetPoints
vtkPoints * GetPoints()
Get the point coordinates for the cell.
Definition: vtkCell.h:127
vtkX3D::index
@ index
Definition: vtkX3D.h:246
vtkIdList.h
vtkCell::vtkCell
vtkCell()
vtkCell::EvaluatePosition
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...
vtkCell::IntersectWithLine
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.