VTK  9.2.6
vtkPointLocator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPointLocator.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 =========================================================================*/
49 #ifndef vtkPointLocator_h
50 #define vtkPointLocator_h
51 
52 #include "vtkCommonDataModelModule.h" // For export macro
54 
55 class vtkCellArray;
56 class vtkIdList;
57 class vtkNeighborPoints;
58 class vtkPoints;
59 
60 class VTKCOMMONDATAMODEL_EXPORT vtkPointLocator : public vtkIncrementalPointLocator
61 {
62 public:
67  static vtkPointLocator* New();
68 
70 
74  void PrintSelf(ostream& os, vtkIndent indent) override;
76 
78 
81  vtkSetVector3Macro(Divisions, int);
82  vtkGetVectorMacro(Divisions, int, 3);
84 
86 
89  vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
90  vtkGetMacro(NumberOfPointsPerBucket, int);
92 
93  // Re-use any superclass signatures that we don't override.
95 
102  vtkIdType FindClosestPoint(const double x[3]) override;
103 
105 
112  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
114  double radius, const double x[3], double inputDataLength, double& dist2);
116 
123  int InitPointInsertion(vtkPoints* newPts, const double bounds[6]) override;
124 
131  int InitPointInsertion(vtkPoints* newPts, const double bounds[6], vtkIdType estNumPts) override;
132 
142  void InsertPoint(vtkIdType ptId, const double x[3]) override;
143 
154  vtkIdType InsertNextPoint(const double x[3]) override;
155 
157 
162  vtkIdType IsInsertedPoint(double x, double y, double z) override
163  {
164  double xyz[3];
165  xyz[0] = x;
166  xyz[1] = y;
167  xyz[2] = z;
168  return this->IsInsertedPoint(xyz);
169  };
170  vtkIdType IsInsertedPoint(const double x[3]) override;
172 
182  int InsertUniquePoint(const double x[3], vtkIdType& ptId) override;
183 
191  vtkIdType FindClosestInsertedPoint(const double x[3]) override;
192 
201  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
202 
204 
211  virtual void FindDistributedPoints(int N, const double x[3], vtkIdList* result, int M);
212  virtual void FindDistributedPoints(int N, double x, double y, double z, vtkIdList* result, int M);
214 
221  void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
222 
229  virtual vtkIdList* GetPointsInBucket(const double x[3], int ijk[3]);
230 
232 
235  vtkGetObjectMacro(Points, vtkPoints);
237 
239 
243  void Initialize() override;
244  void FreeSearchStructure() override;
245  void BuildLocator() override;
246  void ForceBuildLocator() override;
247  void GenerateRepresentation(int level, vtkPolyData* pd) override;
249 
250 protected:
251  vtkPointLocator();
252  ~vtkPointLocator() override;
253 
254  void BuildLocatorInternal() override;
255 
256  // place points in appropriate buckets
257  void GetBucketNeighbors(
258  vtkNeighborPoints* buckets, const int ijk[3], const int ndivs[3], int level);
259  void GetOverlappingBuckets(
260  vtkNeighborPoints* buckets, const double x[3], const int ijk[3], double dist, int level);
261  void GetOverlappingBuckets(vtkNeighborPoints* buckets, const double x[3], double dist,
262  int prevMinLevel[3], int prevMaxLevel[3]);
263  void GenerateFace(int face, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
264  double Distance2ToBucket(const double x[3], const int nei[3]);
265  double Distance2ToBounds(const double x[3], const double bounds[6]);
266 
267  vtkPoints* Points; // Used for merging points
268  int Divisions[3]; // Number of sub-divisions in x-y-z directions
269  int NumberOfPointsPerBucket; // Used with previous boolean to control subdivide
270  vtkIdList** HashTable; // lists of point ids in buckets
271  double H[3]; // width of each bucket in x-y-z directions
272 
276 
277  // These are inlined methods and data members for performance reasons
278  double HX, HY, HZ;
279  double FX, FY, FZ, BX, BY, BZ;
280  vtkIdType XD, YD, ZD, SliceSize;
281 
282  void GetBucketIndices(const double* x, int ijk[3]) const
283  {
284  // Compute point index. Make sure it lies within range of locator.
285  vtkIdType tmp0 = static_cast<vtkIdType>(((x[0] - this->BX) * this->FX));
286  vtkIdType tmp1 = static_cast<vtkIdType>(((x[1] - this->BY) * this->FY));
287  vtkIdType tmp2 = static_cast<vtkIdType>(((x[2] - this->BZ) * this->FZ));
288 
289  ijk[0] = tmp0 < 0 ? 0 : (tmp0 >= this->XD ? this->XD - 1 : tmp0);
290  ijk[1] = tmp1 < 0 ? 0 : (tmp1 >= this->YD ? this->YD - 1 : tmp1);
291  ijk[2] = tmp2 < 0 ? 0 : (tmp2 >= this->ZD ? this->ZD - 1 : tmp2);
292  }
293 
294  vtkIdType GetBucketIndex(const double* x) const
295  {
296  int ijk[3];
297  this->GetBucketIndices(x, ijk);
298  return ijk[0] + ijk[1] * this->XD + ijk[2] * this->SliceSize;
299  }
300 
301  void ComputePerformanceFactors();
302 
303 private:
304  vtkPointLocator(const vtkPointLocator&) = delete;
305  void operator=(const vtkPointLocator&) = delete;
306 };
307 
308 #endif
vtkIdType GetBucketIndex(const double *x) const
vtkIdType IsInsertedPoint(double x, double y, double z) override
Determine whether point given by x[3] has been inserted into points list.
void GetBucketIndices(const double *x, int ijk[3]) const
virtual void BuildLocator()=0
Build the locator from the input dataset.
virtual vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double &dist2)=0
Given a position x and a radius r, return the id of the point closest to the point in that radius...
quickly locate points in 3-space
virtual int InsertUniquePoint(const double x[3], vtkIdType &ptId)=0
Insert a point unless there has been a duplicate in the search structure.
vtkIdType InsertionPointId
virtual vtkIdType IsInsertedPoint(double x, double y, double z)=0
Determine whether or not a given point has been inserted.
#define VTK_INT_MAX
Definition: vtkType.h:155
vtkIdList ** HashTable
Abstract class in support of both point location and point insertion.
int vtkIdType
Definition: vtkType.h:332
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:90
virtual void FreeSearchStructure()=0
Free the memory required for the spatial data structure.
virtual int InitPointInsertion(vtkPoints *newPts, const double bounds[6])=0
Initialize the point insertion process.
a simple class to control print indentation
Definition: vtkIndent.h:39
virtual void BuildLocatorInternal()
This function is not pure virtual to maintain backwards compatibility.
Definition: vtkLocator.h:203
list of point or cell ids
Definition: vtkIdList.h:33
vtkPoints * Points
virtual void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result)=0
Find all points within a specified radius R of position x.
virtual void FindClosestNPoints(int N, const double x[3], vtkIdList *result)=0
Find the closest N points to a position.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard type and print methods.
object to represent cell connectivity
Definition: vtkCellArray.h:186
virtual vtkIdType InsertNextPoint(const double x[3])=0
Insert a given point and return the point index.
virtual vtkIdType FindClosestInsertedPoint(const double x[3])=0
Given a point x assumed to be covered by the search structure, return the index of the closest point ...
virtual vtkIdType FindClosestPoint(const double x[3])=0
Given a position x, return the id of the point closest to it.
virtual void Initialize()
Initialize locator.
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
virtual void InsertPoint(vtkIdType ptId, const double x[3])=0
Insert a given point with a specified point index ptId.
virtual void GenerateRepresentation(int level, vtkPolyData *pd)=0
Method to build a representation at a particular level.
represent and manipulate 3D points
Definition: vtkPoints.h:39
virtual void ForceBuildLocator()
Build the locator from the input dataset (even if UseExistingSearchStructure is on).
Definition: vtkLocator.h:167