xref: /petsc/src/vec/is/utils/kdtree.c (revision 98480730077ed2eb9d6c4a96b667530324fc8426)
1 #include <petsc.h>
2 #include <petscis.h>
3 #include <petsc/private/petscimpl.h>
4 
5 typedef struct {
6   PetscInt   axis;                                // Could make this into a uint8_t to save on memory and bandwidth?
7   PetscBool  is_greater_leaf, is_less_equal_leaf; // Could possibly use PetscBT to save on memory and bandwidth?
8   PetscReal  split;
9   PetscCount greater_handle, less_equal_handle;
10 } KDStem;
11 
12 typedef struct {
13   PetscInt   count;
14   PetscCount indices_handle, coords_handle;
15 } KDLeaf;
16 
17 struct _n_PetscKDTree {
18   PetscInt dim;
19   PetscInt max_bucket_size;
20 
21   PetscBool  is_root_leaf;
22   PetscCount root_handle;
23 
24   PetscCount       num_coords, num_leaves, num_stems, num_bucket_indices;
25   const PetscReal *coords, *coords_owned; // Only free owned on Destroy
26   KDLeaf          *leaves;
27   KDStem          *stems;
28   PetscCount      *bucket_indices;
29 };
30 
31 /*@C
32   PetscKDTreeDestroy - destroy a `PetscKDTree`
33 
34   Not Collective, No Fortran Support
35 
36   Input Parameters:
37 . tree - tree to destroy
38 
39   Level: advanced
40 
41 .seealso: `PetscKDTree`, `PetscKDTreeCreate()`
42 @*/
43 PetscErrorCode PetscKDTreeDestroy(PetscKDTree *tree)
44 {
45   PetscFunctionBeginUser;
46   if (*tree == NULL) PetscFunctionReturn(PETSC_SUCCESS);
47   PetscCall(PetscFree((*tree)->stems));
48   PetscCall(PetscFree((*tree)->leaves));
49   PetscCall(PetscFree((*tree)->bucket_indices));
50   PetscCall(PetscFree((*tree)->coords_owned));
51   PetscCall(PetscFree(*tree));
52   PetscFunctionReturn(PETSC_SUCCESS);
53 }
54 
55 PetscLogEvent         PetscKDTree_Build, PetscKDTree_Query;
56 static PetscErrorCode PetscKDTreeRegisterLogEvents()
57 {
58   static PetscBool is_initialized = PETSC_FALSE;
59 
60   PetscFunctionBeginUser;
61   if (is_initialized) PetscFunctionReturn(PETSC_SUCCESS);
62   PetscCall(PetscLogEventRegister("KDTreeBuild", IS_CLASSID, &PetscKDTree_Build));
63   PetscCall(PetscLogEventRegister("KDTreeQuery", IS_CLASSID, &PetscKDTree_Query));
64   PetscFunctionReturn(PETSC_SUCCESS);
65 }
66 
67 // From http://graphics.stanford.edu/~seander/bithacks.html#RoundUpPowerOf2
68 static inline uint32_t RoundToNextPowerOfTwo(uint32_t v)
69 {
70   v--;
71   v |= v >> 1;
72   v |= v >> 2;
73   v |= v >> 4;
74   v |= v >> 8;
75   v |= v >> 16;
76   v++;
77   return v;
78 }
79 
80 typedef struct {
81   PetscInt    initial_axis;
82   PetscKDTree tree;
83 } *KDTreeSortContext;
84 
85 // Sort nodes based on "superkey"
86 // See "Building a Balanced k-d Tree in O(kn log n) Time" https://jcgt.org/published/0004/01/03/
87 static inline int PetscKDTreeSortFunc(PetscCount left, PetscCount right, PetscKDTree tree, PetscInt axis)
88 {
89   const PetscReal *coords = tree->coords;
90   const PetscInt   dim    = tree->dim;
91 
92   for (PetscInt i = 0; i < dim; i++) {
93     PetscReal diff = coords[left * dim + axis] - coords[right * dim + axis];
94     if (PetscUnlikely(diff == 0)) {
95       axis = (axis + 1) % dim;
96       continue;
97     } else return PetscSign(diff);
98   }
99   return 0; // All components are the same
100 }
101 
102 static int PetscKDTreeTimSort(const void *l, const void *r, void *ctx)
103 {
104   KDTreeSortContext kd_ctx = (KDTreeSortContext)ctx;
105   return PetscKDTreeSortFunc(*(PetscCount *)l, *(PetscCount *)r, kd_ctx->tree, kd_ctx->initial_axis);
106 }
107 
108 static PetscErrorCode PetscKDTreeVerifySortedIndices(PetscKDTree tree, PetscCount sorted_indices[], PetscCount temp[], PetscCount start, PetscCount end)
109 {
110   PetscCount num_coords = tree->num_coords, range_size = end - start, location;
111   PetscBool  has_duplicates;
112 
113   PetscFunctionBeginUser;
114   PetscCall(PetscArraycpy(temp, &sorted_indices[0 * num_coords + start], range_size));
115   PetscCall(PetscSortCount(range_size, temp));
116   PetscCall(PetscSortedCheckDupsCount(range_size, temp, &has_duplicates));
117   PetscCheck(has_duplicates == PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Sorted indices must have unique entries, but found duplicates");
118   for (PetscInt d = 1; d < tree->dim; d++) {
119     for (PetscCount i = start; i < end; i++) {
120       PetscCall(PetscFindCount(sorted_indices[d * num_coords + i], range_size, temp, &location));
121       PetscCheck(location > -1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Sorted indices are not consistent. Could not find %" PetscCount_FMT " from %" PetscInt_FMT " dimensional index in 0th dimension", sorted_indices[d * num_coords + i], d);
122     }
123   }
124   PetscFunctionReturn(PETSC_SUCCESS);
125 }
126 
127 typedef struct {
128   PetscKDTree    tree;
129   PetscSegBuffer stems, leaves, bucket_indices, bucket_coords;
130   PetscBool      debug_build, copy_coords;
131 } *KDTreeBuild;
132 
133 // The range is end exclusive, so [start,end).
134 static PetscErrorCode PetscKDTreeBuildStemAndLeaves(KDTreeBuild kd_build, PetscCount sorted_indices[], PetscCount temp[], PetscCount start, PetscCount end, PetscInt depth, PetscBool *is_node_leaf, PetscCount *node_handle)
135 {
136   PetscKDTree tree = kd_build->tree;
137   PetscInt    dim  = tree->dim;
138 
139   PetscFunctionBeginUser;
140   PetscCheck(start <= end, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Start %" PetscCount_FMT " must be less than or equal to end %" PetscCount_FMT, start, end);
141   if (kd_build->debug_build) PetscCall(PetscKDTreeVerifySortedIndices(tree, sorted_indices, temp, start, end));
142   if (end - start <= tree->max_bucket_size) {
143     KDLeaf     *leaf;
144     PetscCount *bucket_indices;
145 
146     PetscCall(PetscSegBufferGetSize(kd_build->leaves, node_handle));
147     PetscCall(PetscSegBufferGet(kd_build->leaves, 1, &leaf));
148     PetscCall(PetscMemzero(leaf, sizeof(KDLeaf)));
149     *is_node_leaf = PETSC_TRUE;
150 
151     PetscCall(PetscIntCast(end - start, &leaf->count));
152     PetscCall(PetscSegBufferGetSize(kd_build->bucket_indices, &leaf->indices_handle));
153     PetscCall(PetscSegBufferGet(kd_build->bucket_indices, leaf->count, &bucket_indices));
154     PetscCall(PetscArraycpy(bucket_indices, &sorted_indices[start], leaf->count));
155     if (kd_build->copy_coords) {
156       PetscReal *bucket_coords;
157       PetscCall(PetscSegBufferGetSize(kd_build->bucket_coords, &leaf->coords_handle));
158       PetscCall(PetscSegBufferGet(kd_build->bucket_coords, leaf->count * dim, &bucket_coords));
159       // Coords are saved in axis-major ordering for better vectorization
160       for (PetscCount i = 0; i < leaf->count; i++) {
161         for (PetscInt d = 0; d < dim; d++) bucket_coords[d * leaf->count + i] = tree->coords[bucket_indices[i] * dim + d];
162       }
163     } else leaf->coords_handle = -1;
164   } else {
165     KDStem    *stem;
166     PetscCount num_coords = tree->num_coords;
167     PetscInt   axis       = depth % dim;
168     PetscCount median     = start + PetscCeilInt64(end - start, 2) - 1, lower;
169     PetscCount median_idx = sorted_indices[median], medianp1_idx = sorted_indices[median + 1];
170 
171     PetscCall(PetscSegBufferGetSize(kd_build->stems, node_handle));
172     PetscCall(PetscSegBufferGet(kd_build->stems, 1, &stem));
173     PetscCall(PetscMemzero(stem, sizeof(KDStem)));
174     *is_node_leaf = PETSC_FALSE;
175 
176     stem->axis = axis;
177     // Place split halfway between the "boundary" nodes of the partitioning
178     stem->split = (tree->coords[tree->dim * median_idx + axis] + tree->coords[tree->dim * medianp1_idx + axis]) / 2;
179     PetscCall(PetscArraycpy(temp, &sorted_indices[0 * num_coords + start], end - start));
180     lower = median; // Set lower in case dim == 1
181     for (PetscInt d = 1; d < tree->dim; d++) {
182       PetscCount upper = median;
183       lower            = start - 1;
184       for (PetscCount i = start; i < end; i++) {
185         // In case of duplicate coord point equal to the median coord point, limit lower partition to median, ensuring balanced tree
186         if (lower < median && PetscKDTreeSortFunc(sorted_indices[d * num_coords + i], median_idx, tree, axis) <= 0) {
187           sorted_indices[(d - 1) * num_coords + (++lower)] = sorted_indices[d * num_coords + i];
188         } else {
189           sorted_indices[(d - 1) * num_coords + (++upper)] = sorted_indices[d * num_coords + i];
190         }
191       }
192       PetscCheck(lower == median, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Partitioning into less_equal bin failed. Range upper bound should be %" PetscCount_FMT " but partitioning resulted in %" PetscCount_FMT, median, lower);
193       PetscCheck(upper == end - 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Partitioning into greater bin failed. Range upper bound should be %" PetscCount_FMT " but partitioning resulted in %" PetscCount_FMT, upper, end - 1);
194     }
195     PetscCall(PetscArraycpy(&sorted_indices[(tree->dim - 1) * num_coords + start], temp, end - start));
196 
197     PetscCall(PetscKDTreeBuildStemAndLeaves(kd_build, sorted_indices, temp, start, lower + 1, depth + 1, &stem->is_less_equal_leaf, &stem->less_equal_handle));
198     PetscCall(PetscKDTreeBuildStemAndLeaves(kd_build, sorted_indices, temp, lower + 1, end, depth + 1, &stem->is_greater_leaf, &stem->greater_handle));
199   }
200   PetscFunctionReturn(PETSC_SUCCESS);
201 }
202 
203 /*@C
204   PetscKDTreeCreate - create a `PetscKDTree`
205 
206   Not Collective, No Fortran Support
207 
208   Input Parameters:
209 + num_coords      - number of coordinate points to build the `PetscKDTree`
210 . dim             - the dimension of the coordinates
211 . coords          - array of the coordinates, in point-major order
212 . copy_mode       - behavior handling `coords`, `PETSC_COPY_VALUES` generally more performant
213 - max_bucket_size - maximum number of points stored at each leaf
214 
215   Output Parameter:
216 . new_tree - the resulting `PetscKDTree`
217 
218   Level: advanced
219 
220   Note:
221   When `copy_mode == PETSC_COPY_VALUES`, the coordinates are copied and organized to optimize vectorization and cache-coherency.
222   It is recommended to run this way if the extra memory use is not a concern and it has very little impact on the `PetscKDTree` creation time.
223 
224   Developer Note:
225   Building algorithm detailed in 'Building a Balanced k-d Tree in O(kn log n) Time' Brown, 2015
226 
227 .seealso: `PetscKDTree`, `PetscKDTreeDestroy()`, `PetscKDTreeQueryPointsNearestNeighbor()`
228 @*/
229 PetscErrorCode PetscKDTreeCreate(PetscCount num_coords, PetscInt dim, const PetscReal coords[], PetscCopyMode copy_mode, PetscInt max_bucket_size, PetscKDTree *new_tree)
230 {
231   PetscKDTree tree;
232   PetscCount *sorted_indices, *temp;
233 
234   PetscFunctionBeginUser;
235   PetscCall(PetscKDTreeRegisterLogEvents());
236   PetscCall(PetscLogEventBegin(PetscKDTree_Build, 0, 0, 0, 0));
237   PetscCheck(dim > 0, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "Dimension of PetscKDTree must be greater than 0, recieved %" PetscInt_FMT, dim);
238   PetscCheck(num_coords > -1, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "Number of coordinates may not be negative, recieved %" PetscCount_FMT, num_coords);
239   if (num_coords == 0) {
240     *new_tree = NULL;
241     PetscFunctionReturn(PETSC_SUCCESS);
242   }
243   PetscAssertPointer(coords, 3);
244   PetscAssertPointer(new_tree, 6);
245   PetscCall(PetscNew(&tree));
246   tree->dim             = dim;
247   tree->max_bucket_size = max_bucket_size == PETSC_DECIDE ? 32 : max_bucket_size;
248   tree->num_coords      = num_coords;
249 
250   switch (copy_mode) {
251   case PETSC_OWN_POINTER:
252     tree->coords_owned = coords;
253   case PETSC_USE_POINTER:
254     tree->coords = coords;
255     break;
256   case PETSC_COPY_VALUES:
257     PetscCall(PetscMalloc1(num_coords * dim, &tree->coords_owned));
258     PetscCall(PetscArraycpy((PetscReal *)tree->coords_owned, coords, num_coords * dim));
259     tree->coords = tree->coords_owned;
260     break;
261   }
262 
263   KDTreeSortContext kd_ctx;
264   PetscCall(PetscMalloc2(num_coords * dim, &sorted_indices, num_coords, &temp));
265   PetscCall(PetscNew(&kd_ctx));
266   kd_ctx->tree = tree;
267   for (PetscInt j = 0; j < dim; j++) {
268     for (PetscCount i = 0; i < num_coords; i++) sorted_indices[num_coords * j + i] = i;
269     kd_ctx->initial_axis = j;
270     PetscCall(PetscTimSort((PetscInt)num_coords, &sorted_indices[num_coords * j], sizeof(*sorted_indices), PetscKDTreeTimSort, kd_ctx));
271   }
272   PetscCall(PetscFree(kd_ctx));
273 
274   PetscInt    num_leaves = (PetscInt)PetscCeilInt64(num_coords, tree->max_bucket_size);
275   PetscInt    num_stems  = RoundToNextPowerOfTwo((uint32_t)num_leaves);
276   KDTreeBuild kd_build;
277   PetscCall(PetscNew(&kd_build));
278   kd_build->tree        = tree;
279   kd_build->copy_coords = copy_mode == PETSC_COPY_VALUES ? PETSC_TRUE : PETSC_FALSE;
280   PetscCall(PetscOptionsGetBool(NULL, NULL, "-kdtree_debug", &kd_build->debug_build, NULL));
281   PetscCall(PetscSegBufferCreate(sizeof(KDStem), num_stems, &kd_build->stems));
282   PetscCall(PetscSegBufferCreate(sizeof(KDLeaf), num_leaves, &kd_build->leaves));
283   PetscCall(PetscSegBufferCreate(sizeof(PetscCount), num_coords, &kd_build->bucket_indices));
284   if (kd_build->copy_coords) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), num_coords * dim, &kd_build->bucket_coords));
285 
286   PetscCall(PetscKDTreeBuildStemAndLeaves(kd_build, sorted_indices, temp, 0, num_coords, 0, &tree->is_root_leaf, &tree->root_handle));
287 
288   PetscCall(PetscSegBufferGetSize(kd_build->stems, &tree->num_stems));
289   PetscCall(PetscSegBufferGetSize(kd_build->leaves, &tree->num_leaves));
290   PetscCall(PetscSegBufferGetSize(kd_build->bucket_indices, &tree->num_bucket_indices));
291   PetscCall(PetscSegBufferExtractAlloc(kd_build->stems, &tree->stems));
292   PetscCall(PetscSegBufferExtractAlloc(kd_build->leaves, &tree->leaves));
293   PetscCall(PetscSegBufferExtractAlloc(kd_build->bucket_indices, &tree->bucket_indices));
294   if (kd_build->copy_coords) {
295     PetscCall(PetscFree(tree->coords_owned));
296     PetscCall(PetscSegBufferExtractAlloc(kd_build->bucket_coords, &tree->coords_owned));
297     tree->coords = tree->coords_owned;
298     PetscCall(PetscSegBufferDestroy(&kd_build->bucket_coords));
299   }
300   PetscCall(PetscSegBufferDestroy(&kd_build->stems));
301   PetscCall(PetscSegBufferDestroy(&kd_build->leaves));
302   PetscCall(PetscSegBufferDestroy(&kd_build->bucket_indices));
303   PetscCall(PetscFree(kd_build));
304   PetscCall(PetscFree2(sorted_indices, temp));
305   *new_tree = tree;
306   PetscCall(PetscLogEventEnd(PetscKDTree_Build, 0, 0, 0, 0));
307   PetscFunctionReturn(PETSC_SUCCESS);
308 }
309 
310 static inline PetscReal PetscSquareDistance(PetscInt dim, const PetscReal *PETSC_RESTRICT x, const PetscReal *PETSC_RESTRICT y)
311 {
312   PetscReal dist = 0;
313   for (PetscInt j = 0; j < dim; j++) dist += PetscSqr(x[j] - y[j]);
314   return dist;
315 }
316 
317 static inline PetscErrorCode PetscKDTreeQueryLeaf(PetscKDTree tree, KDLeaf leaf, const PetscReal point[], PetscCount *index, PetscReal *distance_sqr)
318 {
319   PetscInt dim = tree->dim;
320 
321   PetscFunctionBeginUser;
322   *distance_sqr = PETSC_MAX_REAL;
323   *index        = -1;
324   for (PetscInt i = 0; i < leaf.count; i++) {
325     PetscCount point_index = tree->bucket_indices[leaf.indices_handle + i];
326     PetscReal  dist        = PetscSquareDistance(dim, point, &tree->coords[point_index * dim]);
327     if (dist < *distance_sqr) {
328       *distance_sqr = dist;
329       *index        = point_index;
330     }
331   }
332   PetscFunctionReturn(PETSC_SUCCESS);
333 }
334 
335 static inline PetscErrorCode PetscKDTreeQueryLeaf_CopyCoords(PetscKDTree tree, KDLeaf leaf, const PetscReal point[], PetscCount *index, PetscReal *distance_sqr)
336 {
337   PetscInt dim = tree->dim;
338 
339   PetscFunctionBeginUser;
340   *distance_sqr = PETSC_MAX_REAL;
341   *index        = -1;
342   for (PetscInt i = 0; i < leaf.count; i++) {
343     // Coord data saved in axis-major ordering for vectorization
344     PetscReal dist = 0.;
345     for (PetscInt d = 0; d < dim; d++) dist += PetscSqr(point[d] - tree->coords[leaf.coords_handle + d * leaf.count + i]);
346     if (dist < *distance_sqr) {
347       *distance_sqr = dist;
348       *index        = tree->bucket_indices[leaf.indices_handle + i];
349     }
350   }
351   PetscFunctionReturn(PETSC_SUCCESS);
352 }
353 
354 // Recursive point query from 'Algorithms for Fast Vector Quantization' by  Sunil Arya and David Mount
355 // Variant also implemented in pykdtree
356 static PetscErrorCode PetscKDTreeQuery_Recurse(PetscKDTree tree, const PetscReal point[], PetscCount node_handle, PetscBool is_node_leaf, PetscReal offset[], PetscReal rd, PetscReal tol_sqr, PetscCount *index, PetscReal *dist_sqr)
357 {
358   PetscFunctionBeginUser;
359   if (*dist_sqr < tol_sqr) PetscFunctionReturn(PETSC_SUCCESS);
360   if (is_node_leaf) {
361     KDLeaf     leaf = tree->leaves[node_handle];
362     PetscReal  dist;
363     PetscCount point_index;
364 
365     if (leaf.coords_handle > -1) PetscCall(PetscKDTreeQueryLeaf_CopyCoords(tree, leaf, point, &point_index, &dist));
366     else PetscCall(PetscKDTreeQueryLeaf(tree, leaf, point, &point_index, &dist));
367     if (dist < *dist_sqr) {
368       *dist_sqr = dist;
369       *index    = point_index;
370     }
371     PetscFunctionReturn(PETSC_SUCCESS);
372   }
373 
374   KDStem    stem       = tree->stems[node_handle];
375   PetscReal old_offset = offset[stem.axis], new_offset = point[stem.axis] - stem.split;
376   if (new_offset <= 0) {
377     PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.less_equal_handle, stem.is_less_equal_leaf, offset, rd, tol_sqr, index, dist_sqr));
378     rd += -PetscSqr(old_offset) + PetscSqr(new_offset);
379     if (rd < *dist_sqr) {
380       offset[stem.axis] = new_offset;
381       PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.greater_handle, stem.is_greater_leaf, offset, rd, tol_sqr, index, dist_sqr));
382       offset[stem.axis] = old_offset;
383     }
384   } else {
385     PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.greater_handle, stem.is_greater_leaf, offset, rd, tol_sqr, index, dist_sqr));
386     rd += -PetscSqr(old_offset) + PetscSqr(new_offset);
387     if (rd < *dist_sqr) {
388       offset[stem.axis] = new_offset;
389       PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.less_equal_handle, stem.is_less_equal_leaf, offset, rd, tol_sqr, index, dist_sqr));
390       offset[stem.axis] = old_offset;
391     }
392   }
393   PetscFunctionReturn(PETSC_SUCCESS);
394 }
395 
396 /*@C
397   PetscKDTreeQueryPointsNearestNeighbor - find the nearest neighbor in a `PetscKDTree`
398 
399   Not Collective, No Fortran Support
400 
401   Input Parameters:
402 + tree       - tree to query
403 . num_points - number of points to query
404 . points     - array of the coordinates, in point-major order
405 - tolerance  - tolerance for nearest neighbor
406 
407   Output Parameter:
408 + indices   - indices of the nearest neighbor to the query point
409 - distances - distance between the queried point and the nearest neighbor
410 
411   Level: advanced
412 
413   Notes:
414   When traversing the tree, if a point has been found to be closer than the `tolerance`, the function short circuits and doesn't check for any closer points.
415 
416   The `indices` and `distances` arrays should be at least of size `num_points`.
417 
418 .seealso: `PetscKDTree`, `PetscKDTreeCreate()`
419 @*/
420 PetscErrorCode PetscKDTreeQueryPointsNearestNeighbor(PetscKDTree tree, PetscCount num_points, const PetscReal points[], PetscReal tolerance, PetscCount indices[], PetscReal distances[])
421 {
422   PetscReal *offsets, rd;
423 
424   PetscFunctionBeginUser;
425   PetscCall(PetscLogEventBegin(PetscKDTree_Query, 0, 0, 0, 0));
426   if (tree == NULL) {
427     PetscCheck(num_points == 0, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "num_points may only be zero, if tree is NULL");
428     PetscFunctionReturn(PETSC_SUCCESS);
429   }
430   PetscAssertPointer(points, 3);
431   PetscAssertPointer(indices, 5);
432   PetscAssertPointer(distances, 6);
433   PetscCall(PetscCalloc1(tree->dim, &offsets));
434 
435   for (PetscCount p = 0; p < num_points; p++) {
436     rd           = 0.;
437     distances[p] = PETSC_MAX_REAL;
438     indices[p]   = -1;
439     PetscCall(PetscKDTreeQuery_Recurse(tree, &points[p * tree->dim], tree->root_handle, tree->is_root_leaf, offsets, rd, PetscSqr(tolerance), &indices[p], &distances[p]));
440     distances[p] = PetscSqrtReal(distances[p]);
441   }
442   PetscCall(PetscFree(offsets));
443   PetscCall(PetscLogEventEnd(PetscKDTree_Query, 0, 0, 0, 0));
444   PetscFunctionReturn(PETSC_SUCCESS);
445 }
446 
447 /*@C
448   PetscKDTreeView - view a `PetscKDTree`
449 
450   Not Collective, No Fortran Support
451 
452   Input Parameters:
453 + tree   - tree to view
454 - viewer - visualization context
455 
456   Level: advanced
457 
458 .seealso: `PetscKDTree`, `PetscKDTreeCreate()`, `PetscViewer`
459 @*/
460 PetscErrorCode PetscKDTreeView(PetscKDTree tree, PetscViewer viewer)
461 {
462   PetscFunctionBeginUser;
463   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
464   else PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &viewer));
465   if (tree == NULL) PetscFunctionReturn(PETSC_SUCCESS);
466 
467   PetscCall(PetscViewerASCIIPrintf(viewer, "KDTree:\n"));
468   PetscCall(PetscViewerASCIIPushTab(viewer)); // KDTree:
469   PetscCall(PetscViewerASCIIPrintf(viewer, "Stems:\n"));
470   PetscCall(PetscViewerASCIIPushTab(viewer)); // Stems:
471   for (PetscCount i = 0; i < tree->num_stems; i++) {
472     KDStem stem = tree->stems[i];
473     PetscCall(PetscViewerASCIIPrintf(viewer, "Stem %" PetscCount_FMT ": Axis=%" PetscInt_FMT " Split=%g Greater_%s=%" PetscCount_FMT " Lesser_Equal_%s=%" PetscCount_FMT "\n", i, stem.axis, (double)stem.split, stem.is_greater_leaf ? "Leaf" : "Stem",
474                                      stem.greater_handle, stem.is_less_equal_leaf ? "Leaf" : "Stem", stem.less_equal_handle));
475   }
476   PetscCall(PetscViewerASCIIPopTab(viewer)); // Stems:
477 
478   PetscCall(PetscViewerASCIIPrintf(viewer, "Leaves:\n"));
479   PetscCall(PetscViewerASCIIPushTab(viewer)); // Leaves:
480   for (PetscCount i = 0; i < tree->num_leaves; i++) {
481     KDLeaf leaf = tree->leaves[i];
482     PetscCall(PetscViewerASCIIPrintf(viewer, "Leaf %" PetscCount_FMT ": Count=%" PetscInt_FMT, i, leaf.count));
483     PetscCall(PetscViewerASCIIPushTab(viewer)); // Coords:
484     for (PetscInt j = 0; j < leaf.count; j++) {
485       PetscInt   tabs;
486       PetscCount bucket_index = tree->bucket_indices[leaf.indices_handle + j];
487       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
488       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscCount_FMT ": ", bucket_index));
489 
490       PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
491       PetscCall(PetscViewerASCIISetTab(viewer, 0));
492       if (leaf.coords_handle > -1) {
493         for (PetscInt k = 0; k < tree->dim; k++) PetscCall(PetscViewerASCIIPrintf(viewer, "%g ", (double)tree->coords[leaf.coords_handle + leaf.count * k + j]));
494         PetscCall(PetscViewerASCIIPrintf(viewer, " (stored at leaf)"));
495       } else {
496         for (PetscInt k = 0; k < tree->dim; k++) PetscCall(PetscViewerASCIIPrintf(viewer, "%g ", (double)tree->coords[bucket_index * tree->dim + k]));
497       }
498       PetscCall(PetscViewerASCIISetTab(viewer, tabs));
499     }
500     PetscCall(PetscViewerASCIIPopTab(viewer)); // Coords:
501     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
502   }
503   PetscCall(PetscViewerASCIIPopTab(viewer)); // Leaves:
504   PetscCall(PetscViewerASCIIPopTab(viewer)); // KDTree:
505   PetscFunctionReturn(PETSC_SUCCESS);
506 }
507