VTK  9.3.1
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
28 #ifndef vtkCell_h
29 #define vtkCell_h
30 
31 #define VTK_CELL_SIZE 512
32 #define VTK_TOL 1.e-05 // Tolerance for geometric calculation
33 
34 #include "vtkCommonDataModelModule.h" // For export macro
35 #include "vtkObject.h"
36 
37 #include "vtkBoundingBox.h" // Needed for IntersectWithCell
38 #include "vtkCellType.h" // Needed to define cell types
39 #include "vtkIdList.h" // Needed for inline methods
40 
41 VTK_ABI_NAMESPACE_BEGIN
42 class vtkCellArray;
43 class vtkCellData;
44 class vtkDataArray;
45 class vtkPointData;
47 class vtkPoints;
48 
49 class VTKCOMMONDATAMODEL_EXPORT vtkCell : public vtkObject
50 {
51 public:
52  vtkTypeMacro(vtkCell, vtkObject);
53  void PrintSelf(ostream& os, vtkIndent indent) override;
54 
59  void Initialize(int npts, const vtkIdType* pts, vtkPoints* p);
60 
67  void Initialize(int npts, vtkPoints* p);
68 
74  virtual void ShallowCopy(vtkCell* c);
75 
80  virtual void DeepCopy(vtkCell* c);
81 
85  virtual int GetCellType() = 0;
86 
90  virtual int GetCellDimension() = 0;
91 
97  virtual int IsLinear() { return 1; }
98 
103  virtual int RequiresInitialization() { return 0; }
104  virtual void Initialize() {}
105 
111  virtual int IsExplicitCell() { return 0; }
112 
118  virtual int RequiresExplicitFaceRepresentation() { return 0; }
119  virtual void SetFaces(vtkIdType* vtkNotUsed(faces)) {}
120  virtual vtkIdType* GetFaces() { return nullptr; }
121 
125  vtkPoints* GetPoints() { return this->Points; }
126 
130  vtkIdType GetNumberOfPoints() const { return this->PointIds->GetNumberOfIds(); }
131 
135  virtual int GetNumberOfEdges() = 0;
136 
140  virtual int GetNumberOfFaces() = 0;
141 
145  vtkIdList* GetPointIds() { return this->PointIds; }
146 
150  vtkIdType GetPointId(int ptId) VTK_EXPECTS(0 <= ptId && ptId < GetPointIds()->GetNumberOfIds())
151  {
152  return this->PointIds->GetId(ptId);
153  }
154 
158  virtual vtkCell* GetEdge(int edgeId) = 0;
159 
171  virtual vtkCell* GetFace(int faceId) = 0;
172 
180  virtual int CellBoundary(int subId, const double pcoords[3], vtkIdList* pts) = 0;
181 
199  virtual int EvaluatePosition(const double x[3], double closestPoint[3], int& subId,
200  double pcoords[3], double& dist2, double weights[]) = 0;
201 
207  virtual void EvaluateLocation(
208  int& subId, const double pcoords[3], double x[3], double* weights) = 0;
209 
223  virtual void Contour(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
224  vtkCellArray* verts, vtkCellArray* lines, vtkCellArray* polys, vtkPointData* inPd,
225  vtkPointData* outPd, vtkCellData* inCd, vtkIdType cellId, vtkCellData* outCd) = 0;
226 
239  virtual void Clip(double value, vtkDataArray* cellScalars, vtkIncrementalPointLocator* locator,
240  vtkCellArray* connectivity, vtkPointData* inPd, vtkPointData* outPd, vtkCellData* inCd,
241  vtkIdType cellId, vtkCellData* outCd, int insideOut) = 0;
242 
256  virtual int Inflate(double dist);
257 
266  virtual double ComputeBoundingSphere(double center[3]) const;
267 
276  virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t,
277  double x[3], double pcoords[3], int& subId) = 0;
278 
280 
287  virtual int IntersectWithCell(vtkCell* other, double tol = 0.0);
288  virtual int IntersectWithCell(vtkCell* other, const vtkBoundingBox& boudingBox,
289  const vtkBoundingBox& otherBoundingBox, double tol = 0.0);
291 
302  virtual int Triangulate(int index, vtkIdList* ptIds, vtkPoints* pts) = 0;
303 
318  virtual void Derivatives(
319  int subId, const double pcoords[3], const double* values, int dim, double* derivs) = 0;
320 
325  void GetBounds(double bounds[6]);
326 
331  double* GetBounds() VTK_SIZEHINT(6);
332 
336  double GetLength2();
337 
344  virtual int GetParametricCenter(double pcoords[3]);
345 
353  virtual double GetParametricDistance(const double pcoords[3]);
354 
362  virtual int IsPrimaryCell() { return 1; }
363 
373  virtual double* GetParametricCoords() VTK_SIZEHINT(3 * GetNumberOfPoints());
374 
380  virtual void InterpolateFunctions(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(weight))
381  {
382  }
383  virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double* vtkNotUsed(derivs)) {}
384 
385  // left public for quick computational access
388 
389 protected:
390  vtkCell();
391  ~vtkCell() override;
392 
393  double Bounds[6];
394 
395 private:
396  vtkCell(const vtkCell&) = delete;
397  void operator=(const vtkCell&) = delete;
398 };
399 
400 VTK_ABI_NAMESPACE_END
401 #endif
vtkIdList * PointIds
Definition: vtkCell.h:387
void GetBounds(T a, double bds[6])
abstract base class for most VTK objects
Definition: vtkObject.h:51
represent and manipulate point attribute data
Definition: vtkPointData.h:29
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkIdType GetNumberOfPoints() const
Return the number of points in the cell.
Definition: vtkCell.h:130
vtkPoints * GetPoints()
Get the point coordinates for the cell.
Definition: vtkCell.h:125
represent and manipulate cell attribute data
Definition: vtkCellData.h:30
Abstract class in support of both point location and point insertion.
int vtkIdType
Definition: vtkType.h:315
vtkIdType GetPointId(int ptId)
For cell point i, return the actual point id.
Definition: vtkCell.h:150
virtual int IsLinear()
Non-linear cells require special treatment beyond the usual cell type and connectivity list informati...
Definition: vtkCell.h:97
virtual void InterpolateDerivs(const double vtkNotUsed(pcoords)[3], double *vtkNotUsed(derivs))
Definition: vtkCell.h:383
virtual void SetFaces(vtkIdType *vtkNotUsed(faces))
Definition: vtkCell.h:119
abstract class to specify cell behavior
Definition: vtkCell.h:49
vtkPoints * Points
Definition: vtkCell.h:386
a simple class to control print indentation
Definition: vtkIndent.h:28
list of point or cell ids
Definition: vtkIdList.h:22
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:44
virtual int RequiresInitialization()
Some cells require initialization prior to access.
Definition: vtkCell.h:103
#define VTK_SIZEHINT(...)
virtual int RequiresExplicitFaceRepresentation()
Determine whether the cell requires explicit face representation, and methods for setting and getting...
Definition: vtkCell.h:118
int GetCellType(const Ioss::ElementTopology *topology)
Returns VTK celltype for a Ioss topology element.
object to represent cell connectivity
Definition: vtkCellArray.h:175
virtual int IsExplicitCell()
Explicit cells require additional representational information beyond the usual cell type and connect...
Definition: vtkCell.h:111
vtkIdList * GetPointIds()
Return the list of point ids defining the cell.
Definition: vtkCell.h:145
virtual vtkIdType * GetFaces()
Definition: vtkCell.h:120
virtual void Initialize()
Definition: vtkCell.h:104
#define VTK_EXPECTS(x)
represent and manipulate 3D points
Definition: vtkPoints.h:28
Fast, simple class for representing and operating on 3D bounds.