VTK  9.2.5
vtkCellLocator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkCellLocator.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 =========================================================================*/
147 #ifndef vtkCellLocator_h
148 #define vtkCellLocator_h
149 
150 #include "vtkAbstractCellLocator.h"
151 #include "vtkCommonDataModelModule.h" // For export macro
152 #include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_2_0
153 #include "vtkNew.h" // For vtkNew
154 
155 class vtkIntArray;
156 
157 class VTKCOMMONDATAMODEL_EXPORT vtkCellLocator : public vtkAbstractCellLocator
158 {
159 public:
161 
165  void PrintSelf(ostream& os, vtkIndent indent) override;
167 
172  static vtkCellLocator* New();
173 
179 
180  // Re-use any superclass signatures that we don't override.
185 
192  int IntersectWithLine(const double p1[3], const double p2[3], double tol, double& t, double x[3],
193  double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
194 
204  int IntersectWithLine(const double p1[3], const double p2[3], const double tol, vtkPoints* points,
205  vtkIdList* cellIds, vtkGenericCell* cell) override;
206 
216  void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
217  vtkIdType& cellId, int& subId, double& dist2) override
218  {
219  this->Superclass::FindClosestPoint(x, closestPoint, cell, cellId, subId, dist2);
220  }
221 
234  vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
235  vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
236 
244  vtkIdType FindCell(double x[3], double vtkNotUsed(tol2), vtkGenericCell* GenCell, int& subId,
245  double pcoords[3], double* weights) override;
246 
251  void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
252 
262  const double p1[3], const double p2[3], double tolerance, vtkIdList* cellsIds) override
263  {
264  this->Superclass::FindCellsAlongLine(p1, p2, tolerance, cellsIds);
265  }
266 
268 
271  void FreeSearchStructure() override;
272  void BuildLocator() override;
273  void ForceBuildLocator() override;
274  void GenerateRepresentation(int level, vtkPolyData* pd) override;
276 
277  VTK_DEPRECATED_IN_9_2_0("This method is deprecated because LazyEvaluation has been deprecated")
278  virtual void BuildLocatorIfNeeded() {}
279 
283  virtual vtkIdList* GetCells(int bucket);
284 
289  virtual int GetNumberOfBuckets(void);
290 
294  void ShallowCopy(vtkAbstractCellLocator* locator) override;
295 
296 protected:
298  ~vtkCellLocator() override;
299 
300  void BuildLocatorInternal() override;
301 
302  //------------------------------------------------------------------------------
304  {
305  public:
307 
308  inline int GetNumberOfNeighbors();
309 
310  inline void Reset();
311 
312  inline int* GetPoint(int i);
313 
314  inline int InsertNextPoint(int* x);
315 
316  protected:
318  };
319 
320  void GetOverlappingBuckets(vtkNeighborCells& buckets, const double x[3], double dist,
321  int prevMinLevel[3], int prevMaxLevel[3]);
322 
323  inline void GetBucketIndices(const double x[3], int ijk[3]);
324 
325  double Distance2ToBucket(const double x[3], int nei[3]);
326  double Distance2ToBounds(const double x[3], double bounds[6]);
327 
328  int NumberOfOctants; // number of octants in tree
329  double Bounds[6]; // bounding box root octant
330  double H[3]; // width of leaf octant in x-y-z directions
331  int NumberOfDivisions; // number of "leaf" octant sub-divisions
332  std::shared_ptr<std::vector<vtkSmartPointer<vtkIdList>>> TreeSharedPtr;
334 
335  void MarkParents(const vtkSmartPointer<vtkIdList>&, int, int, int, int, int);
336  int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType& idx);
338  int face, int numDivs, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
339  void ComputeOctantBounds(double octantBounds[6], int i, int j, int k);
340 
341 private:
342  vtkCellLocator(const vtkCellLocator&) = delete;
343  void operator=(const vtkCellLocator&) = delete;
344 };
345 
346 #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:296
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(void)
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.
int IntersectWithLine(const double p1[3], const double p2[3], const double tol, vtkPoints *points, vtkIdList *cellIds, vtkGenericCell *cell) override
Take the passed line segment and intersect it with the data set.
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)
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:143
a simple class to control print indentation
Definition: vtkIndent.h:119
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:155
represent and manipulate 3D points
Definition: vtkPoints.h:149
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:200
@ points
Definition: vtkX3D.h:452
@ level
Definition: vtkX3D.h:401
@ radius
Definition: vtkX3D.h:258
@ size
Definition: vtkX3D.h:259
@ offset
Definition: vtkX3D.h:444
#define VTK_DEPRECATED_IN_9_2_0(reason)
int vtkIdType
Definition: vtkType.h:332