VTK  9.3.20240419
vtkCellLocator.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
115 #ifndef vtkCellLocator_h
116 #define vtkCellLocator_h
117 
118 #include "vtkAbstractCellLocator.h"
119 #include "vtkCommonDataModelModule.h" // For export macro
120 #include "vtkNew.h" // For vtkNew
121 
122 VTK_ABI_NAMESPACE_BEGIN
123 class vtkIntArray;
124 
125 class VTKCOMMONDATAMODEL_EXPORT vtkCellLocator : public vtkAbstractCellLocator
126 {
127 public:
129 
133  void PrintSelf(ostream& os, vtkIndent indent) override;
135 
140  static vtkCellLocator* New();
141 
147 
148  // Re-use any superclass signatures that we don't override.
153 
160  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
161  double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
162 
172  int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints* points,
173  vtkIdList* cellIds, vtkGenericCell* cell) override;
174 
184  void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
185  vtkIdType& cellId, int& subId, double& dist2) override
186  {
187  this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
188  }
189 
202  vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
203  vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
204 
212  vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell* GenCell, int& subId,
213  double pcoords[3], double* weights) override;
214 
219  void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
220 
230  const double p1[3], const double p2[3], double tolerance, vtkIdList* cellsIds) override
231  {
232  this->Superclass::FindCellsAlongLine(p1, p2, tolerance, cellsIds);
233  }
234 
236 
239  void FreeSearchStructure() override;
240  void BuildLocator() override;
241  void ForceBuildLocator() override;
242  void GenerateRepresentation(int level, vtkPolyData* pd) override;
244 
248  virtual vtkIdList* GetCells(int bucket);
249 
254  virtual int GetNumberOfBuckets();
255 
261  void ShallowCopy(vtkAbstractCellLocator* locator) override;
262 
263 protected:
265  ~vtkCellLocator() override;
266 
267  void BuildLocatorInternal() override;
268 
269  //------------------------------------------------------------------------------
271  {
272  public:
274 
275  inline int GetNumberOfNeighbors();
276 
277  inline void Reset();
278 
279  inline int* GetPoint(int i);
280 
281  inline int InsertNextPoint(int* x);
282 
283  protected:
285  };
286 
287  void GetOverlappingBuckets(vtkNeighborCells& buckets, const double x[3], double dist,
288  int prevMinLevel[3], int prevMaxLevel[3]);
289 
290  inline void GetBucketIndices(const double x[3], int ijk[3]);
291 
292  double Distance2ToBucket(const double x[3], int nei[3]);
293  double Distance2ToBounds(const double x[3], double bounds[6]);
294 
295  int NumberOfOctants; // number of octants in tree
296  double Bounds[6]; // bounding box root octant
297  double H[3]; // width of leaf octant in x-y-z directions
298  int NumberOfDivisions; // number of "leaf" octant sub-divisions
299  std::shared_ptr<std::vector<vtkSmartPointer<vtkIdList>>> TreeSharedPtr;
301 
302  void MarkParents(const vtkSmartPointer<vtkIdList>&, int, int, int, int, int);
303  int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType& idx);
305  int face, int numDivs, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
306  void ComputeOctantBounds(double octantBounds[6], int i, int j, int k);
307 
308 private:
309  vtkCellLocator(const vtkCellLocator&) = delete;
310  void operator=(const vtkCellLocator&) = delete;
311 };
312 
313 VTK_ABI_NAMESPACE_END
314 #endif
an abstract base class for locators which find cells
virtual vtkIdType FindCell(double x[3])
Returns the Id of the cell containing the point, returns -1 if no cell found.
virtual void SetNumberOfCellsPerNode(int)
Specify the preferred/maximum number of cells in each node/bucket.
virtual void FindClosestPoint(const double x[3], double closestPoint[3], vtkIdType &cellId, int &subId, double &dist2)
Return the closest point and the cell which is closest to the point x.
virtual vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3], vtkIdType &cellId, int &subId, double &dist2)
Return the closest point within a specified radius and the cell which is closest to the point x.
virtual int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)
Return intersection point (if any) of finite line with cells contained in cell locator.
object to represent cell connectivity
Definition: vtkCellArray.h:286
vtkNew< vtkIntArray > Points
octree-based spatial search object to quickly locate cells
void MarkParents(const vtkSmartPointer< vtkIdList > &, int, int, int, int, int)
void FreeSearchStructure() override
Satisfy vtkLocator abstract interface.
~vtkCellLocator() override
double Distance2ToBounds(const double x[3], double bounds[6])
int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType &idx)
int GetNumberOfCellsPerBucket()
void GetBucketIndices(const double x[3], int ijk[3])
void GenerateRepresentation(int level, vtkPolyData *pd) override
Satisfy vtkLocator abstract interface.
void GetOverlappingBuckets(vtkNeighborCells &buckets, const double x[3], double dist, int prevMinLevel[3], int prevMaxLevel[3])
virtual int GetNumberOfBuckets()
Return number of buckets available.
void SetNumberOfCellsPerBucket(int N)
Specify the average number of cells in each octant.
void FindCellsWithinBounds(double *bbox, vtkIdList *cells) override
Return a list of unique cell ids inside of a given bounding box.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods to print and obtain type-related information.
void ComputeOctantBounds(double octantBounds[6], int i, int j, int k)
double Distance2ToBucket(const double x[3], int nei[3])
int IntersectWithLine(const double p1[3], const double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId, vtkIdType &cellId, vtkGenericCell *cell) override
Return intersection point (if any) AND the cell which was intersected by the finite line.
std::shared_ptr< std::vector< vtkSmartPointer< vtkIdList > > > TreeSharedPtr
static vtkCellLocator * New()
Construct with automatic computation of divisions, averaging 25 cells per bucket.
void GenerateFace(int face, int numDivs, int i, int j, int k, vtkPoints *pts, vtkCellArray *polys)
int IntersectWithLine(const double p1[3], const double p2[3], double tol, vtkPoints *points, vtkIdList *cellIds, vtkGenericCell *cell) override
Take the passed line segment and intersect it with the data set.
virtual vtkIdList * GetCells(int bucket)
Get the cells in a particular bucket.
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
void FindCellsAlongLine(const double p1[3], const double p2[3], double tolerance, vtkIdList *cellsIds) override
Take the passed line segment and intersect it with the data set.
void ForceBuildLocator() override
Satisfy vtkLocator abstract interface.
void ShallowCopy(vtkAbstractCellLocator *locator) override
Shallow copy of a vtkCellLocator.
vtkSmartPointer< vtkIdList > * Tree
vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2, int &inside) override
Return the closest point within a specified radius and the cell which is closest to the point x.
void BuildLocator() override
Satisfy vtkLocator abstract interface.
void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2) override
Return the closest point and the cell which is closest to the point x.
vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell *GenCell, int &subId, double pcoords[3], double *weights) override
Find the cell containing a given point.
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:133
a simple class to control print indentation
Definition: vtkIndent.h:108
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:144
represent and manipulate 3D points
Definition: vtkPoints.h:139
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:181
@ points
Definition: vtkX3D.h:446
@ level
Definition: vtkX3D.h:395
@ radius
Definition: vtkX3D.h:252
@ size
Definition: vtkX3D.h:253
@ offset
Definition: vtkX3D.h:438
int vtkIdType
Definition: vtkType.h:315