VTK  9.2.6
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=========================================================================*/
58
59#ifndef vtkStaticPointLocator2D_h
60#define vtkStaticPointLocator2D_h
61
63#include "vtkCommonDataModelModule.h" // For export macro
64
65class vtkIdList;
66struct vtkBucketList2D;
67
68class VTKCOMMONDATAMODEL_EXPORT vtkStaticPointLocator2D : public vtkAbstractPointLocator
69{
70public:
76
78
82 void PrintSelf(ostream& os, vtkIndent indent) override;
84
86
91 vtkSetClampMacro(NumberOfPointsPerBucket, int, 1, VTK_INT_MAX);
92 vtkGetMacro(NumberOfPointsPerBucket, int);
94
96
102 vtkSetVector2Macro(Divisions, int);
103 vtkGetVectorMacro(Divisions, int, 2);
105
106 // Re-use any superclass signatures that we don't override.
111
119 vtkIdType FindClosestPoint(const double x[3]) override;
120
122
130 vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double& dist2) override;
132 double radius, const double x[3], double inputDataLength, double& dist2);
134
143 void FindClosestNPoints(int N, const double x[3], vtkIdList* result) override;
144
151 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result) override;
152
162 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
163 double ptX[3], vtkIdType& ptId);
164
166
172 double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList* result);
174
182 void MergePoints(double tol, vtkIdType* mergeMap);
183
185
189 void Initialize() override;
190 void FreeSearchStructure() override;
191 void BuildLocator() override;
192 void ForceBuildLocator() override;
194
202
208 void GetBucketIds(vtkIdType bNum, vtkIdList* bList);
209
211
225 vtkSetClampMacro(MaxNumberOfBuckets, vtkIdType, 1000, VTK_ID_MAX);
228
236 bool GetLargeIds() { return this->LargeIds; }
237
239
242 void GetBounds(double* bounds) override
243 {
244 bounds[0] = this->Bounds[0];
245 bounds[1] = this->Bounds[1];
246 bounds[2] = this->Bounds[2];
247 bounds[3] = this->Bounds[3];
248 }
249
250
252
256 virtual double* GetSpacing() { return this->H; }
257 virtual void GetSpacing(double spacing[3])
258 {
259 spacing[0] = this->H[0];
260 spacing[1] = this->H[1];
261 spacing[2] = 0.0;
262 }
263
264
266
271 void GetBucketIndices(const double* x, int ij[2]) const;
272 vtkIdType GetBucketIndex(const double* x) const;
274
280 void GenerateRepresentation(int level, vtkPolyData* pd) override;
281
282protected:
285
286 void BuildLocatorInternal() override;
287
288 int NumberOfPointsPerBucket; // Used with AutomaticOn to control subdivide
289 int Divisions[2]; // Number of sub-divisions in x-y directions
290 double H[2]; // Width of each bucket in x-y directions
291 vtkBucketList2D* Buckets; // Lists of point ids in each bucket
292 vtkIdType MaxNumberOfBuckets; // Maximum number of buckets in locator
293 bool LargeIds; // indicate whether integer ids are small or large
294
295private:
297 void operator=(const vtkStaticPointLocator2D&) = delete;
298};
299
300#endif
virtual void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result)=0
Find all points within a specified radius R of position x.
virtual vtkIdType FindClosestPoint(const double x[3])=0
Given a position x, return the id of the point closest to it.
virtual double * GetBounds()
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.
list of point or cell ids
Definition vtkIdList.h:34
a simple class to control print indentation
Definition vtkIndent.h:40
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:91
vtkIdType FindClosestPoint(const double x[3]) override
Given a position x, return the id of the point closest to it.
void BuildLocator() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double &dist2) override
Given a position x and a radius r, return the id of the point closest to the point within that radius...
void PrintSelf(ostream &os, vtkIndent indent) override
Standard type and print methods.
void FindClosestNPoints(int N, const double x[3], vtkIdList *result) override
Find the closest N points to a position.
double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList *result)
Special method for 2D operations (e.g., vtkVoronoi2D).
~vtkStaticPointLocator2D() override
bool GetLargeIds()
Inform the user as to whether large ids are being used.
void BuildLocatorInternal() override
This function is not pure virtual to maintain backwards compatibility.
void GetBucketIndices(const double *x, int ij[2]) const
Given a point x[3], return the locator index (i,j) which contains the point.
virtual vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double inputDataLength, double &dist2)
Given a position x and a radius r, return the id of the point closest to the point within that radius...
vtkIdType GetBucketIndex(const double *x) const
Given a point x[3], return the locator index (i,j) which contains the point.
int IntersectWithLine(double a0[3], double a1[3], double tol, double &t, double lineX[3], double ptX[3], vtkIdType &ptId)
Intersect the points contained in the locator with the line defined by (a0,a1).
void Initialize() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result) override
Find all points within a specified radius R of position x.
void GetBounds(double *bounds) override
Provide an accessor to the bounds.
void GenerateRepresentation(int level, vtkPolyData *pd) override
Populate a polydata with the faces of the bins that potentially contain cells.
void GetBucketIds(vtkIdType bNum, vtkIdList *bList)
Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return a list of point ids...
virtual double * GetSpacing()
Provide an accessor to the bucket spacing.
static vtkStaticPointLocator2D * New()
Construct with automatic computation of divisions, averaging 5 points per bucket.
void ForceBuildLocator() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void FreeSearchStructure() override
See vtkLocator and vtkAbstractPointLocator interface documentation.
void MergePoints(double tol, vtkIdType *mergeMap)
Merge points in the locator given a tolerance.
virtual void GetSpacing(double spacing[3])
Provide an accessor to the bucket spacing.
vtkIdType GetNumberOfPointsInBucket(vtkIdType bNum)
Given a bucket number bNum between 0 <= bNum < this->GetNumberOfBuckets(), return the number of point...
int vtkIdType
Definition vtkType.h:332
#define VTK_ID_MAX
Definition vtkType.h:336
#define VTK_INT_MAX
Definition vtkType.h:155