VTK  9.5.20251210
vtkStaticPointLocatorPrivate.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 vtkStaticPointLocatorPrivate_h
14#define vtkStaticPointLocatorPrivate_h
15
16#include "vtkCellArray.h"
17#include "vtkDataSet.h"
18#include "vtkDoubleArray.h"
19#include "vtkMath.h"
20#include "vtkPointSet.h"
21#include "vtkPoints.h"
22#include "vtkPolyData.h"
23#include "vtkSMPThreadLocal.h"
25#include "vtkSMPTools.h"
27#include "vtkStructuredData.h"
28
29VTK_ABI_NAMESPACE_BEGIN
30
31//------------------------------------------------------------------------------
32// The following code supports threaded point locator construction. The locator
33// is assumed to be constructed once (i.e., it does not allow incremental point
34// insertion). The algorithm proceeds in three steps:
35// 1) All points are assigned a bucket index (combined i-j-k bucket location).
36// The index is computed in parallel. This requires a one time allocation of an
37// index array (which is also associated with the originating point ids).
38// 2) vtkSMPTools::Sort() is used to sort the index array. Note that the sort
39// carries along the point ids as well. This creates contiguous runs of points
40// all resident in the same bucket.
41// 3) The bucket offsets are updated to refer to the right entry location into
42// the sorted point ids array. This enables quick access, and an indirect count
43// of the number of points in each bucket.
44
45// Forward declaration of lists of neighboring buckets.
46struct NeighborBuckets;
47
48//------------------------------------------------------------------------------
49// The bucketed points, including the sorted map. This is just a PIMPLd
50// wrapper around the classes that do the real work.
52{
54 vtkIdType NumPts; // the number of points to bucket
57
58 // These are internal data members used for performance reasons
60 int Divisions[3];
61 double Bounds[6];
62 double H[3];
63 double hX, hY, hZ, hX2, hY2, hZ2;
64 double fX, fY, fZ, bX, bY, bZ;
66
67 // Used for accelerated performance for certain methods
68 double* FastPoints; // fast path for accessing points
69 double BinRadius; // circumradius of a single bin/bucket
70 int MaxLevel; // the maximum possible level searches can proceed
71
72 // Construction
73 vtkBucketList(vtkStaticPointLocator* loc, vtkIdType numPts, int numBuckets)
74 {
75 this->Locator = loc;
76 this->NumPts = numPts;
77 this->NumBuckets = numBuckets;
78 this->BatchSize = 10000; // building the offset array
79 this->DataSet = loc->GetDataSet();
80 loc->GetDivisions(this->Divisions);
81
82 // Setup internal data members for more efficient processing.
83 double spacing[3], bounds[6];
84 loc->GetDivisions(this->Divisions);
85 loc->GetSpacing(spacing);
86 loc->GetBounds(bounds);
87 this->hX = this->H[0] = spacing[0];
88 this->hY = this->H[1] = spacing[1];
89 this->hZ = this->H[2] = spacing[2];
90 this->hX2 = this->hX / 2.0;
91 this->hY2 = this->hY / 2.0;
92 this->hZ2 = this->hZ / 2.0;
93 this->fX = 1.0 / spacing[0];
94 this->fY = 1.0 / spacing[1];
95 this->fZ = 1.0 / spacing[2];
96 this->bX = this->Bounds[0] = bounds[0];
97 this->Bounds[1] = bounds[1];
98 this->bY = this->Bounds[2] = bounds[2];
99 this->Bounds[3] = bounds[3];
100 this->bZ = this->Bounds[4] = bounds[4];
101 this->Bounds[5] = bounds[5];
102 this->xD = this->Divisions[0];
103 this->yD = this->Divisions[1];
104 this->zD = this->Divisions[2];
105 this->xyD = this->Divisions[0] * this->Divisions[1];
106
107 this->FastPoints = nullptr;
108 this->BinRadius = sqrt(hX * hX + hY * hY + hZ * hZ) / 2.0;
109 this->MaxLevel = std::max({ this->xD, this->yD, this->zD });
110 }
111
112 // Virtuals for templated subclasses
113 virtual ~vtkBucketList() = default;
114 virtual void BuildLocator() = 0;
115
116 // place points in appropriate buckets
118 NeighborBuckets* buckets, const int ijk[3], const int ndivs[3], int level);
119 void GenerateFace(int face, int i, int j, int k, vtkPoints* pts, vtkCellArray* polys);
120 double Distance2ToBucket(const double x[3], const int nei[3]);
121 double Distance2ToBounds(const double x[3], const double bounds[6]);
122
123 //-----------------------------------------------------------------------------
124 // Inlined for performance. These function invocations must be called after
125 // BuildLocator() is invoked, otherwise the output is indeterminate.
126 void GetBucketIndices(const double* x, int ijk[3]) const
127 {
128 // Compute point index. Make sure it lies within range of locator.
129 vtkIdType tmp0 = static_cast<vtkIdType>(((x[0] - bX) * fX));
130 vtkIdType tmp1 = static_cast<vtkIdType>(((x[1] - bY) * fY));
131 vtkIdType tmp2 = static_cast<vtkIdType>(((x[2] - bZ) * fZ));
132
133 ijk[0] = tmp0 < 0 ? 0 : std::min(xD - 1, tmp0);
134 ijk[1] = tmp1 < 0 ? 0 : std::min(yD - 1, tmp1);
135 ijk[2] = tmp2 < 0 ? 0 : std::min(zD - 1, tmp2);
136 }
137
138 //-----------------------------------------------------------------------------
139 vtkIdType GetBucketIndex(const double* x) const
140 {
141 int ijk[3];
142 this->GetBucketIndices(x, ijk);
143 return ijk[0] + ijk[1] * xD + ijk[2] * xyD;
144 }
145
146 //-----------------------------------------------------------------------------
147 // Return the center of the bucket/bin at (i,j,k).
148 void GetBucketCenter(int i, int j, int k, double center[3])
149 {
150 center[0] = this->bX + this->hX2 + i * this->hX;
151 center[1] = this->bY + this->hY2 + j * this->hY;
152 center[2] = this->bZ + this->hZ2 + k * this->hZ;
153 }
154
155 //-----------------------------------------------------------------------------
156 // Return the bounding box (min,max) of a specified bucket (i,j,k).
157 void GetBucketBounds(int i, int j, int k, double min[3], double max[3])
158 {
159 min[0] = this->bX + i * this->hX;
160 min[1] = this->bY + j * this->hY;
161 min[2] = this->bZ + k * this->hZ;
162 max[0] = min[0] + this->hX;
163 max[1] = min[1] + this->hY;
164 max[2] = min[2] + this->hZ;
165 }
166
167 //-----------------------------------------------------------------------------
168 // Determine whether a bin/bucket specified by i,j,k is completely contained
169 // inside the sphere (center,r2). Return true if contained; false otherwise.
170 bool BucketInsideSphere(int i, int j, int k, double center[3], double r2)
171 {
172 double min[3], max[3];
173 min[0] = this->bX + i * this->hX;
174 min[1] = this->bY + j * this->hY;
175 min[2] = this->bZ + k * this->hZ;
176 max[0] += this->hX;
177 max[1] += this->hY;
178 max[2] += this->hZ;
179 return vtkBoundingBox::InsideSphere(min, max, center, r2);
180 }
181}; // vtkBucketList
182
183//------------------------------------------------------------------------------
184// This templates class manages the creation of the static locator
185// structures. It also implements the operator() functors which are supplied
186// to vtkSMPTools for threaded processesing.
187template <typename TIds>
189{
190 // Okay the various ivars
191 vtkLocatorTuple<TIds>* Map; // the map to be sorted
192 TIds* Offsets; // offsets for each bucket into the map
193
194 // Construction
195 BucketList(vtkStaticPointLocator* loc, vtkIdType numPts, int numBuckets)
196 : vtkBucketList(loc, numPts, numBuckets)
197 {
198 // one extra to simplify traversal
199 this->Map = new vtkLocatorTuple<TIds>[numPts + 1];
200 this->Map[numPts].Bucket = numBuckets;
201 this->Offsets = new TIds[numBuckets + 1];
202 this->Offsets[numBuckets] = numPts;
203 }
204
205 // Release allocated memory
206 ~BucketList() override
207 {
208 delete[] this->Map;
209 delete[] this->Offsets;
210 }
211
212 // The number of point ids in a bucket is determined by computing the
213 // difference between the offsets into the sorted points array.
215 {
216 return (this->Offsets[bucketNum + 1] - this->Offsets[bucketNum]);
217 }
218
219 // Given a bucket number, return the point ids in that bucket.
221 {
222 return this->Map + this->Offsets[bucketNum];
223 }
224
225 // Given a bucket number, return the point ids in that bucket.
226 void GetIds(vtkIdType bucketNum, vtkIdList* bList)
227 {
228 const vtkLocatorTuple<TIds>* ids = this->GetIds(bucketNum);
229 vtkIdType numIds = this->GetNumberOfIds(bucketNum);
230 bList->SetNumberOfIds(numIds);
231 for (int i = 0; i < numIds; i++)
232 {
233 bList->SetId(i, ids[i].PtId);
234 }
235 }
236
237 // Templated implementations of the locator
238 vtkIdType FindClosestPoint(const double x[3]);
240 double radius, const double x[3], double inputDataLength, double& dist2);
241 void FindClosestNPoints(int N, const double x[3], vtkIdList* result);
242 double FindNPointsInShell(int N, const double x[3], vtkDist2TupleArray& results,
243 double minDist2 = (-0.1), bool sort = true, vtkDoubleArray* petals = nullptr);
244 void FindPointsWithinRadius(double R, const double x[3], vtkIdList* result);
245 int IntersectWithLine(double a0[3], double a1[3], double tol, double& t, double lineX[3],
246 double ptX[3], vtkIdType& ptId);
247 void MergePoints(double tol, vtkIdType* pointMap, int orderingMode);
249 void GenerateRepresentation(int vtkNotUsed(level), vtkPolyData* pd);
250
251 // Internal methods
253 NeighborBuckets* buckets, const double x[3], const int ijk[3], double dist, int level);
254 void GetOverlappingBuckets(NeighborBuckets* buckets, const double x[3], double dist,
255 int prevMinLevel[3], int prevMaxLevel[3]);
256
257 // Implicit point representation, slower path
258 template <typename T>
260 {
263
265 : BList(blist)
266 , DataSet(ds)
267 {
268 }
269
271 {
272 double p[3];
273 vtkLocatorTuple<T>* t = this->BList->Map + ptId;
274
275 for (; ptId < end; ++ptId, ++t)
276 {
277 this->DataSet->GetPoint(ptId, p);
278 t->Bucket = this->BList->GetBucketIndex(p);
279 t->PtId = ptId;
280 } // for all points in this batch
281 }
282 };
283
284 // Explicit point representation (e.g., vtkPointSet), faster path
285 template <typename T, typename TPts>
287 {
289 const TPts* Points;
290
291 MapPointsArray(BucketList<T>* blist, const TPts* pts)
292 : BList(blist)
293 , Points(pts)
294 {
295 }
296
298 {
299 double p[3];
300 const TPts* x = this->Points + 3 * ptId;
301 vtkLocatorTuple<T>* t = this->BList->Map + ptId;
302
303 for (; ptId < end; ++ptId, x += 3, ++t)
304 {
305 p[0] = static_cast<double>(x[0]);
306 p[1] = static_cast<double>(x[1]);
307 p[2] = static_cast<double>(x[2]);
308 t->Bucket = this->BList->GetBucketIndex(p);
309 t->PtId = ptId;
310 } // for all points in this batch
311 }
312 };
313
314 // A clever way to build offsets in parallel. Basically each thread builds
315 // offsets across a range of the sorted map. Recall that offsets are an
316 // integral value referring to the locations of the sorted points that
317 // reside in each bucket.
318 template <typename T>
320 {
324
326 : BList(blist)
327 {
328 this->NumBuckets = this->BList->NumBuckets;
329 this->NumPts = this->BList->NumPts;
330 }
331
332 // Traverse sorted points (i.e., tuples) and update bucket offsets.
333 void operator()(vtkIdType batch, vtkIdType batchEnd)
334 {
335 T* offsets = this->BList->Offsets;
336 const vtkLocatorTuple<T>* curPt = this->BList->Map + batch * this->BList->BatchSize;
337 const vtkLocatorTuple<T>* endBatchPt = this->BList->Map + batchEnd * this->BList->BatchSize;
338 const vtkLocatorTuple<T>* endPt = this->BList->Map + this->NumPts;
339 const vtkLocatorTuple<T>* prevPt;
340 endBatchPt = (endBatchPt > endPt ? endPt : endBatchPt);
341
342 // Special case at the very beginning of the mapped points array. If
343 // the first point is in bucket# N, then all buckets up and including
344 // N must refer to the first point.
345 if (curPt == this->BList->Map)
346 {
347 prevPt = this->BList->Map;
348 std::fill_n(offsets, curPt->Bucket + 1, 0); // point to the first points
349 } // at the very beginning of the map (sorted points array)
350
351 // We are entering this functor somewhere in the interior of the
352 // mapped points array. All we need to do is point to the entry
353 // position because we are interested only in prevPt->Bucket.
354 else
355 {
356 prevPt = curPt;
357 } // else in the middle of a batch
358
359 // Okay we have a starting point for a bucket run. Now we can begin
360 // filling in the offsets in this batch. A previous thread should
361 // have/will have completed the previous and subsequent runs outside
362 // of the [batch,batchEnd) range
363 for (curPt = prevPt; curPt < endBatchPt;)
364 {
365 for (; curPt->Bucket == prevPt->Bucket && curPt <= endBatchPt; ++curPt)
366 {
367 // advance
368 }
369 // Fill in any gaps in the offset array
370 std::fill_n(
371 offsets + prevPt->Bucket + 1, curPt->Bucket - prevPt->Bucket, curPt - this->BList->Map);
372 prevPt = curPt;
373 } // for all batches in this range
374 } // operator()
375 };
376
377 // Merge points that are pecisely coincident. Operates in parallel on
378 // locator buckets. Does not need to check neighbor buckets.
379 template <typename T>
381 {
385
387 : BList(blist)
388 , MergeMap(mergeMap)
389 {
390 this->DataSet = blist->DataSet;
391 }
392
393 void operator()(vtkIdType bucket, vtkIdType endBucket)
394 {
395 BucketList<T>* bList = this->BList;
396 vtkIdType* mergeMap = this->MergeMap;
397 int i, j;
398 const vtkLocatorTuple<TIds>* ids;
399 double p[3], p2[3];
400 vtkIdType ptId, ptId2, numIds;
401
402 for (; bucket < endBucket; ++bucket)
403 {
404 if ((numIds = bList->GetNumberOfIds(bucket)) > 0)
405 {
406 ids = bList->GetIds(bucket);
407 for (i = 0; i < numIds; i++)
408 {
409 ptId = ids[i].PtId;
410 if (mergeMap[ptId] < 0)
411 {
412 mergeMap[ptId] = ptId;
413 this->DataSet->GetPoint(ptId, p);
414 for (j = i + 1; j < numIds; j++)
415 {
416 ptId2 = ids[j].PtId;
417 if (mergeMap[ptId2] < 0)
418 {
419 this->DataSet->GetPoint(ptId2, p2);
420 if (p[0] == p2[0] && p[1] == p2[1] && p[2] == p2[2])
421 {
422 mergeMap[ptId2] = ptId;
423 }
424 }
425 }
426 } // if point not yet visited
427 }
428 }
429 }
430 }
431 };
432
433 // Merge points that are coincident within a specified tolerance. Depending
434 // on the orderingMode, either a serialized ordering process is used (i.e.,
435 // POINT_ORDER) or a threaded ordering process is used (i.e.,
436 // BIN_ORDER). Note that due to the tolerance, the merging tolerance
437 // needs to check neighbor buckets which slows the algorithm down
438 // considerably. Note that merging is in one direction: larger ids are
439 // merged to lower ids.
440 template <typename T>
442 {
446 double Tol;
447
449
450 MergeClose(BucketList<T>* blist, double tol, vtkIdType* mergeMap)
451 : BList(blist)
452 , MergeMap(mergeMap)
453 , Tol(tol)
454 {
455 this->DataSet = blist->DataSet;
456 }
457
458 // The core merging process around the point ptId.
459 void MergePoint(vtkIdType ptId, vtkIdList* nearby)
460 {
461 vtkIdType* mergeMap = this->MergeMap;
462
463 // Make sure the point is not already merged
464 if (mergeMap[ptId] < 0)
465 {
466 mergeMap[ptId] = ptId;
467 double p[3];
468 this->DataSet->GetPoint(ptId, p);
469 this->BList->FindPointsWithinRadius(this->Tol, p, nearby);
470 vtkIdType numIds = nearby->GetNumberOfIds();
471 if (numIds > 0)
472 {
473 for (auto i = 0; i < numIds; ++i)
474 {
475 vtkIdType nearId = nearby->GetId(i);
476 if (mergeMap[nearId] < 0)
477 {
478 mergeMap[nearId] = ptId;
479 } // if eligible for merging and not yet merged
480 } // for all nearby points
481 } // if nearby points exist
482 } // if point not yet merged
483 } // MergePoint
484
485 // Just allocate a little bit of memory to get started.
487 {
488 vtkIdList*& pIds = this->PIds.Local();
489 pIds->Allocate(128); // allocate some memory
490 }
491
492 void Reduce() {}
493 };
494
495 // Merge points with non-zero tolerance. Order of point merging guarantees
496 // that any two merged point ids (p0,p1) are such that p0<p1. Consequently
497 // this is a completely serial algorithm.
498 template <typename T>
499 struct MergePointOrder : public MergeClose<T>
500 {
501 MergePointOrder(BucketList<T>* blist, double tol, vtkIdType* mergeMap)
502 : MergeClose<T>(blist, tol, mergeMap)
503 {
504 }
505
507
508 // Process serially, point by point.
510 {
511 vtkIdList*& nearby = this->PIds.Local();
512
513 // Serial operation over all points in the locator.
514 for (vtkIdType ptId = 0; ptId < numPts; ++ptId)
515 {
516 this->MergePoint(ptId, nearby);
517 } // for all points in the locator
518 } // operator()
519
520 void Reduce() { this->MergeClose<T>::Reduce(); }
521 }; // Merge points in point ordering
522
523 // Merge points with non-zero tolerance. The order of point merging depends
524 // on the order in which the bins are traversed (using a checkerboard
525 // pattern). While the algorithm is threaded, the checkerboarding acts as
526 // a barrier to full threading so the performance is not optimal (but at
527 // least deterministic / reproducible).
528 //
529 // Checkerboarding works as follows. The locator bin volume of dimensions
530 // Divisions[3] is divided into a collection of "blocks" which are
531 // subvolumes of bins of dimensions d^3. The algorithm makes multiple,
532 // threaded passes over the blocks (a total of d^3 threaded traversals),
533 // choosing one of the bins in each block to process via the current
534 // checkerboard index. The dimension d of the blocks is determined by the
535 // tolerance and locator bin size, and is chosen in such a way as to
536 // separate the point merging computation so as to avoid threading data
537 // races / write contention.
538 template <typename T>
539 struct MergeBinOrder : public MergeClose<T>
540 {
541 int CheckerboardDimension; // the dimension of the checkerboard block/subvolume
542 int NumBlocks; // how many blocks/subvolumes are in the binned locator
543 int BlockDims[3]; // the number of blocks in each coordinate direction
544 int CheckerboardIndex[3]; // which bin is being processed in the blocks
545
546 // The main function of the constructor is the setup the checkerboard
547 // traversal. This means configuring the checkerboard subvolume, and
548 // set up the traversal indices.
549 MergeBinOrder(BucketList<T>* blist, double tol, vtkIdType* mergeMap)
550 : MergeClose<T>(blist, tol, mergeMap)
551 {
552 BucketList<T>* bl = this->BList;
553 double hMin = std::min({ bl->hX, bl->hY, bl->hZ });
554 this->CheckerboardDimension =
555 1 + (hMin == 0.0 ? 1 : (1 + vtkMath::Floor(tol / (hMin / 2.0))));
556
557 // Determine how many blocks there are in the locater, and determine the
558 // dimensions of the blocks.
559 this->NumBlocks = 1;
560 for (auto i = 0; i < 3; ++i)
561 {
562 double numBlocks =
563 static_cast<double>(bl->Divisions[i]) / static_cast<double>(this->CheckerboardDimension);
564 this->BlockDims[i] = (bl->Divisions[i] <= 1 ? 1 : vtkMath::Ceil(numBlocks));
565 this->NumBlocks *= this->BlockDims[i];
566 }
568 }
569
570 // Initialize the checkerboard traversal process. A pointer to the
571 // current traversal state (within the checkerboard region) is returned.
573 {
574 // Control checkerboard traversal
575 this->CheckerboardIndex[0] = 0;
576 this->CheckerboardIndex[1] = 0;
577 this->CheckerboardIndex[2] = 0;
578 return this->CheckerboardIndex;
579 }
580
581 // Given a blockId and the current checkerboard index, compute the
582 // current locator bin/bucket id. May return <0 if no bin exists.
583 vtkIdType GetCurrentBin(int blockId, int cIdx[3])
584 {
585 // Which checkerboard block are we in?
586 int ijk[3];
587 vtkStructuredData::ComputePointStructuredCoords(blockId, this->BlockDims, ijk);
588
589 // Combine the block index with the checkerboard index. Make sure that
590 // we are still inside the locator bins (partial blocks may exist at
591 // the boundary). Recall that the blocks are composed of d^3 bins.
592 for (auto i = 0; i < 3; ++i)
593 {
594 ijk[i] = ijk[i] * this->CheckerboardDimension + cIdx[i];
595 if (ijk[i] >= this->BList->Divisions[i])
596 {
597 return (-1);
598 }
599 }
600
601 // Okay return the bin index
602 return (ijk[0] + ijk[1] * this->BList->Divisions[0] +
603 ijk[2] * this->BList->Divisions[0] * this->BList->Divisions[1]);
604 }
605
607
608 // Process locator blocks/subvolumes.
609 void operator()(vtkIdType blockId, vtkIdType endBlockId)
610 {
611 vtkIdList*& nearby = this->PIds.Local();
612 for (; blockId < endBlockId; ++blockId)
613 {
614 vtkIdType bin = this->GetCurrentBin(blockId, this->CheckerboardIndex);
615 vtkIdType numIds;
616
617 if (bin >= 0 && (numIds = this->BList->GetNumberOfIds(bin)) > 0)
618 {
619 const vtkLocatorTuple<TIds>* ids = this->BList->GetIds(bin);
620 for (auto i = 0; i < numIds; ++i)
621 {
622 vtkIdType ptId = ids[i].PtId;
623 this->MergePoint(ptId, nearby);
624 } // for all points in bin/bucket
625 } // if points exist in bin/bucket
626 } // for all blocks
627 } // operator()
628
629 void Reduce() { this->MergeClose<T>::Reduce(); }
630
631 // Coordinate the checkerboard threading process. Checkerboarding simply
632 // processes a subset of the locator bins to avoid write contention. The
633 // checkerboard footprint (its subvolume size) is a function of the
634 // tolerance, and is effectively a d^3 subvolume that is traversed
635 // (across all subvolumes) in a synchronized fashion. Hence there are d^3
636 // separate SMP traversals - if d becomes too large, the fallback is
637 // simply a serial (MergePointOrder()) to avoid thread thrashing.
638 void Execute()
639 {
640 int cDim = this->CheckerboardDimension;
641 int* cIdx = this->InitializeCheckerboardIndex();
642
643 // Coordinate the checkerboarding by synchronized traversal of the
644 // the checkerboard subblocks.
645 for (cIdx[2] = 0; cIdx[2] < cDim; ++cIdx[2])
646 {
647 for (cIdx[1] = 0; cIdx[1] < cDim; ++cIdx[1])
648 {
649 for (cIdx[0] = 0; cIdx[0] < cDim; ++cIdx[0])
650 {
651 vtkSMPTools::For(0, this->NumBlocks, *this);
652 }
653 }
654 }
655 } // Execute()
656
657 }; // MergeBinOrder
658
659 // Merge points that are geometrically coincident and have matching data
660 // values. Operates in parallel on locator buckets. Does not need to check
661 // neighbor buckets.
662 template <typename T>
664 {
671
673 : BList(blist)
674 , DataArray(da)
675 , MergeMap(mergeMap)
676 {
677 this->DataSet = blist->DataSet;
678 }
679
680 bool TuplesEqual(int tupleSize, double* t1, double* t2)
681 {
682 for (auto i = 0; i < tupleSize; ++i)
683 {
684 if (t1[i] != t2[i])
685 {
686 return false;
687 }
688 }
689 return true;
690 }
691
693 {
694 vtkIdType numComp = this->DataArray->GetNumberOfComponents();
695 this->Tuple.Local().resize(numComp);
696 this->Tuple2.Local().resize(numComp);
697 }
698
699 void operator()(vtkIdType bucket, vtkIdType endBucket)
700 {
701 BucketList<T>* bList = this->BList;
702 vtkIdType* mergeMap = this->MergeMap;
703 int i, j;
704 const vtkLocatorTuple<TIds>* ids;
705 double p[3], p2[3];
706 vtkIdType ptId, ptId2, numIds;
707 int tupleSize = static_cast<int>(this->Tuple.Local().size());
708 double* t = this->Tuple.Local().data();
709 double* t2 = this->Tuple2.Local().data();
710
711 for (; bucket < endBucket; ++bucket)
712 {
713 if ((numIds = bList->GetNumberOfIds(bucket)) > 0)
714 {
715 ids = bList->GetIds(bucket);
716 for (i = 0; i < numIds; i++)
717 {
718 ptId = ids[i].PtId;
719 if (mergeMap[ptId] < 0)
720 {
721 mergeMap[ptId] = ptId;
722 this->DataSet->GetPoint(ptId, p);
723 this->DataArray->GetTuple(ptId, t);
724 for (j = i + 1; j < numIds; j++)
725 {
726 ptId2 = ids[j].PtId;
727 if (mergeMap[ptId2] < 0)
728 {
729 this->DataSet->GetPoint(ptId2, p2);
730 if (p[0] == p2[0] && p[1] == p2[1] && p[2] == p2[2])
731 {
732 this->DataArray->GetTuple(ptId2, t2);
733 if (this->TuplesEqual(tupleSize, t, t2))
734 {
735 mergeMap[ptId2] = ptId;
736 } // if point's data match
737 } // if points geometrically coincident
738 } // if point not yet visited
739 } // for the remaining points in the bin
740 } // if point not yet merged
741 } // for all points in bucket
742 } // if bucket contains points
743 } // for all buckets
744 } // operator()
745
746 void Reduce() {}
747 }; // MergePointsWithData
748
749 // Build the map and other structures to support locator operations
750 void BuildLocator() override
751 {
752 // Place each point in a bucket
753 //
755 if (ps)
756 { // map points array: explicit points representation of float or double
757 int dataType = ps->GetPoints()->GetDataType();
758 void* pts = ps->GetPoints()->GetVoidPointer(0);
759 if (dataType == VTK_FLOAT)
760 {
761 MapPointsArray<TIds, float> mapper(this, static_cast<float*>(pts));
762 vtkSMPTools::For(0, this->NumPts, mapper);
763 }
764 else if (dataType == VTK_DOUBLE)
765 {
766 this->FastPoints = static_cast<double*>(pts);
767 MapPointsArray<TIds, double> mapper(this, this->FastPoints);
768 vtkSMPTools::For(0, this->NumPts, mapper);
769 }
770 }
771 else // if (!mapped)
772 { // map dataset points: non-float points or implicit points representation
773 MapDataSet<TIds> mapper(this, this->DataSet);
774 vtkSMPTools::For(0, this->NumPts, mapper);
775 }
776
777 // Now group the points into contiguous runs within buckets (recall that
778 // sorting is occurring based on bin/bucket id).
779 vtkSMPTools::Sort(this->Map, this->Map + this->NumPts);
780
781 // Build the offsets into the Map. The offsets are the positions of
782 // each bucket into the sorted list. They mark the beginning of the
783 // list of points in each bucket. Amazingly, this can be done in
784 // parallel.
785 int numBatches = static_cast<int>(ceil(static_cast<double>(this->NumPts) / this->BatchSize));
786 MapOffsets<TIds> offMapper(this);
787 vtkSMPTools::For(0, numBatches, offMapper);
788 }
789}; // BucketList
790
791VTK_ABI_NAMESPACE_END
792#endif // vtkStaticPointLocatorPrivate_h
793// VTK-HeaderTest-Exclude: vtkStaticPointLocatorPrivate.h
RealT r2
Definition PyrC2Basis.h:20
RealT t2
Definition PyrC2Basis.h:22
int GetNumberOfComponents() const
Set/Get the dimension (n) of the components.
static bool InsideSphere(const double min[3], const double max[3], const double center[3], double r2)
Performant method to determine if box if fully inside a sphere.
object to represent cell connectivity
virtual double * GetTuple(vtkIdType tupleIdx)=0
Get the data tuple at tupleIdx.
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
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.
static int Ceil(double x)
Rounds a double to the nearest integer not less than itself.
Definition vtkMath.h:1945
static int Floor(double x)
Rounds a double to the nearest integer not greater than itself.
Definition vtkMath.h:1936
concrete class for storing a set of points
Definition vtkPointSet.h:98
static vtkPointSet * SafeDownCast(vtkObjectBase *o)
Standard methods for type information and printing.
vtkPoints * GetPoints() override
Specify point array to define point coordinates.
represent and manipulate 3D points
Definition vtkPoints.h:139
void * GetVoidPointer(const int id)
Return a void pointer.
Definition vtkPoints.h:196
virtual int GetDataType() const
Return the underlying data type.
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.
Thread local storage for VTK objects.
T & Local()
This needs to be called mainly within a threaded execution path.
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 3-space
virtual double * GetBounds()
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-z directions.
static void ComputePointStructuredCoords(vtkIdType ptId, const int dim[3], int ijk[3], int dataDescription=vtkStructuredData::VTK_STRUCTURED_EMPTY)
Given a pointId and grid dimensions 'dim', get the structured coordinates (i-j-k).
MapDataSet(BucketList< T > *blist, vtkDataSet *ds)
void operator()(vtkIdType ptId, vtkIdType end)
void operator()(vtkIdType batch, vtkIdType batchEnd)
MapPointsArray(BucketList< T > *blist, const TPts *pts)
void operator()(vtkIdType ptId, vtkIdType end)
void operator()(vtkIdType blockId, vtkIdType endBlockId)
MergeBinOrder(BucketList< T > *blist, double tol, vtkIdType *mergeMap)
vtkIdType GetCurrentBin(int blockId, int cIdx[3])
vtkSMPThreadLocalObject< vtkIdList > PIds
void MergePoint(vtkIdType ptId, vtkIdList *nearby)
MergeClose(BucketList< T > *blist, double tol, vtkIdType *mergeMap)
MergePointOrder(BucketList< T > *blist, double tol, vtkIdType *mergeMap)
MergePointsAndData(BucketList< T > *blist, vtkDataArray *da, vtkIdType *mergeMap)
bool TuplesEqual(int tupleSize, double *t1, double *t2)
vtkSMPThreadLocal< std::vector< double > > Tuple
void operator()(vtkIdType bucket, vtkIdType endBucket)
vtkSMPThreadLocal< std::vector< double > > Tuple2
MergePrecise(BucketList< T > *blist, vtkIdType *mergeMap)
void operator()(vtkIdType bucket, vtkIdType endBucket)
void GetIds(vtkIdType bucketNum, vtkIdList *bList)
BucketList(vtkStaticPointLocator *loc, vtkIdType numPts, int numBuckets)
void GetOverlappingBuckets(NeighborBuckets *buckets, const double x[3], double dist, int prevMinLevel[3], int prevMaxLevel[3])
double FindNPointsInShell(int N, const double x[3], vtkDist2TupleArray &results, double minDist2=(-0.1), bool sort=true, vtkDoubleArray *petals=nullptr)
void GenerateRepresentation(int level, vtkPolyData *pd)
void MergePointsWithData(vtkDataArray *data, vtkIdType *pointMap)
vtkIdType FindClosestPoint(const double x[3])
vtkIdType GetNumberOfIds(vtkIdType bucketNum)
void FindPointsWithinRadius(double R, const double x[3], vtkIdList *result)
vtkLocatorTuple< TIds > * Map
void BuildLocator() override
const vtkLocatorTuple< TIds > * GetIds(vtkIdType bucketNum)
vtkIdType FindClosestPointWithinRadius(double radius, const double x[3], double inputDataLength, double &dist2)
int IntersectWithLine(double a0[3], double a1[3], double tol, double &t, double lineX[3], double ptX[3], vtkIdType &ptId)
void GetOverlappingBuckets(NeighborBuckets *buckets, const double x[3], const int ijk[3], double dist, int level)
void MergePoints(double tol, vtkIdType *pointMap, int orderingMode)
void FindClosestNPoints(int N, const double x[3], vtkIdList *result)
void GenerateFace(int face, int i, int j, int k, vtkPoints *pts, vtkCellArray *polys)
vtkIdType GetBucketIndex(const double *x) const
virtual void BuildLocator()=0
void GetBucketNeighbors(NeighborBuckets *buckets, const int ijk[3], const int ndivs[3], int level)
vtkBucketList(vtkStaticPointLocator *loc, vtkIdType numPts, int numBuckets)
double Distance2ToBucket(const double x[3], const int nei[3])
virtual ~vtkBucketList()=default
vtkStaticPointLocator * Locator
double Distance2ToBounds(const double x[3], const double bounds[6])
bool BucketInsideSphere(int i, int j, int k, double center[3], double r2)
void GetBucketIndices(const double *x, int ijk[3]) const
void GetBucketBounds(int i, int j, int k, double min[3], double max[3])
void GetBucketCenter(int i, int j, int k, double center[3])
Represent an array of vtkDist2Tuples.
#define vtkDataArray
int vtkIdType
Definition vtkType.h:367
#define VTK_DOUBLE
Definition vtkType.h:44
#define VTK_FLOAT
Definition vtkType.h:43
#define max(a, b)