VTK  9.2.20230323
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 =========================================================================*/
129 #ifndef vtkCell_h
130 #define vtkCell_h
131 
132 #define VTK_CELL_SIZE 512
133 #define VTK_TOL 1.e-05 // Tolerance for geometric calculation
134 
135 #include "vtkCommonDataModelModule.h" // For export macro
136 #include "vtkObject.h"
137 
138 #include "vtkBoundingBox.h" // Needed for IntersectWithCell
139 #include "vtkCellType.h" // Needed to define cell types
140 #include "vtkIdList.h" // Needed for inline methods
141 
142 VTK_ABI_NAMESPACE_BEGIN
143 class vtkCellArray;
144 class vtkCellData;
145 class vtkDataArray;
146 class vtkPointData;
148 class vtkPoints;
149 
150 class VTKCOMMONDATAMODEL_EXPORT vtkCell : public vtkObject
151 {
152 public:
153  vtkTypeMacro(vtkCell, vtkObject);
154  void PrintSelf(ostream& os, vtkIndent indent) override;
155 
160  void Initialize(int npts, const vtkIdType* pts, vtkPoints* p);
161 
168  void Initialize(int npts, vtkPoints* p);
169 
175  virtual void ShallowCopy(vtkCell* c);
176 
181  virtual void DeepCopy(vtkCell* c);
182 
186  virtual int GetCellType() = 0;
187 
191  virtual int GetCellDimension() = 0;
192 
198  virtual int IsLinear() { return 1; }
199 
204  virtual int RequiresInitialization() { return 0; }
205  virtual void Initialize() {}
206 
212  virtual int IsExplicitCell() { return 0; }
213 
219  virtual int RequiresExplicitFaceRepresentation() { return 0; }
220  virtual void SetFaces(vtkIdType* vtkNotUsed(faces)) {}
221  virtual vtkIdType* GetFaces() { return nullptr; }
222 
226  vtkPoints* GetPoints() { return this->Points; }
227 
231  vtkIdType GetNumberOfPoints() const { return this->PointIds->GetNumberOfIds(); }
232 
236  virtual int GetNumberOfEdges() = 0;
237 
241  virtual int GetNumberOfFaces() = 0;
242 
246  vtkIdList* GetPointIds() { return this->PointIds; }
247 
251  vtkIdType GetPointId(int ptId) VTK_EXPECTS(0 <= ptId && ptId < GetPointIds()->GetNumberOfIds())
252  {
253  return this->PointIds->GetId(ptId);
254  }
255 
259  virtual vtkCell* GetEdge(int edgeId) = 0;
260 
272  virtual vtkCell* GetFace(int faceId) = 0;
273 
281  virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) = 0;
282 
300  virtual int EvaluatePosition(const double x[3], double closestPoint[3], int& subId,
301  double pcoords[3], double& dist2, double weights[]) = 0;
302 
308  virtual void EvaluateLocation(
309  int& subId, const double pcoords[3], double x[3], double* weights) = 0;
310 
324  virtual void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
325  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
326  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) = 0;
327 
340  virtual void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
341  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
342  vtkIdType cellId, vtkCellData* outCd, int insideOut) = 0;
343 
357  virtual int Inflate(double dist);
358 
367  virtual double ComputeBoundingSphere(double center[3]) const;
368 
377  virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
378  double x[3], double pcoords[3], int& subId) = 0;
379 
381 
388  virtual int IntersectWithCell(vtkCell* other, double tol = 0.0);
389  virtual int IntersectWithCell(vtkCell* other, const vtkBoundingBox& boudingBox,
390  const vtkBoundingBox& otherBoundingBox, double tol = 0.0);
392 
403  virtual int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) = 0;
404 
419  virtual void Derivatives(
420  int subId, const double pcoords[3], const double* values, int dim, double* derivs) = 0;
421 
426  void GetBounds(double bounds[6]);
427 
432  double* GetBounds() VTK_SIZEHINT(6);
433 
437  double GetLength2();
438 
445  virtual int GetParametricCenter(double pcoords[3]);
446 
454  virtual double GetParametricDistance(const double pcoords[3]);
455 
463  virtual int IsPrimaryCell() { return 1; }
464 
474  virtual double* GetParametricCoords() VTK_SIZEHINT(3 * GetNumberOfPoints());
475 
481  virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(weight))
482  {
483  }
484  virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(derivs)) {}
485 
486  // left public for quick computational access
489 
490 protected:
492  ~vtkCell() override;
493 
494  double Bounds[6];
495 
496 private:
497  vtkCell(const vtkCell&) = delete;
498  void operator=(const vtkCell&) = delete;
499 };
500 
501 VTK_ABI_NAMESPACE_END
502 #endif
Fast, simple class for representing and operating on 3D bounds.
object to represent cell connectivity
Definition: vtkCellArray.h:297
represent and manipulate cell attribute data
Definition: vtkCellData.h:152
abstract class to specify cell behavior
Definition: vtkCell.h:151
virtual vtkCell * GetEdge(int edgeId)=0
Return the edge cell from the edgeId of the cell.
virtual void SetFaces(vtkIdType *vtkNotUsed(faces))
Definition: vtkCell.h:220
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:251
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 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:219
virtual void Initialize()
Definition: vtkCell.h:205
void Initialize(int npts, vtkPoints *p)
Initialize the cell with point coordinates specified.
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:246
~vtkCell() override
virtual int RequiresInitialization()
Some cells require initialization prior to access.
Definition: vtkCell.h:204
vtkPoints * Points
Definition: vtkCell.h:487
virtual int IsLinear()
Non-linear cells require special treatment beyond the usual cell type and connectivity list informati...
Definition: vtkCell.h:198
vtkPoints * GetPoints()
Get the point coordinates for the cell.
Definition: vtkCell.h:226
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).
virtual vtkIdType * GetFaces()
Definition: vtkCell.h:221
vtkIdList * PointIds
Definition: vtkCell.h:488
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:231
virtual int IsExplicitCell()
Explicit cells require additional representational information beyond the usual cell type and connect...
Definition: vtkCell.h:212
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:484
virtual int Triangulate(int index, vtkIdList *ptIds, vtkPoints *pts)=0
Generate simplices of proper dimension.
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:166
list of point or cell ids
Definition: vtkIdList.h:144
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:120
abstract base class for most VTK objects
Definition: vtkObject.h:83
represent and manipulate point attribute data
Definition: vtkPointData.h:151
represent and manipulate 3D points
Definition: vtkPoints.h:150
@ value
Definition: vtkX3D.h:232
@ center
Definition: vtkX3D.h:242
@ weight
Definition: vtkX3D.h:544
@ index
Definition: vtkX3D.h:258
int vtkIdType
Definition: vtkType.h:327
#define VTK_SIZEHINT(...)
#define VTK_EXPECTS(x)