VTK  9.0.1
vtkStaticPointLocator2D.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkStaticPointLocator2D.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 =========================================================================*/
56 #ifndef vtkStaticPointLocator2D_h
57 #define vtkStaticPointLocator2D_h
58 
60 #include "vtkCommonDataModelModule.h" // For export macro
61 
62 class vtkIdList;
63 struct vtkBucketList2D;
64 
65 class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator2D : public vtkAbstractPointLocator
66 {
67 public:
72  static vtkStaticPointLocator2D* New();
73 
75 
79  void PrintSelf(ostream& os, vtkIndent indent) override;
81 
83 
88  vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
89  vtkGetMacro(NumberOfPointsPerBucket, int);
91 
93 
99  vtkSetVector2Macro(Divisions, int);
100  vtkGetVectorMacro(Divisions, int, 2);
102 
103  // Re-use any superclass signatures that we don't override.
108 
116  vtkIdType FindClosestPoint(const double x[3]) override;
117 
119 
127  vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
129  double radius, const double x[3], double inputDataLength, double& dist2);
131 
140  void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
141 
148  void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
149 
159  int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
160  double ptX[3], vtkIdType& ptId);
161 
163 
169  double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList* result);
171 
179  void MergePoints(double tol, vtkIdType* mergeMap);
180 
182 
186  void Initialize() override;
187  void FreeSearchStructure() override;
188  void BuildLocator() override;
190 
197  vtkIdType GetNumberOfPointsInBucket(vtkIdType bNum);
198 
204  void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
205 
207 
221  vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
222  vtkGetMacro(MaxNumberOfBuckets, vtkIdType);
224 
232  bool GetLargeIds() { return this->LargeIds; }
233 
235 
238  void GetBounds(double* bounds) override
239  {
240  bounds[0] = this->Bounds[0];
241  bounds[1] = this->Bounds[1];
242  bounds[2] = this->Bounds[2];
243  bounds[3] = this->Bounds[3];
244  }
246 
248 
252  virtual double* GetSpacing() { return this->H; }
253  virtual void GetSpacing(double spacing[3])
254  {
255  spacing[0] = this->H[0];
256  spacing[1] = this->H[1];
257  spacing[2] = 0.0;
258  }
260 
262 
267  void GetBucketIndices(const double* x, int ij[2]) const;
268  vtkIdType GetBucketIndex(const double* x) const;
270 
276  void GenerateRepresentation(int level, vtkPolyData* pd) override;
277 
278 protected:
280  ~vtkStaticPointLocator2D() override;
281 
282  int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
283  int Divisions[2]; // Number of sub-divisions in x-y directions
284  double H[2]; // Width of each bucket in x-y directions
285  vtkBucketList2D* Buckets; // Lists of point ids in each bucket
286  vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
287  bool LargeIds; // indicate whether integer ids are small or large
288 
289 private:
291  void operator=(const vtkStaticPointLocator2D&) = delete;
292 };
293 
294 #endif
virtual void BuildLocator()=0
Build the locator from the input dataset.
virtual void GetSpacing(double spacing[3])
Provide an accessor to the bucket spacing.
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...
#define VTK_INT_MAX
Definition: vtkType.h:155
virtual double * GetBounds()
Provide an accessor to the bounds.
virtual double * GetSpacing()
Provide an accessor to the bucket spacing.
int vtkIdType
Definition: vtkType.h:338
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:84
virtual void FreeSearchStructure()=0
Free the memory required for the spatial data structure.
a simple class to control print indentation
Definition: vtkIndent.h:33
void PrintSelf(ostream &os, vtkIndent indent) override
Standard type and print methods.
abstract class to quickly locate points in 3-space
list of point or cell ids
Definition: vtkIdList.h:30
virtual void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result)=0
Find all points within a specified radius R of position x.
void GetBounds(double *bounds) override
Provide an accessor to the bounds.
virtual void FindClosestNPoints(int N, const double x[3], vtkIdList *result)=0
Find the closest N points to a position.
#define VTK_ID_MAX
Definition: vtkType.h:342
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...
bool GetLargeIds()
Inform the user as to whether large ids are being used.
virtual void GenerateRepresentation(int level, vtkPolyData *pd)=0
Method to build a representation at a particular level.
quickly locate points in 2-space