VTK  9.0.1
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 =========================================================================*/
37 #ifndef vtkCell_h
38 #define vtkCell_h
39 
40 #define VTK_CELL_SIZE 512
41 #define VTK_TOL 1.e-05 // Tolerance for geometric calculation
42 
43 #include "vtkCommonDataModelModule.h" // For export macro
44 #include "vtkObject.h"
45 
46 #include "vtkCellType.h" // Needed to define cell types
47 #include "vtkIdList.h" // Needed for inline methods
48 
49 class vtkCellArray;
50 class vtkCellData;
51 class vtkDataArray;
52 class vtkPointData;
54 class vtkPoints;
55 
56 class VTKCOMMONDATAMODEL_EXPORT vtkCell : public vtkObject
57 {
58 public:
59  vtkTypeMacro(vtkCell, vtkObject);
60  void PrintSelf(ostream& os, vtkIndent indent) override;
61 
66  void Initialize(int npts, const vtkIdType* pts, vtkPoints* p);
67 
74  void Initialize(int npts, vtkPoints* p);
75 
81  virtual void ShallowCopy(vtkCell* c);
82 
87  virtual void DeepCopy(vtkCell* c);
88 
92  virtual int GetCellType() = 0;
93 
97  virtual int GetCellDimension() = 0;
98 
104  virtual int IsLinear() { return 1; }
105 
110  virtual int RequiresInitialization() { return 0; }
111  virtual void Initialize() {}
112 
118  virtual int IsExplicitCell() { return 0; }
119 
125  virtual int RequiresExplicitFaceRepresentation() { return 0; }
126  virtual void SetFaces(vtkIdType* vtkNotUsed(faces)) {}
127  virtual vtkIdType* GetFaces() { return nullptr; }
128 
132  vtkPoints* GetPoints() { return this->Points; }
133 
137  vtkIdType GetNumberOfPoints() { return this->PointIds->GetNumberOfIds(); }
138 
142  virtual int GetNumberOfEdges() = 0;
143 
147  virtual int GetNumberOfFaces() = 0;
148 
152  vtkIdList* GetPointIds() { return this->PointIds; }
153 
157  vtkIdType GetPointId(int ptId) VTK_EXPECTS(0 <= ptId && ptId < GetPointIds()->GetNumberOfIds())
158  {
159  return this->PointIds->GetId(ptId);
160  }
161 
165  virtual vtkCell* GetEdge(int edgeId) = 0;
166 
170  virtual vtkCell* GetFace(int faceId) = 0;
171 
179  virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) = 0;
180 
198  virtual int EvaluatePosition(const double x[3], double closestPoint[3], int& subId,
199  double pcoords[3], double& dist2, double weights[]) = 0;
200 
206  virtual void EvaluateLocation(
207  int& subId, const double pcoords[3], double x[3], double* weights) = 0;
208 
222  virtual void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
223  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
224  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) = 0;
225 
238  virtual void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
239  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
240  vtkIdType cellId, vtkCellData* outCd, int insideOut) = 0;
241 
250  virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
251  double x[3], double pcoords[3], int& subId) = 0;
252 
263  virtual int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) = 0;
264 
279  virtual void Derivatives(
280  int subId, const double pcoords[3], const double* values, int dim, double* derivs) = 0;
281 
286  void GetBounds(double bounds[6]);
287 
292  double* GetBounds() VTK_SIZEHINT(6);
293 
297  double GetLength2();
298 
305  virtual int GetParametricCenter(double pcoords[3]);
306 
314  virtual double GetParametricDistance(const double pcoords[3]);
315 
323  virtual int IsPrimaryCell() { return 1; }
324 
334  virtual double* GetParametricCoords() VTK_SIZEHINT(3 * GetNumberOfPoints());
335 
341  virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(weight))
342  {
343  }
344  virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(derivs)) {}
345 
346  // left public for quick computational access
349 
350 protected:
351  vtkCell();
352  ~vtkCell() override;
353 
354  double Bounds[6];
355 
356 private:
357  vtkCell(const vtkCell&) = delete;
358  void operator=(const vtkCell&) = delete;
359 };
360 
361 #endif
vtkIdList * PointIds
Definition: vtkCell.h:348
vtkIdType GetNumberOfPoints()
Return the number of points in the cell.
Definition: vtkCell.h:137
abstract base class for most VTK objects
Definition: vtkObject.h:62
represent and manipulate point attribute data
Definition: vtkPointData.h:31
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkPoints * GetPoints()
Get the point coordinates for the cell.
Definition: vtkCell.h:132
represent and manipulate cell attribute data
Definition: vtkCellData.h:32
Abstract class in support of both point location and point insertion.
int vtkIdType
Definition: vtkType.h:338
vtkIdType GetPointId(int ptId)
For cell point i, return the actual point id.
Definition: vtkCell.h:157
virtual int IsLinear()
Non-linear cells require special treatment beyond the usual cell type and connectivity list informati...
Definition: vtkCell.h:104
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:344
virtual void SetFaces(vtkIdType *vtkNotUsed(faces))
Definition: vtkCell.h:126
abstract class to specify cell behavior
Definition: vtkCell.h:56
vtkPoints * Points
Definition: vtkCell.h:347
a simple class to control print indentation
Definition: vtkIndent.h:33
list of point or cell ids
Definition: vtkIdList.h:30
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:49
virtual int RequiresInitialization()
Some cells require initialization prior to access.
Definition: vtkCell.h:110
#define VTK_SIZEHINT(...)
virtual int RequiresExplicitFaceRepresentation()
Determine whether the cell requires explicit face representation, and methods for setting and getting...
Definition: vtkCell.h:125
object to represent cell connectivity
Definition: vtkCellArray.h:179
virtual int IsExplicitCell()
Explicit cells require additional representational information beyond the usual cell type and connect...
Definition: vtkCell.h:118
vtkIdList * GetPointIds()
Return the list of point ids defining the cell.
Definition: vtkCell.h:152
virtual vtkIdType * GetFaces()
Definition: vtkCell.h:127
virtual void Initialize()
Definition: vtkCell.h:111
#define VTK_EXPECTS(x)
represent and manipulate 3D points
Definition: vtkPoints.h:33