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 =========================================================================*/
42 #ifndef vtkCell_h
43 #define vtkCell_h
44 
45 #define VTK_CELL_SIZE 512
46 #define VTK_TOL 1.e-05 // Tolerance for geometric calculation
47 
48 #include "vtkCommonDataModelModule.h" // For export macro
49 #include "vtkObject.h"
50 
51 #include "vtkIdList.h" // Needed for inline methods
52 #include "vtkCellType.h" // Needed to define cell types
53 
54 class vtkCellArray;
55 class vtkCellData;
56 class vtkDataArray;
57 class vtkPointData;
59 class vtkPoints;
60 
62 {
63 public:
64  vtkTypeMacro(vtkCell,vtkObject);
65  void PrintSelf(ostream& os, vtkIndent indent);
66 
69  void Initialize(int npts, vtkIdType *pts, vtkPoints *p);
70 
74  virtual void ShallowCopy(vtkCell *c);
75 
78  virtual void DeepCopy(vtkCell *c);
79 
81  virtual int GetCellType() = 0;
82 
84  virtual int GetCellDimension() = 0;
85 
89  virtual int IsLinear() {return 1;}
90 
92 
95  virtual int RequiresInitialization() {return 0;}
96  virtual void Initialize() {}
98 
102  virtual int IsExplicitCell() {return 0;}
103 
105 
108  virtual int RequiresExplicitFaceRepresentation() {return 0;}
109  virtual void SetFaces(vtkIdType *vtkNotUsed(faces)) {}
110  virtual vtkIdType *GetFaces() {return NULL;}
112 
114  vtkPoints *GetPoints() {return this->Points;}
115 
117  vtkIdType GetNumberOfPoints() {return this->PointIds->GetNumberOfIds();}
118 
120  virtual int GetNumberOfEdges() = 0;
121 
123  virtual int GetNumberOfFaces() = 0;
124 
126  vtkIdList *GetPointIds() {return this->PointIds;}
127 
129  vtkIdType GetPointId(int ptId) {return this->PointIds->GetId(ptId);}
130 
132  virtual vtkCell *GetEdge(int edgeId) = 0;
133 
135  virtual vtkCell *GetFace(int faceId) = 0;
136 
142  virtual int CellBoundary(int subId, double pcoords[3], vtkIdList *pts) = 0;
143 
145 
160  virtual int EvaluatePosition(double x[3], double* closestPoint,
161  int& subId, double pcoords[3],
162  double& dist2, double *weights) = 0;
164 
166 
169  virtual void EvaluateLocation(int& subId, double pcoords[3],
170  double x[3], double *weights) = 0;
172 
174 
185  virtual void Contour(double value, vtkDataArray *cellScalars,
186  vtkIncrementalPointLocator *locator, vtkCellArray *verts,
187  vtkCellArray *lines, vtkCellArray *polys,
188  vtkPointData *inPd, vtkPointData *outPd,
189  vtkCellData *inCd, vtkIdType cellId,
190  vtkCellData *outCd) = 0;
192 
194 
205  virtual void Clip(double value, vtkDataArray *cellScalars,
206  vtkIncrementalPointLocator *locator, vtkCellArray *connectivity,
207  vtkPointData *inPd, vtkPointData *outPd,
208  vtkCellData *inCd, vtkIdType cellId, vtkCellData *outCd,
209  int insideOut) = 0;
211 
213 
221  virtual int IntersectWithLine(double p1[3], double p2[3],
222  double tol, double& t, double x[3],
223  double pcoords[3], int& subId) = 0;
225 
235  virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts) = 0;
236 
238 
250  virtual void Derivatives(int subId, double pcoords[3], double *values,
251  int dim, double *derivs) = 0;
253 
254 
257  void GetBounds(double bounds[6]);
258 
259 
262  double *GetBounds();
263 
264 
266  double GetLength2();
267 
268 
274  virtual int GetParametricCenter(double pcoords[3]);
275 
276 
282  virtual double GetParametricDistance(double pcoords[3]);
283 
284 
290  virtual int IsPrimaryCell() {return 1;}
291 
292 
300  virtual double *GetParametricCoords();
301 
303 
306  virtual void InterpolateFunctions(double vtkNotUsed(pcoords)[3], double* vtkNotUsed(weight))
307  {
308  }
309  virtual void InterpolateDerivs(double vtkNotUsed(pcoords)[3], double* vtkNotUsed(derivs))
310  {
311  }
313 
314  // left public for quick computational access
317 
318 protected:
319  vtkCell();
320  ~vtkCell();
321 
322  double Bounds[6];
323 
324 private:
325  vtkCell(const vtkCell&); // Not implemented.
326  void operator=(const vtkCell&); // Not implemented.
327 };
328 
329 #endif
330 
331 
vtkIdList * PointIds
Definition: vtkCell.h:316
vtkIdType GetNumberOfPoints()
Definition: vtkCell.h:117
abstract base class for most VTK objects
Definition: vtkObject.h:61
represent and manipulate point attribute data
Definition: vtkPointData.h:36
virtual void InterpolateFunctions(double vtkNotUsed(pcoords)[3], double *vtkNotUsed(weight))
Definition: vtkCell.h:306
vtkPoints * GetPoints()
Definition: vtkCell.h:114
represent and manipulate cell attribute data
Definition: vtkCellData.h:37
Abstract class in support of both point location and point insertion.
void DeepCopy(vtkPistonReference *self, vtkPistonReference *other)
int vtkIdType
Definition: vtkType.h:275
vtkIdType GetPointId(int ptId)
Definition: vtkCell.h:129
virtual int IsLinear()
Definition: vtkCell.h:89
virtual void SetFaces(vtkIdType *vtkNotUsed(faces))
Definition: vtkCell.h:109
abstract class to specify cell behavior
Definition: vtkCell.h:61
virtual void PrintSelf(ostream &os, vtkIndent indent)
vtkPoints * Points
Definition: vtkCell.h:315
a simple class to control print indentation
Definition: vtkIndent.h:38
list of point or cell ids
Definition: vtkIdList.h:35
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:54
virtual int RequiresInitialization()
Definition: vtkCell.h:95
virtual int RequiresExplicitFaceRepresentation()
Definition: vtkCell.h:108
virtual void InterpolateDerivs(double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:309
object to represent cell connectivity
Definition: vtkCellArray.h:49
virtual int IsExplicitCell()
Definition: vtkCell.h:102
virtual int IsPrimaryCell()
Definition: vtkCell.h:290
int Contour(vtkDataSet *input, vtkPolyData *output, vtkDataArray *field, float isoValue, bool computeScalars)
vtkIdList * GetPointIds()
Definition: vtkCell.h:126
virtual vtkIdType * GetFaces()
Definition: vtkCell.h:110
virtual void Initialize()
Definition: vtkCell.h:96
#define VTKCOMMONDATAMODEL_EXPORT
represent and manipulate 3D points
Definition: vtkPoints.h:38