VTK  9.1.0
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 =========================================================================*/
42 #ifndef vtkCellLocator_h
43 #define vtkCellLocator_h
44 
45 #include "vtkAbstractCellLocator.h"
46 #include "vtkCommonDataModelModule.h" // For export macro
47 
48 class vtkNeighborCells;
49 
50 class VTKCOMMONDATAMODEL_EXPORT vtkCellLocator : public vtkAbstractCellLocator
51 {
52 public:
54  void PrintSelf(ostream& os, vtkIndent indent) override;
55 
60  static vtkCellLocator* New();
61 
67 
68  // Re-use any superclass signatures that we don't override.
73 
81  int IntersectWithLine(const double a0[3], const double a1[3], double tol, double& t, double x[3],
82  double pcoords[3], int& subId, vtkIdType& cellId, vtkGenericCell* cell) override;
83 
96  void FindClosestPoint(const double x[3], double closestPoint[3], vtkGenericCell* cell,
97  vtkIdType& cellId, int& subId, double& dist2) override;
98 
117  vtkIdType FindClosestPointWithinRadius(double x[3], double radius, double closestPoint[3],
118  vtkGenericCell* cell, vtkIdType& cellId, int& subId, double& dist2, int& inside) override;
119 
123  virtual vtkIdList* GetCells(int bucket);
124 
129  virtual int GetNumberOfBuckets(void);
130 
137  double x[3], double tol2, vtkGenericCell* GenCell, double pcoords[3], double* weights) override;
138 
144  void FindCellsWithinBounds(double* bbox, vtkIdList* cells) override;
145 
155  const double p1[3], const double p2[3], double tolerance, vtkIdList* cells) override;
156 
158 
161  void FreeSearchStructure() override;
162  void BuildLocator() override;
163  virtual void BuildLocatorIfNeeded();
164  virtual void ForceBuildLocator();
165  virtual void BuildLocatorInternal();
166  void GenerateRepresentation(int level, vtkPolyData* pd) override;
168 
169 protected:
171  ~vtkCellLocator() override;
172 
173  void GetBucketNeighbors(int ijk[3], int ndivs, int level);
175  const double x[3], int ijk[3], double dist, int prevMinLevel[3], int prevMaxLevel[3]);
176 
179 
180  double Distance2ToBucket(const double x[3], int nei[3]);
181  double Distance2ToBounds(const double x[3], double bounds[6]);
182 
183  int NumberOfOctants; // number of octants in tree
184  double Bounds[6]; // bounding box root octant
185  int NumberOfParents; // number of parent octants
186  double H[3]; // width of leaf octant in x-y-z directions
187  int NumberOfDivisions; // number of "leaf" octant sub-divisions
188  vtkIdList** Tree; // octree
189 
190  void MarkParents(void*, int, int, int, int, int);
191  void GetChildren(int idx, int level, int children[8]);
192  int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType& idx);
194  int face, int numDivs, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
195 
196  vtkNeighborCells* Buckets;
197  unsigned char* CellHasBeenVisited;
198  unsigned char QueryNumber;
199 
200  void ComputeOctantBounds(int i, int j, int k);
201  double OctantBounds[6]; // the bounds of the current octant
202  int IsInOctantBounds(const double x[3], double tol = 0.0)
203  {
204  if (this->OctantBounds[0] - tol <= x[0] && x[0] <= this->OctantBounds[1] + tol &&
205  this->OctantBounds[2] - tol <= x[1] && x[1] <= this->OctantBounds[3] + tol &&
206  this->OctantBounds[4] - tol <= x[2] && x[2] <= this->OctantBounds[5] + tol)
207  {
208  return 1;
209  }
210  else
211  {
212  return 0;
213  }
214  }
215 
216 private:
217  vtkCellLocator(const vtkCellLocator&) = delete;
218  void operator=(const vtkCellLocator&) = delete;
219 };
220 
221 #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:181
octree-based spatial search object to quickly locate cells
void FreeSearchStructure() override
Satisfy vtkLocator abstract interface.
~vtkCellLocator() override
double Distance2ToBounds(const double x[3], double bounds[6])
void ClearCellHasBeenVisited()
int GenerateIndex(int offset, int numDivs, int i, int j, int k, vtkIdType &idx)
int GetNumberOfCellsPerBucket()
int IntersectWithLine(const double a0[3], const double a1[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.
vtkNeighborCells * Buckets
void GenerateRepresentation(int level, vtkPolyData *pd) override
Satisfy vtkLocator abstract interface.
unsigned char QueryNumber
virtual int GetNumberOfBuckets(void)
Return number of buckets available.
void FindCellsAlongLine(const double p1[3], const double p2[3], double tolerance, vtkIdList *cells) override
Given a finite line defined by the two points (p1,p2), return the list of unique cell ids in the buck...
virtual void BuildLocatorInternal()
Satisfy vtkLocator abstract interface.
void GetOverlappingBuckets(const double x[3], int ijk[3], double dist, int prevMinLevel[3], int prevMaxLevel[3])
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.
vtkIdType FindCell(double x[3], double tol2, vtkGenericCell *GenCell, double pcoords[3], double *weights) override
Find the cell containing a given point.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkIdList ** Tree
void GetChildren(int idx, int level, int children[8])
void MarkParents(void *, int, int, int, int, int)
double Distance2ToBucket(const double x[3], int nei[3])
static vtkCellLocator * New()
Construct with automatic computation of divisions, averaging 25 cells per bucket.
void ComputeOctantBounds(int i, int j, int k)
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.
virtual void ForceBuildLocator()
Satisfy vtkLocator abstract interface.
void ClearCellHasBeenVisited(vtkIdType id)
void GetBucketNeighbors(int ijk[3], int ndivs, int level)
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.
int IsInOctantBounds(const double x[3], double tol=0.0)
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.
virtual void BuildLocatorIfNeeded()
Satisfy vtkLocator abstract interface.
unsigned char * CellHasBeenVisited
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:31
a simple class to control print indentation
Definition: vtkIndent.h:34
represent and manipulate 3D points
Definition: vtkPoints.h:34
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:86
@ level
Definition: vtkX3D.h:401
@ radius
Definition: vtkX3D.h:258
@ offset
Definition: vtkX3D.h:444
int vtkIdType
Definition: vtkType.h:332