VTK  9.3.20240327
vtkCell.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
107 #ifndef vtkCell_h
108 #define vtkCell_h
109 
110 #define VTK_CELL_SIZE 512
111 #define VTK_TOL 1.e-05 // Tolerance for geometric calculation
112 
113 #include "vtkCommonDataModelModule.h" // For export macro
114 #include "vtkObject.h"
115 
116 #include "vtkBoundingBox.h" // Needed for IntersectWithCell
117 #include "vtkCellType.h" // Needed to define cell types
118 #include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_4_0
119 #include "vtkIdList.h" // Needed for inline methods
120 
121 VTK_ABI_NAMESPACE_BEGIN
122 class vtkCellArray;
123 class vtkCellData;
124 class vtkDataArray;
125 class vtkPointData;
127 class vtkPoints;
128 
129 class VTKCOMMONDATAMODEL_EXPORT vtkCell : public vtkObject
130 {
131 public:
132  vtkTypeMacro(vtkCell, vtkObject);
133  void PrintSelf(ostream& os, vtkIndent indent) override;
134 
139  void Initialize(int npts, const vtkIdType* pts, vtkPoints* p);
140 
147  void Initialize(int npts, vtkPoints* p);
148 
154  virtual void ShallowCopy(vtkCell* c);
155 
160  virtual void DeepCopy(vtkCell* c);
161 
165  virtual int GetCellType() = 0;
166 
170  virtual int GetCellDimension() = 0;
171 
177  virtual int IsLinear() { return 1; }
178 
183  virtual int RequiresInitialization() { return 0; }
184  virtual void Initialize() {}
185 
191  virtual int IsExplicitCell() { return 0; }
192 
198  virtual int RequiresExplicitFaceRepresentation() { return 0; }
199 
200  VTK_DEPRECATED_IN_9_4_0("Use SetCellFaces() after casting the cell to vtkPolyhedron.")
201  virtual void SetFaces(vtkIdType* vtkNotUsed(faces)) {}
202  VTK_DEPRECATED_IN_9_4_0("Use GetCellFaces() after casting the cell to vtkPolyhedron.")
203  virtual vtkIdType* GetFaces() { return nullptr; }
204 
208  vtkPoints* GetPoints() { return this->Points; }
209 
213  vtkIdType GetNumberOfPoints() const { return this->PointIds->GetNumberOfIds(); }
214 
218  virtual int GetNumberOfEdges() = 0;
219 
223  virtual int GetNumberOfFaces() = 0;
224 
228  vtkIdList* GetPointIds() { return this->PointIds; }
229 
233  vtkIdType GetPointId(int ptId) VTK_EXPECTS(0 <= ptId && ptId < GetPointIds()->GetNumberOfIds())
234  {
235  return this->PointIds->GetId(ptId);
236  }
237 
241  virtual vtkCell* GetEdge(int edgeId) = 0;
242 
254  virtual vtkCell* GetFace(int faceId) = 0;
255 
263  virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) = 0;
264 
282  virtual int EvaluatePosition(const double x[3], double closestPoint[3], int& subId,
283  double pcoords[3], double& dist2, double weights[]) = 0;
284 
290  virtual void EvaluateLocation(
291  int& subId, const double pcoords[3], double x[3], double* weights) = 0;
292 
306  virtual void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
307  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
308  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) = 0;
309 
322  virtual void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
323  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
324  vtkIdType cellId, vtkCellData* outCd, int insideOut) = 0;
325 
339  virtual int Inflate(double dist);
340 
349  virtual double ComputeBoundingSphere(double center[3]) const;
350 
359  virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
360  double x[3], double pcoords[3], int& subId) = 0;
361 
363 
370  virtual int IntersectWithCell(vtkCell* other, double tol = 0.0);
371  virtual int IntersectWithCell(vtkCell* other, const vtkBoundingBox& boudingBox,
372  const vtkBoundingBox& otherBoundingBox, double tol = 0.0);
374 
385  virtual int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts);
386 
397  virtual int TriangulateIds(int index, vtkIdList* ptIds);
398 
410  virtual int TriangulateLocalIds(int index, vtkIdList* ptIds) = 0;
411 
426  virtual void Derivatives(
427  int subId, const double pcoords[3], const double* values, int dim, double* derivs) = 0;
428 
433  void GetBounds(double bounds[6]);
434 
439  double* GetBounds() VTK_SIZEHINT(6);
440 
444  double GetLength2();
445 
452  virtual int GetParametricCenter(double pcoords[3]);
453 
461  virtual double GetParametricDistance(const double pcoords[3]);
462 
470  virtual int IsPrimaryCell() { return 1; }
471 
481  virtual double* GetParametricCoords() VTK_SIZEHINT(3 * GetNumberOfPoints());
482 
488  virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(weight))
489  {
490  }
491  virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(derivs)) {}
492 
493  // left public for quick computational access
496 
497 protected:
499  ~vtkCell() override;
500 
501  double Bounds[6];
502 
503 private:
504  vtkCell(const vtkCell&) = delete;
505  void operator=(const vtkCell&) = delete;
506 };
507 
508 VTK_ABI_NAMESPACE_END
509 #endif
Fast, simple class for representing and operating on 3D bounds.
object to represent cell connectivity
Definition: vtkCellArray.h:285
represent and manipulate cell attribute data
Definition: vtkCellData.h:140
abstract class to specify cell behavior
Definition: vtkCell.h:130
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
virtual int GetNumberOfEdges()=0
Return the number of edges in the cell.
virtual vtkCell * GetFace(int faceId)=0
Return the face cell from the faceId of the cell.
virtual void Derivatives(int subId, const double pcoords[3], const double *values, int dim, double *derivs)=0
Compute derivatives given cell subId and parametric coordinates.
virtual int EvaluatePosition(const double x[3], double closestPoint[3], 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...
vtkIdType GetPointId(int ptId)
For cell point i, return the actual point id.
Definition: vtkCell.h:233
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.
virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList *pts)=0
Given parametric coordinates of a point, return the closest cell boundary, and whether the point is i...
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual int Inflate(double dist)
Inflates the cell.
virtual int GetCellDimension()=0
Return the topological dimensional of the cell (0,1,2, or 3).
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.
virtual int GetNumberOfFaces()=0
Return the number of faces in the cell.
virtual void EvaluateLocation(int &subId, const double pcoords[3], double x[3], double *weights)=0
Determine global coordinate (x[3]) from subId and parametric coordinates.
void GetBounds(double bounds[6])
Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax).
virtual void ShallowCopy(vtkCell *c)
Copy this cell by reference counting the internal data structures.
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)
Generate simplices of proper dimension.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)=0
Intersect with a ray.
virtual double ComputeBoundingSphere(double center[3]) const
Computes the bounding sphere of the cell.
virtual int GetCellType()=0
Return the type of cell.
virtual void DeepCopy(vtkCell *c)
Copy this cell by completely copying internal data structures.
virtual int RequiresExplicitFaceRepresentation()
Determine whether the cell requires explicit face representation, and methods for setting and getting...
Definition: vtkCell.h:198
virtual void Initialize()
Definition: vtkCell.h:184
virtual int TriangulateLocalIds(int index, vtkIdList *ptIds)=0
Generate simplices of proper dimension.
void Initialize(int npts, vtkPoints *p)
Initialize the cell with point coordinates specified.
virtual int TriangulateIds(int index, vtkIdList *ptIds)
Generate simplices of proper dimension.
virtual int IntersectWithCell(vtkCell *other, const vtkBoundingBox &boudingBox, const vtkBoundingBox &otherBoundingBox, double tol=0.0)
Intersects with an other cell.
virtual double * GetParametricCoords()
Return a contiguous array of parametric coordinates of the points defining this cell.
vtkIdList * GetPointIds()
Return the list of point ids defining the cell.
Definition: vtkCell.h:228
~vtkCell() override
virtual int RequiresInitialization()
Some cells require initialization prior to access.
Definition: vtkCell.h:183
vtkPoints * Points
Definition: vtkCell.h:494
virtual int IsLinear()
Non-linear cells require special treatment beyond the usual cell type and connectivity list informati...
Definition: vtkCell.h:177
vtkPoints * GetPoints()
Get the point coordinates for the cell.
Definition: vtkCell.h:208
void Initialize(int npts, const vtkIdType *pts, vtkPoints *p)
Initialize cell from outside with point ids and point coordinates specified.
double * GetBounds()
Compute cell bounding box (xmin,xmax,ymin,ymax,zmin,zmax).
vtkIdList * PointIds
Definition: vtkCell.h:495
virtual int IntersectWithCell(vtkCell *other, double tol=0.0)
Intersects with an other cell.
vtkIdType GetNumberOfPoints() const
Return the number of points in the cell.
Definition: vtkCell.h:213
virtual int IsExplicitCell()
Explicit cells require additional representational information beyond the usual cell type and connect...
Definition: vtkCell.h:191
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:491
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:154
list of point or cell ids
Definition: vtkIdList.h:132
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:108
abstract base class for most VTK objects
Definition: vtkObject.h:161
represent and manipulate point attribute data
Definition: vtkPointData.h:139
represent and manipulate 3D points
Definition: vtkPoints.h:138
@ value
Definition: vtkX3D.h:220
@ center
Definition: vtkX3D.h:230
@ weight
Definition: vtkX3D.h:532
@ index
Definition: vtkX3D.h:246
#define VTK_DEPRECATED_IN_9_4_0(reason)
int vtkIdType
Definition: vtkType.h:315
#define VTK_SIZEHINT(...)
#define VTK_EXPECTS(x)