VTK  9.6.20260307
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}; // vtkBucketList2D
207
208//------------------------------------------------------------------------------
209// This templated class manages the creation of the static locator
210// structures. It also implements the operator() functors which are supplied
211// to vtkSMPTools for threaded processesing.
212template <typename TIds>
214{
215 // Okay the various ivars
216 vtkLocatorTuple<TIds>* Map; // the map to be sorted
217 TIds* Offsets; // offsets for each bucket into the map
218
219 // Construction
220 BucketList2D(vtkStaticPointLocator2D* loc, vtkIdType numPts, int numBuckets)
221 : vtkBucketList2D(loc, numPts, numBuckets)
222 {
223 // one extra to simplify traversal
224 this->Map = new vtkLocatorTuple<TIds>[numPts + 1];
225 this->Map[numPts].Bucket = numBuckets;
226 this->Offsets = new TIds[numBuckets + 1];
227 this->Offsets[numBuckets] = numPts;
228 }
229
230 // Release allocated memory
231 ~BucketList2D() override
232 {
233 delete[] this->Map;
234 delete[] this->Offsets;
235 }
236
237 // The number of point ids in a bucket is determined by computing the
238 // difference between the offsets into the sorted points array.
240 {
241 if (bucketNum < 0 || bucketNum >= this->NumBuckets)
242 {
243 return 0;
244 }
245 return (this->Offsets[bucketNum + 1] - this->Offsets[bucketNum]);
246 }
247
248 // Given a bucket number, return the point ids in that bucket.
250 {
251 return this->Map + this->Offsets[bucketNum];
252 }
253
254 // Given a bucket number, return the point ids in that bucket.
255 void GetIds(vtkIdType bucketNum, vtkIdList* bList)
256 {
257 const vtkLocatorTuple<TIds>* ids = this->GetIds(bucketNum);
258 vtkIdType numIds = this->GetNumberOfIds(bucketNum);
259 bList->SetNumberOfIds(numIds);
260 for (int i = 0; i < numIds; i++)
261 {
262 bList->SetId(i, ids[i].PtId);
263 }
264 }
265
266 // Templated implementations of the locator
267 vtkIdType FindClosestPoint(const double x[3]);
269 double radius, const double x[3], double inputDataLength, double& dist2);
270 void FindClosestNPoints(int N, const double x[3], vtkIdList* result);
271 double FindNPointsInAnnulus(int N, const double x[3], vtkDist2TupleArray& results,
272 double minDist2 = (-0.1), bool sort = true, vtkDoubleArray* petals = nullptr);
273 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result);
274 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
275 double ptX[3], vtkIdType& ptId);
276 double FindCloseNBoundedPoints(int N, const double x[3], vtkIdList* result);
277
278 void MergePoints(double tol, vtkIdType* pointMap);
279 void GenerateRepresentation(int vtkNotUsed(level), vtkPolyData* pd);
280
281 // Internal methods
282 bool BucketIntersectsCircle(int i, int j, const double center[3], double R2);
284 NeighborBuckets2D* buckets, const double x[3], const int ij[2], double dist, int level);
285 void GetOverlappingBuckets(NeighborBuckets2D* buckets, const double x[3], double dist,
286 int prevMinLevel[2], int prevMaxLevel[2]);
287
288 // Implicit point representation, slower path
289 template <typename T>
291 {
294
296 : BList(blist)
297 , DataSet(ds)
298 {
299 }
300
302 {
303 double p[3];
304 vtkLocatorTuple<T>* t = this->BList->Map + ptId;
305 for (; ptId < end; ++ptId, ++t)
306 {
307 this->DataSet->GetPoint(ptId, p);
308 t->PtId = ptId;
309 t->Bucket = this->BList->GetBucketIndex(p);
310 } // for all points in this batch
311 }
312 };
313
314 template <typename T, typename TPointsArray>
316 {
318 TPointsArray* Points;
319
320 MapPointsArray(BucketList2D<T>* blist, TPointsArray* pts)
321 : BList(blist)
322 , Points(pts)
323 {
324 }
325
327 {
328 double p[3];
329 auto x = vtk::DataArrayTupleRange<3>(this->Points, ptId, end).begin();
330 vtkLocatorTuple<T>* t = this->BList->Map + ptId;
331 for (; ptId < end; ++ptId, ++x, ++t)
332 {
333 x->GetTuple(p);
334 t->PtId = ptId;
335 t->Bucket = this->BList->GetBucketIndex(p);
336 } // for all points in this batch
337 }
338 };
339
341 {
342 template <typename TPointsArray>
343 void operator()(TPointsArray* points, BucketList2D* blist)
344 {
345 MapPointsArray<TIds, TPointsArray> mapper(blist, points);
346 vtkSMPTools::For(0, blist->NumPts, mapper);
347 }
348 };
349
350 // A clever way to build offsets in parallel. Basically each thread builds
351 // offsets across a range of the sorted map. Recall that offsets are an
352 // integral value referring to the locations of the sorted points that
353 // reside in each bucket.
354 template <typename T>
356 {
360
362 : BList(blist)
363 {
364 this->NumPts = this->BList->NumPts;
365 this->NumBuckets = this->BList->NumBuckets;
366 }
367
368 // Traverse sorted points (i.e., tuples) and update bucket offsets.
369 void operator()(vtkIdType batch, vtkIdType batchEnd)
370 {
371 T* offsets = this->BList->Offsets;
372 const vtkLocatorTuple<T>* curPt = this->BList->Map + batch * this->BList->BatchSize;
373 const vtkLocatorTuple<T>* endBatchPt = this->BList->Map + batchEnd * this->BList->BatchSize;
374 const vtkLocatorTuple<T>* endPt = this->BList->Map + this->NumPts;
375 const vtkLocatorTuple<T>* prevPt;
376 endBatchPt = (endBatchPt > endPt ? endPt : endBatchPt);
377
378 // Special case at the very beginning of the mapped points array. If
379 // the first point is in bucket# N, then all buckets up and including
380 // N must refer to the first point.
381 if (curPt == this->BList->Map)
382 {
383 prevPt = this->BList->Map;
384 std::fill_n(offsets, curPt->Bucket + 1, 0); // point to the first points
385 } // at the very beginning of the map (sorted points array)
386
387 // We are entering this functor somewhere in the interior of the
388 // mapped points array. All we need to do is point to the entry
389 // position because we are interested only in prevPt->Bucket.
390 else
391 {
392 prevPt = curPt;
393 } // else in the middle of a batch
394
395 // Okay we have a starting point for a bucket run. Now we can begin
396 // filling in the offsets in this batch. A previous thread should
397 // have/will have completed the previous and subsequent runs outside
398 // of the [batch,batchEnd) range
399 for (curPt = prevPt; curPt < endBatchPt;)
400 {
401 for (; curPt->Bucket == prevPt->Bucket && curPt <= endBatchPt; ++curPt)
402 {
403 // advance
404 }
405 // Fill in any gaps in the offset array
406 std::fill_n(
407 offsets + prevPt->Bucket + 1, curPt->Bucket - prevPt->Bucket, curPt - this->BList->Map);
408 prevPt = curPt;
409 } // for all batches in this range
410 } // operator()
411 };
412
413 // Merge points that are pecisely coincident. Operates in parallel on
414 // locator buckets. Does not need to check neighbor buckets.
415 template <typename T>
417 {
421
423 : BList(blist)
424 , MergeMap(mergeMap)
425 {
426 this->DataSet = blist->DataSet;
427 }
428
429 void operator()(vtkIdType bucket, vtkIdType endBucket)
430 {
431 BucketList2D<T>* bList = this->BList;
432 vtkIdType* mergeMap = this->MergeMap;
433 int i, j;
434 const vtkLocatorTuple<TIds>* ids;
435 double p[3], p2[3];
436 vtkIdType ptId, ptId2, numIds;
437
438 for (; bucket < endBucket; ++bucket)
439 {
440 if ((numIds = bList->GetNumberOfIds(bucket)) > 0)
441 {
442 ids = bList->GetIds(bucket);
443 for (i = 0; i < numIds; i++)
444 {
445 ptId = ids[i].PtId;
446 if (mergeMap[ptId] < 0)
447 {
448 mergeMap[ptId] = ptId;
449 this->DataSet->GetPoint(ptId, p);
450 for (j = i + 1; j < numIds; j++)
451 {
452 ptId2 = ids[j].PtId;
453 if (mergeMap[ptId2] < 0)
454 {
455 this->DataSet->GetPoint(ptId2, p2);
456 if (p[0] == p2[0] && p[1] == p2[1])
457 {
458 mergeMap[ptId2] = ptId;
459 }
460 }
461 }
462 } // if point not yet visited
463 }
464 }
465 }
466 }
467 };
468
469 // Merge points that are coincident within a tolerance. Operates in
470 // parallel on points. Needs to check neighbor buckets which slows it down
471 // considerably. Note that merging is one direction: larger ids are merged
472 // to lower.
473 template <typename T>
475 {
479 double Tol;
480
482
483 MergeClose(BucketList2D<T>* blist, double tol, vtkIdType* mergeMap)
484 : BList(blist)
485 , MergeMap(mergeMap)
486 , Tol(tol)
487 {
488 this->DataSet = blist->DataSet;
489 }
490
491 // Just allocate a little bit of memory to get started.
493 {
494 vtkIdList*& pIds = this->PIds.Local();
495 pIds->Allocate(128); // allocate some memory
496 }
497
498 void operator()(vtkIdType ptId, vtkIdType endPtId)
499 {
500 BucketList2D<T>* bList = this->BList;
501 vtkIdType* mergeMap = this->MergeMap;
502 int i;
503 double p[3];
504 vtkIdType nearId, numIds;
505 vtkIdList*& nearby = this->PIds.Local();
506
507 for (; ptId < endPtId; ++ptId)
508 {
509 if (mergeMap[ptId] < 0)
510 {
511 mergeMap[ptId] = ptId;
512 this->DataSet->GetPoint(ptId, p);
513 bList->FindPointsWithinRadius(this->Tol, p, nearby);
514 if ((numIds = nearby->GetNumberOfIds()) > 0)
515 {
516 for (i = 0; i < numIds; i++)
517 {
518 nearId = nearby->GetId(i);
519 if (ptId < nearId && (mergeMap[nearId] < 0 || ptId < mergeMap[nearId]))
520 {
521 mergeMap[nearId] = ptId;
522 }
523 }
524 }
525 } // if point not yet processed
526 } // for all points in this batch
527 }
528
529 void Reduce() {}
530 };
531
532 // Build the map and other structures to support locator operations
533 void BuildLocator() override
534 {
535 // Place each point in a bucket
536 auto points = this->DataSet->GetPoints()->GetData();
539 points, worker, this))
540 {
541 worker(points, this);
542 }
543
544 // Provide accelerated access to points. Needed for Voronoi bin iterators.
546 {
547 this->FastPoints = static_cast<double*>(vtkDoubleArray::SafeDownCast(points)->GetPointer(0));
548 }
549
550 // Now gather the points into contiguous runs in buckets
551 //
552 vtkSMPTools::Sort(this->Map, this->Map + this->NumPts);
553
554 // Build the offsets into the Map. The offsets are the positions of
555 // each bucket into the sorted list. They mark the beginning of the
556 // list of points in each bucket. Amazingly, this can be done in
557 // parallel.
558 //
559 int numBatches = static_cast<int>(ceil(static_cast<double>(this->NumPts) / this->BatchSize));
560 MapOffsets<TIds> offMapper(this);
561 vtkSMPTools::For(0, numBatches, offMapper);
562 }
563}; // BucketList2D
564
565VTK_ABI_NAMESPACE_END
566#endif // vtkStaticPointLocator2DPrivate_h
567// 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:133
void SetNumberOfIds(vtkIdType number)
Specify the number of ids for this object to hold.
int Allocate(vtkIdType sz, int strategy=0)
Allocate a capacity for sz ids in the list and set the number of stored ids in the list to 0.
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition vtkIdList.h:159
void SetId(vtkIdType i, vtkIdType vtkid)
Set the id at location i.
Definition vtkIdList.h:188
vtkIdType GetId(vtkIdType i)
Return the id at location i.
Definition vtkIdList.h:164
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
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)