VTK  9.6.20260617
vtkStaticPointLocator2DPrivate.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright 2011 Sandia Corporation
3// SPDX-License-Identifier: LicenseRef-BSD-3-Clause-Sandia-USGov
12
13#ifndef vtkStaticPointLocator2DPrivate_h
14#define vtkStaticPointLocator2DPrivate_h
15
16#include "vtkArrayDispatch.h"
18#include "vtkCellArray.h"
19#include "vtkDataArrayRange.h"
20#include "vtkDataSet.h"
21#include "vtkDoubleArray.h"
22#include "vtkMath.h"
23#include "vtkPoints.h"
25#include "vtkSMPTools.h"
27#include "vtkStructuredData.h"
28
29VTK_ABI_NAMESPACE_BEGIN
30
31#define Distance2BetweenPoints2D(p1, p2) \
32 ((p1[0] - p2[0]) * (p1[0] - p2[0]) + (p1[1] - p2[1]) * (p1[1] - p2[1]))
33
39inline bool IntersectsCircle(
40 const double min[2], const double max[2], const double center[2], double r2)
41{
42 double d2 = 0.0;
43 for (int i = 0; i < 2; ++i)
44 {
45 if (center[i] < min[i])
46 {
47 d2 += (center[i] - min[i]) * (center[i] - min[i]);
48 }
49 else if (center[i] > max[i])
50 {
51 d2 += (center[i] - max[i]) * (center[i] - max[i]);
52 }
53 }
54 return (d2 <= r2);
55}
56
62inline bool InsideCircle(
63 const double min[2], const double max[2], const double center[2], double r2)
64{
65 double dmin = 0.0, dmax = 0.0;
66 for (int i = 0; i < 2; ++i)
67 {
68 double a = (center[i] - min[i]) * (center[i] - min[i]);
69 double b = (center[i] - max[i]) * (center[i] - max[i]);
70 dmax += std::max(a, b);
71 if (min[i] <= center[i] && center[i] <= max[i])
72 {
73 dmin += std::min(a, b);
74 }
75 }
76 return (!(dmin <= r2 && r2 <= dmax));
77}
78
79//------------------------------------------------------------------------------
80// The following code supports threaded point locator construction. The locator
81// is assumed to be constructed once (i.e., it does not allow incremental point
82// insertion). The algorithm proceeds in three steps:
83// 1) All points are assigned a bucket index (combined i-j bucket location).
84// The index is computed in parallel. This requires a one time allocation of an
85// index array (which is also associated with the originating point ids).
86// 2) vtkSMPTools::Sort() is used to sort the index array. Note that the sort
87// carries along the point ids as well. This creates contiguous runs of points
88// all resident in the same bucket.
89// 3) The bucket offsets are updated to refer to the right entry location into
90// the sorted point ids array. This enables quick access, and an indirect count
91// of the number of points in each bucket.
92
93struct NeighborBuckets2D;
94
95//------------------------------------------------------------------------------
96// The bucketed points, including the sorted map. This is just a PIMPLd
97// wrapper around the classes that do the real work.
99{
101 vtkIdType NumPts; // the number of points to bucket
104
105 // These are internal data members used for performance reasons
107 int Divisions[3];
108 double Bounds[6];
109 double H[3];
110 double hX, hY, hX2, hY2;
111 double fX, fY, bX, bY;
113
114 // Used for accelerated performance for certain methods
115 double* FastPoints; // fast path for accessing points
116 double BinRadius; // circumradius of a single bin/bucket
117 int MaxLevel; // the maximum possible level searches can proceed
118
119 // Construction
121 {
122 this->Locator = loc;
123 this->NumPts = numPts;
124 this->NumBuckets = numBuckets;
125 this->BatchSize = 10000; // building the offset array
126 this->DataSet = loc->GetDataSet();
127
128 // Setup internal data members for more efficient processing. Remember this is
129 // a 2D locator so just processing (x,y) points.
130 double spacing[3], bounds[6];
131 loc->GetDivisions(this->Divisions);
132 loc->GetSpacing(spacing);
133 loc->GetBounds(bounds);
134 this->hX = this->H[0] = spacing[0];
135 this->hY = this->H[1] = spacing[1];
136 this->hX2 = this->hX / 2.0;
137 this->hY2 = this->hY / 2.0;
138 this->fX = 1.0 / spacing[0];
139 this->fY = 1.0 / spacing[1];
140 this->bX = this->Bounds[0] = bounds[0];
141 this->Bounds[1] = bounds[1];
142 this->bY = this->Bounds[2] = bounds[2];
143 this->Bounds[3] = bounds[3];
144 this->xD = this->Divisions[0];
145 this->yD = this->Divisions[1];
146 this->zD = 1;
147
148 this->FastPoints = nullptr;
149 this->BinRadius = sqrt(hX * hX + hY * hY) / 2.0;
150 this->MaxLevel = std::max({ this->xD, this->yD, this->zD });
151 }
152
153 // Virtuals for templated subclasses
154 virtual ~vtkBucketList2D() = default;
155 virtual void BuildLocator() = 0;
156
157 // place points in appropriate buckets
159 NeighborBuckets2D* buckets, const int ij[2], const int ndivs[2], int level);
160 void GenerateFace(int face, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
161 double Distance2ToBucket(const double x[3], const int nei[3]);
162 double Distance2ToBounds(const double x[3], const double bounds[6]);
163
164 //-----------------------------------------------------------------------------
165 // Inlined for performance. These function invocations must be called after
166 // BuildLocator() is invoked, otherwise the output is indeterminate.
167 void GetBucketIndices(const double* x, int ij[2]) const
168 {
169 // Compute point index. Make sure it lies within range of locator.
170 vtkIdType tmp0 = static_cast<vtkIdType>(((x[0] - bX) * fX));
171 vtkIdType tmp1 = static_cast<vtkIdType>(((x[1] - bY) * fY));
172
173 ij[0] = std::min(std::max<vtkIdType>(tmp0, 0), xD - 1);
174 ij[1] = std::min(std::max<vtkIdType>(tmp1, 0), yD - 1);
175 }
176
177 //-----------------------------------------------------------------------------
178 vtkIdType GetBucketIndex(const double* x) const
179 {
180 int ij[2];
181 this->GetBucketIndices(x, ij);
182 return ij[0] + ij[1] * xD;
183 }
184
185 //-----------------------------------------------------------------------------
186 // Return the center of the bucket/bin at (i,j). Note, returns a 2D point
187 // center[2].
188 void GetBucketCenter(int i, int j, double center[3])
189 {
190 center[0] = this->bX + this->hX2 + i * this->hX;
191 center[1] = this->bY + this->hY2 + j * this->hY;
192 center[2] = 0.0;
193 }
194
195 //-----------------------------------------------------------------------------
196 // Return the bounding box (min,max) of a specified bucket (i,j,k).
197 void GetBucketBounds(int i, int j, double min[3], double max[3])
198 {
199 min[0] = this->bX + i * this->hX;
200 min[1] = this->bY + j * this->hY;
201 min[2] = 0.0;
202 max[0] = min[0] + this->hX;
203 max[1] = min[1] + this->hY;
204 max[2] = 0.0;
205 }
206
207 //-----------------------------------------------------------------------------
208 // Determine whether a bin/bucket specified by i,j is completely contained
209 // inside the circle (center,r2). Return true if contained; false otherwise.
210 bool BucketInsideCircle(int i, int j, const double center[3], double r2)
211 {
212 double xMin = this->bX + i * this->hX;
213 double xMax = xMin + this->hX;
214 double yMin = this->bY + j * this->hY;
215 double yMax = yMin + this->hY;
216 double dx0 = center[0] - xMin, dx1 = center[0] - xMax;
217 double dy0 = center[1] - yMin, dy1 = center[1] - yMax;
218 double dmax2 = std::max(dx0 * dx0, dx1 * dx1) + std::max(dy0 * dy0, dy1 * dy1);
219 return dmax2 <= r2;
220 }
221}; // vtkBucketList2D
222
223//------------------------------------------------------------------------------
224// This templated class manages the creation of the static locator
225// structures. It also implements the operator() functors which are supplied
226// to vtkSMPTools for threaded processesing.
227template <typename TIds>
229{
230 // Okay the various ivars
231 vtkLocatorTuple<TIds>* Map; // the map to be sorted
232 TIds* Offsets; // offsets for each bucket into the map
233
234 // Construction
235 BucketList2D(vtkStaticPointLocator2D* loc, vtkIdType numPts, int numBuckets)
236 : vtkBucketList2D(loc, numPts, numBuckets)
237 {
238 // one extra to simplify traversal
239 this->Map = new vtkLocatorTuple<TIds>[numPts + 1];
240 this->Map[numPts].Bucket = numBuckets;
241 this->Offsets = new TIds[numBuckets + 1];
242 this->Offsets[numBuckets] = numPts;
243 }
244
245 // Release allocated memory
246 ~BucketList2D() override
247 {
248 delete[] this->Map;
249 delete[] this->Offsets;
250 }
251
252 // The number of point ids in a bucket is determined by computing the
253 // difference between the offsets into the sorted points array.
255 {
256 if (bucketNum < 0 || bucketNum >= this->NumBuckets)
257 {
258 return 0;
259 }
260 return (this->Offsets[bucketNum + 1] - this->Offsets[bucketNum]);
261 }
262
263 // Given a bucket number, return the point ids in that bucket.
265 {
266 return this->Map + this->Offsets[bucketNum];
267 }
268
269 // Given a bucket number, return the point ids in that bucket.
270 void GetIds(vtkIdType bucketNum, vtkIdList* bList)
271 {
272 const vtkLocatorTuple<TIds>* ids = this->GetIds(bucketNum);
273 vtkIdType numIds = this->GetNumberOfIds(bucketNum);
274 bList->SetNumberOfIds(numIds);
275 for (int i = 0; i < numIds; i++)
276 {
277 bList->SetId(i, ids[i].PtId);
278 }
279 }
280
281 // Templated implementations of the locator
282 vtkIdType FindClosestPoint(const double x[3]);
284 double radius, const double x[3], double inputDataLength, double& dist2);
285 void FindClosestNPoints(int N, const double x[3], vtkIdList* result);
286 double FindNPointsInAnnulus(int N, const double x[3], vtkDist2TupleArray& results,
287 double minDist2 = (-0.1), bool sort = true, vtkDoubleArray* petals = nullptr);
288 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result);
289 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
290 double ptX[3], vtkIdType& ptId);
291 double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList* result);
292
293 void MergePoints(double tol, vtkIdType* pointMap);
294 void GenerateRepresentation(int vtkNotUsed(level), vtkPolyData* pd);
295
296 // Internal methods
297 bool BucketIntersectsCircle(int i, int j, const double center[3], double R2);
299 NeighborBuckets2D* buckets, const double x[3], const int ij[2], double dist, int level);
300 void GetOverlappingBuckets(NeighborBuckets2D* buckets, const double x[3], double dist,
301 int prevMinLevel[2], int prevMaxLevel[2]);
302
303 // Implicit point representation, slower path
304 template <typename T>
306 {
309
311 : BList(blist)
312 , DataSet(ds)
313 {
314 }
315
317 {
318 double p[3];
319 vtkLocatorTuple<T>* t = this->BList->Map + ptId;
320 for (; ptId < end; ++ptId, ++t)
321 {
322 this->DataSet->GetPoint(ptId, p);
323 t->PtId = ptId;
324 t->Bucket = this->BList->GetBucketIndex(p);
325 } // for all points in this batch
326 }
327 };
328
329 template <typename T, typename TPointsArray>
331 {
333 TPointsArray* Points;
334
335 MapPointsArray(BucketList2D<T>* blist, TPointsArray* pts)
336 : BList(blist)
337 , Points(pts)
338 {
339 }
340
342 {
343 double p[3];
344 auto x = vtk::DataArrayTupleRange<3>(this->Points, ptId, end).begin();
345 vtkLocatorTuple<T>* t = this->BList->Map + ptId;
346 for (; ptId < end; ++ptId, ++x, ++t)
347 {
348 x->GetTuple(p);
349 t->PtId = ptId;
350 t->Bucket = this->BList->GetBucketIndex(p);
351 } // for all points in this batch
352 }
353 };
354
356 {
357 template <typename TPointsArray>
358 void operator()(TPointsArray* points, BucketList2D* blist)
359 {
360 MapPointsArray<TIds, TPointsArray> mapper(blist, points);
361 vtkSMPTools::For(0, blist->NumPts, mapper);
362 }
363 };
364
365 // A clever way to build offsets in parallel. Basically each thread builds
366 // offsets across a range of the sorted map. Recall that offsets are an
367 // integral value referring to the locations of the sorted points that
368 // reside in each bucket.
369 template <typename T>
371 {
375
377 : BList(blist)
378 {
379 this->NumPts = this->BList->NumPts;
380 this->NumBuckets = this->BList->NumBuckets;
381 }
382
383 // Traverse sorted points (i.e., tuples) and update bucket offsets.
384 void operator()(vtkIdType batch, vtkIdType batchEnd)
385 {
386 T* offsets = this->BList->Offsets;
387 const vtkLocatorTuple<T>* curPt = this->BList->Map + batch * this->BList->BatchSize;
388 const vtkLocatorTuple<T>* endBatchPt = this->BList->Map + batchEnd * this->BList->BatchSize;
389 const vtkLocatorTuple<T>* endPt = this->BList->Map + this->NumPts;
390 const vtkLocatorTuple<T>* prevPt;
391 endBatchPt = (endBatchPt > endPt ? endPt : endBatchPt);
392
393 // Special case at the very beginning of the mapped points array. If
394 // the first point is in bucket# N, then all buckets up and including
395 // N must refer to the first point.
396 if (curPt == this->BList->Map)
397 {
398 prevPt = this->BList->Map;
399 std::fill_n(offsets, curPt->Bucket + 1, 0); // point to the first points
400 } // at the very beginning of the map (sorted points array)
401
402 // We are entering this functor somewhere in the interior of the
403 // mapped points array. All we need to do is point to the entry
404 // position because we are interested only in prevPt->Bucket.
405 else
406 {
407 prevPt = curPt;
408 } // else in the middle of a batch
409
410 // Okay we have a starting point for a bucket run. Now we can begin
411 // filling in the offsets in this batch. A previous thread should
412 // have/will have completed the previous and subsequent runs outside
413 // of the [batch,batchEnd) range
414 for (curPt = prevPt; curPt < endBatchPt;)
415 {
416 for (; curPt->Bucket == prevPt->Bucket && curPt <= endBatchPt; ++curPt)
417 {
418 // advance
419 }
420 // Fill in any gaps in the offset array
421 std::fill_n(
422 offsets + prevPt->Bucket + 1, curPt->Bucket - prevPt->Bucket, curPt - this->BList->Map);
423 prevPt = curPt;
424 } // for all batches in this range
425 } // operator()
426 };
427
428 // Merge points that are pecisely coincident. Operates in parallel on
429 // locator buckets. Does not need to check neighbor buckets.
430 template <typename T>
432 {
436
438 : BList(blist)
439 , MergeMap(mergeMap)
440 {
441 this->DataSet = blist->DataSet;
442 }
443
444 void operator()(vtkIdType bucket, vtkIdType endBucket)
445 {
446 BucketList2D<T>* bList = this->BList;
447 vtkIdType* mergeMap = this->MergeMap;
448 int i, j;
449 const vtkLocatorTuple<TIds>* ids;
450 double p[3], p2[3];
451 vtkIdType ptId, ptId2, numIds;
452
453 for (; bucket < endBucket; ++bucket)
454 {
455 if ((numIds = bList->GetNumberOfIds(bucket)) > 0)
456 {
457 ids = bList->GetIds(bucket);
458 for (i = 0; i < numIds; i++)
459 {
460 ptId = ids[i].PtId;
461 if (mergeMap[ptId] < 0)
462 {
463 mergeMap[ptId] = ptId;
464 this->DataSet->GetPoint(ptId, p);
465 for (j = i + 1; j < numIds; j++)
466 {
467 ptId2 = ids[j].PtId;
468 if (mergeMap[ptId2] < 0)
469 {
470 this->DataSet->GetPoint(ptId2, p2);
471 if (p[0] == p2[0] && p[1] == p2[1])
472 {
473 mergeMap[ptId2] = ptId;
474 }
475 }
476 }
477 } // if point not yet visited
478 }
479 }
480 }
481 }
482 };
483
484 // Merge points that are coincident within a tolerance. Operates in
485 // parallel on points. Needs to check neighbor buckets which slows it down
486 // considerably. Note that merging is one direction: larger ids are merged
487 // to lower.
488 template <typename T>
490 {
494 double Tol;
495
497
498 MergeClose(BucketList2D<T>* blist, double tol, vtkIdType* mergeMap)
499 : BList(blist)
500 , MergeMap(mergeMap)
501 , Tol(tol)
502 {
503 this->DataSet = blist->DataSet;
504 }
505
506 // Just allocate a little bit of memory to get started.
508 {
509 vtkIdList*& pIds = this->PIds.Local();
510 pIds->Reserve(128); // allocate some memory
511 }
512
513 void operator()(vtkIdType ptId, vtkIdType endPtId)
514 {
515 BucketList2D<T>* bList = this->BList;
516 vtkIdType* mergeMap = this->MergeMap;
517 int i;
518 double p[3];
519 vtkIdType nearId, numIds;
520 vtkIdList*& nearby = this->PIds.Local();
521
522 for (; ptId < endPtId; ++ptId)
523 {
524 if (mergeMap[ptId] < 0)
525 {
526 mergeMap[ptId] = ptId;
527 this->DataSet->GetPoint(ptId, p);
528 bList->FindPointsWithinRadius(this->Tol, p, nearby);
529 if ((numIds = nearby->GetNumberOfIds()) > 0)
530 {
531 for (i = 0; i < numIds; i++)
532 {
533 nearId = nearby->GetId(i);
534 if (ptId < nearId && (mergeMap[nearId] < 0 || ptId < mergeMap[nearId]))
535 {
536 mergeMap[nearId] = ptId;
537 }
538 }
539 }
540 } // if point not yet processed
541 } // for all points in this batch
542 }
543
544 void Reduce() {}
545 };
546
547 // Build the map and other structures to support locator operations
548 void BuildLocator() override
549 {
550 // Place each point in a bucket
551 auto points = this->DataSet->GetPoints()->GetData();
554 points, worker, this))
555 {
556 worker(points, this);
557 }
558
559 // Provide accelerated access to points. Needed for Voronoi bin iterators.
561 {
562 this->FastPoints = static_cast<double*>(vtkDoubleArray::SafeDownCast(points)->GetPointer(0));
563 }
564
565 // Now gather the points into contiguous runs in buckets
566 //
567 vtkSMPTools::Sort(this->Map, this->Map + this->NumPts);
568
569 // Build the offsets into the Map. The offsets are the positions of
570 // each bucket into the sorted list. They mark the beginning of the
571 // list of points in each bucket. Amazingly, this can be done in
572 // parallel.
573 //
574 int numBatches = static_cast<int>(ceil(static_cast<double>(this->NumPts) / this->BatchSize));
575 MapOffsets<TIds> offMapper(this);
576 vtkSMPTools::For(0, numBatches, offMapper);
577 }
578}; // BucketList2D
579
580VTK_ABI_NAMESPACE_END
581#endif // vtkStaticPointLocator2DPrivate_h
582// VTK-HeaderTest-Exclude: vtkStaticPointLocator2DPrivate.h
RealT r2
Definition PyrC2Basis.h:20
ValueType * GetPointer(vtkIdType valueIdx)
Get the address of a particular data index.
object to represent cell connectivity
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
virtual double * GetPoint(vtkIdType ptId)=0
Get point coordinates with ptId such that: 0 <= ptId < NumberOfPoints.
dynamic, self-adjusting array of double
static vtkDoubleArray * SafeDownCast(vtkObjectBase *o)
list of point or cell ids
Definition vtkIdList.h:135
void SetNumberOfIds(vtkIdType number)
Specify the number of ids for this object to hold.
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition vtkIdList.h:185
vtkIdType GetId(vtkIdType i)
Return the id at location i.
Definition vtkIdList.h:190
void SetId(vtkIdType i, vtkIdType id)
Set the id at location i.
Definition vtkIdList.h:222
vtkTypeBool Reserve(vtkIdType size)
Reserve the id list to the requested number of ids and preserve data.
virtual vtkDataSet * GetDataSet()
Build the locator from the points/cells defining this dataset.
represent and manipulate 3D points
Definition vtkPoints.h:140
concrete dataset represents vertices, lines, polygons, and triangle strips
Thread local storage for VTK objects.
T *& Local()
Returns an object local to the current thread.
static void Sort(RandomAccessIterator begin, RandomAccessIterator end)
A convenience method for sorting data.
static void For(vtkIdType first, vtkIdType last, vtkIdType grain, Functor &f)
Execute a for operation in parallel.
quickly locate points in 2-space
void GetBounds(double *bounds) override
Provide an accessor to the bounds.
virtual double * GetSpacing()
Provide an accessor to the bucket spacing.
virtual int * GetDivisions()
Set the number of divisions in x-y directions.
VTK_ITER_INLINE auto DataArrayTupleRange(const ArrayTypePtr &array, TupleIdType start=-1, TupleIdType end=-1) -> typename detail::SelectTupleRange< ArrayTypePtr, TupleSize >::type
Generate an stl and for-range compatible range of tuple iterators from a vtkDataArray.
MapDataSet(BucketList2D< T > *blist, vtkDataSet *ds)
void operator()(vtkIdType ptId, vtkIdType end)
void operator()(vtkIdType batch, vtkIdType batchEnd)
void operator()(TPointsArray *points, BucketList2D *blist)
void operator()(vtkIdType ptId, vtkIdType end)
MapPointsArray(BucketList2D< T > *blist, TPointsArray *pts)
MergeClose(BucketList2D< T > *blist, double tol, vtkIdType *mergeMap)
vtkSMPThreadLocalObject< vtkIdList > PIds
void operator()(vtkIdType ptId, vtkIdType endPtId)
void operator()(vtkIdType bucket, vtkIdType endBucket)
MergePrecise(BucketList2D< T > *blist, vtkIdType *mergeMap)
void FindClosestNPoints(int N, const double x[3], vtkIdList *result)
double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList *result)
const vtkLocatorTuple< TIds > * GetIds(vtkIdType bucketNum)
BucketList2D(vtkStaticPointLocator2D *loc, vtkIdType numPts, int numBuckets)
void MergePoints(double tol, vtkIdType *pointMap)
void GetOverlappingBuckets(NeighborBuckets2D *buckets, const double x[3], const int ij[2], double dist, int level)
void GetOverlappingBuckets(NeighborBuckets2D *buckets, const double x[3], double dist, int prevMinLevel[2], int prevMaxLevel[2])
vtkLocatorTuple< TIds > * Map
void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result)
int IntersectWithLine(double a0[3], double a1[3], double tol, double &t, double lineX[3], double ptX[3], vtkIdType &ptId)
vtkIdType GetNumberOfIds(vtkIdType bucketNum)
bool BucketIntersectsCircle(int i, int j, const double center[3], double R2)
void GetIds(vtkIdType bucketNum, vtkIdList *bList)
vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double inputDataLength, double &dist2)
vtkIdType FindClosestPoint(const double x[3])
double FindNPointsInAnnulus(int N, const double x[3], vtkDist2TupleArray &results, double minDist2=(-0.1), bool sort=true, vtkDoubleArray *petals=nullptr)
void GenerateRepresentation(int level, vtkPolyData *pd)
Dispatch a single array against all array types mentioned in the ArrayList template parameter.
double Distance2ToBucket(const double x[3], const int nei[3])
void GetBucketBounds(int i, int j, double min[3], double max[3])
double Distance2ToBounds(const double x[3], const double bounds[6])
void GetBucketIndices(const double *x, int ij[2]) const
void GenerateFace(int face, int i, int j, int k, vtkPoints *pts, vtkCellArray *polys)
vtkIdType GetBucketIndex(const double *x) const
bool BucketInsideCircle(int i, int j, const double center[3], double r2)
virtual ~vtkBucketList2D()=default
vtkStaticPointLocator2D * Locator
virtual void BuildLocator()=0
void GetBucketCenter(int i, int j, double center[3])
void GetBucketNeighbors(NeighborBuckets2D *buckets, const int ij[2], const int ndivs[2], int level)
vtkBucketList2D(vtkStaticPointLocator2D *loc, vtkIdType numPts, int numBuckets)
Represent an array of vtkDist2Tuples.
STL-compatible iterable ranges that provide access to vtkDataArray elements.
bool InsideCircle(const double min[2], const double max[2], const double center[2], double r2)
Performant method to determine if a box if fully inside a circle.
bool IntersectsCircle(const double min[2], const double max[2], const double center[2], double r2)
Performant method to intersect a box with a circle.
int vtkIdType
Definition vtkType.h:363
#define max(a, b)