16c6e6638SPierre Jolivet #include <petscviewer.h> 298480730SJames Wright #include <petscis.h> 398480730SJames Wright #include <petsc/private/petscimpl.h> 498480730SJames Wright 56c6e6638SPierre Jolivet // For accessing bitwise Boolean values in are_handles_leaves 68279f6d9SJames Wright #define GREATER_BIT 0 78279f6d9SJames Wright #define LESS_EQUAL_BIT 1 88279f6d9SJames Wright 998480730SJames Wright typedef struct { 10c87b6f5fSJames Wright uint8_t axis; // Coordinate direction that stem splits on 11*73bc78fdSLisandro Dalcin PetscByte are_handles_leaves; // Bit-wise boolean for whether the greater_handle and less_equal_handle index kdtree->stems or kstree->leaves 12c87b6f5fSJames Wright PetscReal split; // Coordinate value that stem splits on 13c87b6f5fSJames Wright PetscCount greater_handle, less_equal_handle; // Handle (index) of the child nodes of the tree 1498480730SJames Wright } KDStem; 1598480730SJames Wright 1698480730SJames Wright typedef struct { 17c87b6f5fSJames Wright PetscInt count; // Number of points in leaf 18c87b6f5fSJames Wright PetscCount indices_handle, coords_handle; // Index into the indices or coordinates array for the points in leaf 1998480730SJames Wright } KDLeaf; 2098480730SJames Wright 2198480730SJames Wright struct _n_PetscKDTree { 22c87b6f5fSJames Wright PetscInt dim; // Coordinate dimension of the tree 23c87b6f5fSJames Wright PetscInt max_bucket_size; // Maximum number of points stored at each leaf 2498480730SJames Wright 2598480730SJames Wright PetscBool is_root_leaf; 2698480730SJames Wright PetscCount root_handle; 2798480730SJames Wright 2898480730SJames Wright PetscCount num_coords, num_leaves, num_stems, num_bucket_indices; 2998480730SJames Wright const PetscReal *coords, *coords_owned; // Only free owned on Destroy 3098480730SJames Wright KDLeaf *leaves; 3198480730SJames Wright KDStem *stems; 3298480730SJames Wright PetscCount *bucket_indices; 3398480730SJames Wright }; 3498480730SJames Wright 3598480730SJames Wright /*@C 3698480730SJames Wright PetscKDTreeDestroy - destroy a `PetscKDTree` 3798480730SJames Wright 3898480730SJames Wright Not Collective, No Fortran Support 3998480730SJames Wright 4098480730SJames Wright Input Parameters: 4198480730SJames Wright . tree - tree to destroy 4298480730SJames Wright 4398480730SJames Wright Level: advanced 4498480730SJames Wright 4598480730SJames Wright .seealso: `PetscKDTree`, `PetscKDTreeCreate()` 4698480730SJames Wright @*/ 4798480730SJames Wright PetscErrorCode PetscKDTreeDestroy(PetscKDTree *tree) 4898480730SJames Wright { 4998480730SJames Wright PetscFunctionBeginUser; 5098480730SJames Wright if (*tree == NULL) PetscFunctionReturn(PETSC_SUCCESS); 5198480730SJames Wright PetscCall(PetscFree((*tree)->stems)); 5298480730SJames Wright PetscCall(PetscFree((*tree)->leaves)); 5398480730SJames Wright PetscCall(PetscFree((*tree)->bucket_indices)); 5498480730SJames Wright PetscCall(PetscFree((*tree)->coords_owned)); 5598480730SJames Wright PetscCall(PetscFree(*tree)); 5698480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5798480730SJames Wright } 5898480730SJames Wright 5998480730SJames Wright PetscLogEvent PetscKDTree_Build, PetscKDTree_Query; 60eba205daSPierre Jolivet static PetscErrorCode PetscKDTreeRegisterLogEvents(void) 6198480730SJames Wright { 6298480730SJames Wright static PetscBool is_initialized = PETSC_FALSE; 6398480730SJames Wright 6498480730SJames Wright PetscFunctionBeginUser; 6598480730SJames Wright if (is_initialized) PetscFunctionReturn(PETSC_SUCCESS); 6698480730SJames Wright PetscCall(PetscLogEventRegister("KDTreeBuild", IS_CLASSID, &PetscKDTree_Build)); 6798480730SJames Wright PetscCall(PetscLogEventRegister("KDTreeQuery", IS_CLASSID, &PetscKDTree_Query)); 6898480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 6998480730SJames Wright } 7098480730SJames Wright 7198480730SJames Wright // From http://graphics.stanford.edu/~seander/bithacks.html#RoundUpPowerOf2 7298480730SJames Wright static inline uint32_t RoundToNextPowerOfTwo(uint32_t v) 7398480730SJames Wright { 7498480730SJames Wright v--; 7598480730SJames Wright v |= v >> 1; 7698480730SJames Wright v |= v >> 2; 7798480730SJames Wright v |= v >> 4; 7898480730SJames Wright v |= v >> 8; 7998480730SJames Wright v |= v >> 16; 8098480730SJames Wright v++; 8198480730SJames Wright return v; 8298480730SJames Wright } 8398480730SJames Wright 8498480730SJames Wright typedef struct { 858279f6d9SJames Wright uint8_t initial_axis; 8698480730SJames Wright PetscKDTree tree; 8798480730SJames Wright } *KDTreeSortContext; 8898480730SJames Wright 8998480730SJames Wright // Sort nodes based on "superkey" 9098480730SJames Wright // See "Building a Balanced k-d Tree in O(kn log n) Time" https://jcgt.org/published/0004/01/03/ 918279f6d9SJames Wright static inline int PetscKDTreeSortFunc(PetscCount left, PetscCount right, PetscKDTree tree, uint8_t axis) 9298480730SJames Wright { 9398480730SJames Wright const PetscReal *coords = tree->coords; 9498480730SJames Wright const PetscInt dim = tree->dim; 9598480730SJames Wright 9698480730SJames Wright for (PetscInt i = 0; i < dim; i++) { 9798480730SJames Wright PetscReal diff = coords[left * dim + axis] - coords[right * dim + axis]; 9898480730SJames Wright if (PetscUnlikely(diff == 0)) { 9998480730SJames Wright axis = (axis + 1) % dim; 10098480730SJames Wright continue; 10198480730SJames Wright } else return PetscSign(diff); 10298480730SJames Wright } 10398480730SJames Wright return 0; // All components are the same 10498480730SJames Wright } 10598480730SJames Wright 1062a8381b2SBarry Smith static int PetscKDTreeTimSort(const void *l, const void *r, PetscCtx ctx) 10798480730SJames Wright { 10898480730SJames Wright KDTreeSortContext kd_ctx = (KDTreeSortContext)ctx; 10998480730SJames Wright return PetscKDTreeSortFunc(*(PetscCount *)l, *(PetscCount *)r, kd_ctx->tree, kd_ctx->initial_axis); 11098480730SJames Wright } 11198480730SJames Wright 11298480730SJames Wright static PetscErrorCode PetscKDTreeVerifySortedIndices(PetscKDTree tree, PetscCount sorted_indices[], PetscCount temp[], PetscCount start, PetscCount end) 11398480730SJames Wright { 11498480730SJames Wright PetscCount num_coords = tree->num_coords, range_size = end - start, location; 11598480730SJames Wright PetscBool has_duplicates; 11698480730SJames Wright 11798480730SJames Wright PetscFunctionBeginUser; 11898480730SJames Wright PetscCall(PetscArraycpy(temp, &sorted_indices[0 * num_coords + start], range_size)); 11998480730SJames Wright PetscCall(PetscSortCount(range_size, temp)); 12098480730SJames Wright PetscCall(PetscSortedCheckDupsCount(range_size, temp, &has_duplicates)); 12198480730SJames Wright PetscCheck(has_duplicates == PETSC_FALSE, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Sorted indices must have unique entries, but found duplicates"); 12298480730SJames Wright for (PetscInt d = 1; d < tree->dim; d++) { 12398480730SJames Wright for (PetscCount i = start; i < end; i++) { 12498480730SJames Wright PetscCall(PetscFindCount(sorted_indices[d * num_coords + i], range_size, temp, &location)); 12598480730SJames Wright 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); 12698480730SJames Wright } 12798480730SJames Wright } 12898480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 12998480730SJames Wright } 13098480730SJames Wright 13198480730SJames Wright typedef struct { 13298480730SJames Wright PetscKDTree tree; 13398480730SJames Wright PetscSegBuffer stems, leaves, bucket_indices, bucket_coords; 13498480730SJames Wright PetscBool debug_build, copy_coords; 13598480730SJames Wright } *KDTreeBuild; 13698480730SJames Wright 13798480730SJames Wright // The range is end exclusive, so [start,end). 13898480730SJames Wright static PetscErrorCode PetscKDTreeBuildStemAndLeaves(KDTreeBuild kd_build, PetscCount sorted_indices[], PetscCount temp[], PetscCount start, PetscCount end, PetscInt depth, PetscBool *is_node_leaf, PetscCount *node_handle) 13998480730SJames Wright { 14098480730SJames Wright PetscKDTree tree = kd_build->tree; 14198480730SJames Wright PetscInt dim = tree->dim; 14298480730SJames Wright 14398480730SJames Wright PetscFunctionBeginUser; 14498480730SJames Wright PetscCheck(start <= end, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Start %" PetscCount_FMT " must be less than or equal to end %" PetscCount_FMT, start, end); 14598480730SJames Wright if (kd_build->debug_build) PetscCall(PetscKDTreeVerifySortedIndices(tree, sorted_indices, temp, start, end)); 14698480730SJames Wright if (end - start <= tree->max_bucket_size) { 14798480730SJames Wright KDLeaf *leaf; 14898480730SJames Wright PetscCount *bucket_indices; 14998480730SJames Wright 15098480730SJames Wright PetscCall(PetscSegBufferGetSize(kd_build->leaves, node_handle)); 15198480730SJames Wright PetscCall(PetscSegBufferGet(kd_build->leaves, 1, &leaf)); 15298480730SJames Wright PetscCall(PetscMemzero(leaf, sizeof(KDLeaf))); 15398480730SJames Wright *is_node_leaf = PETSC_TRUE; 15498480730SJames Wright 15598480730SJames Wright PetscCall(PetscIntCast(end - start, &leaf->count)); 15698480730SJames Wright PetscCall(PetscSegBufferGetSize(kd_build->bucket_indices, &leaf->indices_handle)); 15798480730SJames Wright PetscCall(PetscSegBufferGet(kd_build->bucket_indices, leaf->count, &bucket_indices)); 15898480730SJames Wright PetscCall(PetscArraycpy(bucket_indices, &sorted_indices[start], leaf->count)); 15998480730SJames Wright if (kd_build->copy_coords) { 16098480730SJames Wright PetscReal *bucket_coords; 16198480730SJames Wright PetscCall(PetscSegBufferGetSize(kd_build->bucket_coords, &leaf->coords_handle)); 16298480730SJames Wright PetscCall(PetscSegBufferGet(kd_build->bucket_coords, leaf->count * dim, &bucket_coords)); 16398480730SJames Wright // Coords are saved in axis-major ordering for better vectorization 16498480730SJames Wright for (PetscCount i = 0; i < leaf->count; i++) { 16598480730SJames Wright for (PetscInt d = 0; d < dim; d++) bucket_coords[d * leaf->count + i] = tree->coords[bucket_indices[i] * dim + d]; 16698480730SJames Wright } 16798480730SJames Wright } else leaf->coords_handle = -1; 16898480730SJames Wright } else { 16998480730SJames Wright KDStem *stem; 17098480730SJames Wright PetscCount num_coords = tree->num_coords; 1718279f6d9SJames Wright uint8_t axis = (uint8_t)(depth % dim); 1728279f6d9SJames Wright PetscBool is_greater_leaf, is_less_equal_leaf; 17398480730SJames Wright PetscCount median = start + PetscCeilInt64(end - start, 2) - 1, lower; 17498480730SJames Wright PetscCount median_idx = sorted_indices[median], medianp1_idx = sorted_indices[median + 1]; 17598480730SJames Wright 17698480730SJames Wright PetscCall(PetscSegBufferGetSize(kd_build->stems, node_handle)); 17798480730SJames Wright PetscCall(PetscSegBufferGet(kd_build->stems, 1, &stem)); 17898480730SJames Wright PetscCall(PetscMemzero(stem, sizeof(KDStem))); 17998480730SJames Wright *is_node_leaf = PETSC_FALSE; 18098480730SJames Wright 18198480730SJames Wright stem->axis = axis; 18298480730SJames Wright // Place split halfway between the "boundary" nodes of the partitioning 18398480730SJames Wright stem->split = (tree->coords[tree->dim * median_idx + axis] + tree->coords[tree->dim * medianp1_idx + axis]) / 2; 18498480730SJames Wright PetscCall(PetscArraycpy(temp, &sorted_indices[0 * num_coords + start], end - start)); 18598480730SJames Wright lower = median; // Set lower in case dim == 1 18698480730SJames Wright for (PetscInt d = 1; d < tree->dim; d++) { 18798480730SJames Wright PetscCount upper = median; 18898480730SJames Wright lower = start - 1; 18998480730SJames Wright for (PetscCount i = start; i < end; i++) { 19098480730SJames Wright // In case of duplicate coord point equal to the median coord point, limit lower partition to median, ensuring balanced tree 19198480730SJames Wright if (lower < median && PetscKDTreeSortFunc(sorted_indices[d * num_coords + i], median_idx, tree, axis) <= 0) { 19298480730SJames Wright sorted_indices[(d - 1) * num_coords + (++lower)] = sorted_indices[d * num_coords + i]; 19398480730SJames Wright } else { 19498480730SJames Wright sorted_indices[(d - 1) * num_coords + (++upper)] = sorted_indices[d * num_coords + i]; 19598480730SJames Wright } 19698480730SJames Wright } 19798480730SJames Wright 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); 19898480730SJames Wright 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); 19998480730SJames Wright } 20098480730SJames Wright PetscCall(PetscArraycpy(&sorted_indices[(tree->dim - 1) * num_coords + start], temp, end - start)); 20198480730SJames Wright 2028279f6d9SJames Wright PetscCall(PetscKDTreeBuildStemAndLeaves(kd_build, sorted_indices, temp, start, lower + 1, depth + 1, &is_less_equal_leaf, &stem->less_equal_handle)); 2038279f6d9SJames Wright if (is_less_equal_leaf) PetscCall(PetscBTSet(&stem->are_handles_leaves, LESS_EQUAL_BIT)); 2048279f6d9SJames Wright PetscCall(PetscKDTreeBuildStemAndLeaves(kd_build, sorted_indices, temp, lower + 1, end, depth + 1, &is_greater_leaf, &stem->greater_handle)); 2058279f6d9SJames Wright if (is_greater_leaf) PetscCall(PetscBTSet(&stem->are_handles_leaves, GREATER_BIT)); 20698480730SJames Wright } 20798480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 20898480730SJames Wright } 20998480730SJames Wright 21098480730SJames Wright /*@C 21198480730SJames Wright PetscKDTreeCreate - create a `PetscKDTree` 21298480730SJames Wright 21398480730SJames Wright Not Collective, No Fortran Support 21498480730SJames Wright 21598480730SJames Wright Input Parameters: 21698480730SJames Wright + num_coords - number of coordinate points to build the `PetscKDTree` 21798480730SJames Wright . dim - the dimension of the coordinates 21898480730SJames Wright . coords - array of the coordinates, in point-major order 21998480730SJames Wright . copy_mode - behavior handling `coords`, `PETSC_COPY_VALUES` generally more performant 22098480730SJames Wright - max_bucket_size - maximum number of points stored at each leaf 22198480730SJames Wright 22298480730SJames Wright Output Parameter: 22398480730SJames Wright . new_tree - the resulting `PetscKDTree` 22498480730SJames Wright 22598480730SJames Wright Level: advanced 22698480730SJames Wright 22798480730SJames Wright Note: 22898480730SJames Wright When `copy_mode == PETSC_COPY_VALUES`, the coordinates are copied and organized to optimize vectorization and cache-coherency. 22998480730SJames Wright 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. 23098480730SJames Wright 23198480730SJames Wright Developer Note: 23298480730SJames Wright Building algorithm detailed in 'Building a Balanced k-d Tree in O(kn log n) Time' Brown, 2015 23398480730SJames Wright 23498480730SJames Wright .seealso: `PetscKDTree`, `PetscKDTreeDestroy()`, `PetscKDTreeQueryPointsNearestNeighbor()` 23598480730SJames Wright @*/ 23698480730SJames Wright PetscErrorCode PetscKDTreeCreate(PetscCount num_coords, PetscInt dim, const PetscReal coords[], PetscCopyMode copy_mode, PetscInt max_bucket_size, PetscKDTree *new_tree) 23798480730SJames Wright { 23898480730SJames Wright PetscKDTree tree; 23998480730SJames Wright PetscCount *sorted_indices, *temp; 24098480730SJames Wright 24198480730SJames Wright PetscFunctionBeginUser; 24298480730SJames Wright PetscCall(PetscKDTreeRegisterLogEvents()); 24398480730SJames Wright PetscCall(PetscLogEventBegin(PetscKDTree_Build, 0, 0, 0, 0)); 244bfe80ac4SPierre Jolivet PetscCheck(dim > 0, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "Dimension of PetscKDTree must be greater than 0, received %" PetscInt_FMT, dim); 245bfe80ac4SPierre Jolivet PetscCheck(num_coords > -1, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "Number of coordinates may not be negative, received %" PetscCount_FMT, num_coords); 24698480730SJames Wright if (num_coords == 0) { 24798480730SJames Wright *new_tree = NULL; 24898480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 24998480730SJames Wright } 25098480730SJames Wright PetscAssertPointer(coords, 3); 25198480730SJames Wright PetscAssertPointer(new_tree, 6); 25298480730SJames Wright PetscCall(PetscNew(&tree)); 25398480730SJames Wright tree->dim = dim; 25498480730SJames Wright tree->max_bucket_size = max_bucket_size == PETSC_DECIDE ? 32 : max_bucket_size; 25598480730SJames Wright tree->num_coords = num_coords; 25698480730SJames Wright 25798480730SJames Wright switch (copy_mode) { 25898480730SJames Wright case PETSC_OWN_POINTER: 259565897f2SJames Wright tree->coords_owned = coords; // fallthrough 26098480730SJames Wright case PETSC_USE_POINTER: 26198480730SJames Wright tree->coords = coords; 26298480730SJames Wright break; 26398480730SJames Wright case PETSC_COPY_VALUES: 26498480730SJames Wright PetscCall(PetscMalloc1(num_coords * dim, &tree->coords_owned)); 26598480730SJames Wright PetscCall(PetscArraycpy((PetscReal *)tree->coords_owned, coords, num_coords * dim)); 26698480730SJames Wright tree->coords = tree->coords_owned; 26798480730SJames Wright break; 26898480730SJames Wright } 26998480730SJames Wright 27098480730SJames Wright KDTreeSortContext kd_ctx; 27198480730SJames Wright PetscCall(PetscMalloc2(num_coords * dim, &sorted_indices, num_coords, &temp)); 27298480730SJames Wright PetscCall(PetscNew(&kd_ctx)); 27398480730SJames Wright kd_ctx->tree = tree; 27498480730SJames Wright for (PetscInt j = 0; j < dim; j++) { 27598480730SJames Wright for (PetscCount i = 0; i < num_coords; i++) sorted_indices[num_coords * j + i] = i; 2768279f6d9SJames Wright kd_ctx->initial_axis = (uint8_t)j; 27798480730SJames Wright PetscCall(PetscTimSort((PetscInt)num_coords, &sorted_indices[num_coords * j], sizeof(*sorted_indices), PetscKDTreeTimSort, kd_ctx)); 27898480730SJames Wright } 27998480730SJames Wright PetscCall(PetscFree(kd_ctx)); 28098480730SJames Wright 28198480730SJames Wright PetscInt num_leaves = (PetscInt)PetscCeilInt64(num_coords, tree->max_bucket_size); 28298480730SJames Wright PetscInt num_stems = RoundToNextPowerOfTwo((uint32_t)num_leaves); 28398480730SJames Wright KDTreeBuild kd_build; 28498480730SJames Wright PetscCall(PetscNew(&kd_build)); 28598480730SJames Wright kd_build->tree = tree; 28698480730SJames Wright kd_build->copy_coords = copy_mode == PETSC_COPY_VALUES ? PETSC_TRUE : PETSC_FALSE; 28798480730SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-kdtree_debug", &kd_build->debug_build, NULL)); 28898480730SJames Wright PetscCall(PetscSegBufferCreate(sizeof(KDStem), num_stems, &kd_build->stems)); 28998480730SJames Wright PetscCall(PetscSegBufferCreate(sizeof(KDLeaf), num_leaves, &kd_build->leaves)); 29098480730SJames Wright PetscCall(PetscSegBufferCreate(sizeof(PetscCount), num_coords, &kd_build->bucket_indices)); 29198480730SJames Wright if (kd_build->copy_coords) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), num_coords * dim, &kd_build->bucket_coords)); 29298480730SJames Wright 29398480730SJames Wright PetscCall(PetscKDTreeBuildStemAndLeaves(kd_build, sorted_indices, temp, 0, num_coords, 0, &tree->is_root_leaf, &tree->root_handle)); 29498480730SJames Wright 29598480730SJames Wright PetscCall(PetscSegBufferGetSize(kd_build->stems, &tree->num_stems)); 29698480730SJames Wright PetscCall(PetscSegBufferGetSize(kd_build->leaves, &tree->num_leaves)); 29798480730SJames Wright PetscCall(PetscSegBufferGetSize(kd_build->bucket_indices, &tree->num_bucket_indices)); 29898480730SJames Wright PetscCall(PetscSegBufferExtractAlloc(kd_build->stems, &tree->stems)); 29998480730SJames Wright PetscCall(PetscSegBufferExtractAlloc(kd_build->leaves, &tree->leaves)); 30098480730SJames Wright PetscCall(PetscSegBufferExtractAlloc(kd_build->bucket_indices, &tree->bucket_indices)); 30198480730SJames Wright if (kd_build->copy_coords) { 30298480730SJames Wright PetscCall(PetscFree(tree->coords_owned)); 30398480730SJames Wright PetscCall(PetscSegBufferExtractAlloc(kd_build->bucket_coords, &tree->coords_owned)); 30498480730SJames Wright tree->coords = tree->coords_owned; 30598480730SJames Wright PetscCall(PetscSegBufferDestroy(&kd_build->bucket_coords)); 30698480730SJames Wright } 30798480730SJames Wright PetscCall(PetscSegBufferDestroy(&kd_build->stems)); 30898480730SJames Wright PetscCall(PetscSegBufferDestroy(&kd_build->leaves)); 30998480730SJames Wright PetscCall(PetscSegBufferDestroy(&kd_build->bucket_indices)); 31098480730SJames Wright PetscCall(PetscFree(kd_build)); 31198480730SJames Wright PetscCall(PetscFree2(sorted_indices, temp)); 31298480730SJames Wright *new_tree = tree; 31398480730SJames Wright PetscCall(PetscLogEventEnd(PetscKDTree_Build, 0, 0, 0, 0)); 31498480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 31598480730SJames Wright } 31698480730SJames Wright 31798480730SJames Wright static inline PetscReal PetscSquareDistance(PetscInt dim, const PetscReal *PETSC_RESTRICT x, const PetscReal *PETSC_RESTRICT y) 31898480730SJames Wright { 31998480730SJames Wright PetscReal dist = 0; 32098480730SJames Wright for (PetscInt j = 0; j < dim; j++) dist += PetscSqr(x[j] - y[j]); 32198480730SJames Wright return dist; 32298480730SJames Wright } 32398480730SJames Wright 32498480730SJames Wright static inline PetscErrorCode PetscKDTreeQueryLeaf(PetscKDTree tree, KDLeaf leaf, const PetscReal point[], PetscCount *index, PetscReal *distance_sqr) 32598480730SJames Wright { 32698480730SJames Wright PetscInt dim = tree->dim; 32798480730SJames Wright 32898480730SJames Wright PetscFunctionBeginUser; 32998480730SJames Wright *distance_sqr = PETSC_MAX_REAL; 33098480730SJames Wright *index = -1; 33198480730SJames Wright for (PetscInt i = 0; i < leaf.count; i++) { 33298480730SJames Wright PetscCount point_index = tree->bucket_indices[leaf.indices_handle + i]; 33398480730SJames Wright PetscReal dist = PetscSquareDistance(dim, point, &tree->coords[point_index * dim]); 33498480730SJames Wright if (dist < *distance_sqr) { 33598480730SJames Wright *distance_sqr = dist; 33698480730SJames Wright *index = point_index; 33798480730SJames Wright } 33898480730SJames Wright } 33998480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 34098480730SJames Wright } 34198480730SJames Wright 34298480730SJames Wright static inline PetscErrorCode PetscKDTreeQueryLeaf_CopyCoords(PetscKDTree tree, KDLeaf leaf, const PetscReal point[], PetscCount *index, PetscReal *distance_sqr) 34398480730SJames Wright { 34498480730SJames Wright PetscInt dim = tree->dim; 34598480730SJames Wright 34698480730SJames Wright PetscFunctionBeginUser; 34798480730SJames Wright *distance_sqr = PETSC_MAX_REAL; 34898480730SJames Wright *index = -1; 34998480730SJames Wright for (PetscInt i = 0; i < leaf.count; i++) { 35098480730SJames Wright // Coord data saved in axis-major ordering for vectorization 35198480730SJames Wright PetscReal dist = 0.; 35298480730SJames Wright for (PetscInt d = 0; d < dim; d++) dist += PetscSqr(point[d] - tree->coords[leaf.coords_handle + d * leaf.count + i]); 35398480730SJames Wright if (dist < *distance_sqr) { 35498480730SJames Wright *distance_sqr = dist; 35598480730SJames Wright *index = tree->bucket_indices[leaf.indices_handle + i]; 35698480730SJames Wright } 35798480730SJames Wright } 35898480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 35998480730SJames Wright } 36098480730SJames Wright 36198480730SJames Wright // Recursive point query from 'Algorithms for Fast Vector Quantization' by Sunil Arya and David Mount 36298480730SJames Wright // Variant also implemented in pykdtree 3638279f6d9SJames Wright static PetscErrorCode PetscKDTreeQuery_Recurse(PetscKDTree tree, const PetscReal point[], PetscCount node_handle, char is_node_leaf, PetscReal offset[], PetscReal rd, PetscReal tol_sqr, PetscCount *index, PetscReal *dist_sqr) 36498480730SJames Wright { 36598480730SJames Wright PetscFunctionBeginUser; 36698480730SJames Wright if (*dist_sqr < tol_sqr) PetscFunctionReturn(PETSC_SUCCESS); 36798480730SJames Wright if (is_node_leaf) { 36898480730SJames Wright KDLeaf leaf = tree->leaves[node_handle]; 36998480730SJames Wright PetscReal dist; 37098480730SJames Wright PetscCount point_index; 37198480730SJames Wright 37298480730SJames Wright if (leaf.coords_handle > -1) PetscCall(PetscKDTreeQueryLeaf_CopyCoords(tree, leaf, point, &point_index, &dist)); 37398480730SJames Wright else PetscCall(PetscKDTreeQueryLeaf(tree, leaf, point, &point_index, &dist)); 37498480730SJames Wright if (dist < *dist_sqr) { 37598480730SJames Wright *dist_sqr = dist; 37698480730SJames Wright *index = point_index; 37798480730SJames Wright } 37898480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 37998480730SJames Wright } 38098480730SJames Wright 38198480730SJames Wright KDStem stem = tree->stems[node_handle]; 38298480730SJames Wright PetscReal old_offset = offset[stem.axis], new_offset = point[stem.axis] - stem.split; 38398480730SJames Wright if (new_offset <= 0) { 3848279f6d9SJames Wright PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.less_equal_handle, PetscBTLookup(&stem.are_handles_leaves, LESS_EQUAL_BIT), offset, rd, tol_sqr, index, dist_sqr)); 38598480730SJames Wright rd += -PetscSqr(old_offset) + PetscSqr(new_offset); 38698480730SJames Wright if (rd < *dist_sqr) { 38798480730SJames Wright offset[stem.axis] = new_offset; 3888279f6d9SJames Wright PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.greater_handle, PetscBTLookup(&stem.are_handles_leaves, GREATER_BIT), offset, rd, tol_sqr, index, dist_sqr)); 38998480730SJames Wright offset[stem.axis] = old_offset; 39098480730SJames Wright } 39198480730SJames Wright } else { 3928279f6d9SJames Wright PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.greater_handle, PetscBTLookup(&stem.are_handles_leaves, GREATER_BIT), offset, rd, tol_sqr, index, dist_sqr)); 39398480730SJames Wright rd += -PetscSqr(old_offset) + PetscSqr(new_offset); 39498480730SJames Wright if (rd < *dist_sqr) { 39598480730SJames Wright offset[stem.axis] = new_offset; 3968279f6d9SJames Wright PetscCall(PetscKDTreeQuery_Recurse(tree, point, stem.less_equal_handle, PetscBTLookup(&stem.are_handles_leaves, LESS_EQUAL_BIT), offset, rd, tol_sqr, index, dist_sqr)); 39798480730SJames Wright offset[stem.axis] = old_offset; 39898480730SJames Wright } 39998480730SJames Wright } 40098480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 40198480730SJames Wright } 40298480730SJames Wright 40398480730SJames Wright /*@C 40498480730SJames Wright PetscKDTreeQueryPointsNearestNeighbor - find the nearest neighbor in a `PetscKDTree` 40598480730SJames Wright 40698480730SJames Wright Not Collective, No Fortran Support 40798480730SJames Wright 40898480730SJames Wright Input Parameters: 40998480730SJames Wright + tree - tree to query 41098480730SJames Wright . num_points - number of points to query 41198480730SJames Wright . points - array of the coordinates, in point-major order 41298480730SJames Wright - tolerance - tolerance for nearest neighbor 41398480730SJames Wright 414d6ae5217SJose E. Roman Output Parameters: 41598480730SJames Wright + indices - indices of the nearest neighbor to the query point 41698480730SJames Wright - distances - distance between the queried point and the nearest neighbor 41798480730SJames Wright 41898480730SJames Wright Level: advanced 41998480730SJames Wright 42098480730SJames Wright Notes: 42198480730SJames Wright 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. 42298480730SJames Wright 42398480730SJames Wright The `indices` and `distances` arrays should be at least of size `num_points`. 42498480730SJames Wright 42598480730SJames Wright .seealso: `PetscKDTree`, `PetscKDTreeCreate()` 42698480730SJames Wright @*/ 42798480730SJames Wright PetscErrorCode PetscKDTreeQueryPointsNearestNeighbor(PetscKDTree tree, PetscCount num_points, const PetscReal points[], PetscReal tolerance, PetscCount indices[], PetscReal distances[]) 42898480730SJames Wright { 42998480730SJames Wright PetscReal *offsets, rd; 43098480730SJames Wright 43198480730SJames Wright PetscFunctionBeginUser; 43298480730SJames Wright PetscCall(PetscLogEventBegin(PetscKDTree_Query, 0, 0, 0, 0)); 43398480730SJames Wright if (tree == NULL) { 43498480730SJames Wright PetscCheck(num_points == 0, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "num_points may only be zero, if tree is NULL"); 43598480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 43698480730SJames Wright } 43798480730SJames Wright PetscAssertPointer(points, 3); 43898480730SJames Wright PetscAssertPointer(indices, 5); 43998480730SJames Wright PetscAssertPointer(distances, 6); 44098480730SJames Wright PetscCall(PetscCalloc1(tree->dim, &offsets)); 44198480730SJames Wright 44298480730SJames Wright for (PetscCount p = 0; p < num_points; p++) { 44398480730SJames Wright rd = 0.; 44498480730SJames Wright distances[p] = PETSC_MAX_REAL; 44598480730SJames Wright indices[p] = -1; 4468279f6d9SJames Wright PetscCall(PetscKDTreeQuery_Recurse(tree, &points[p * tree->dim], tree->root_handle, (char)tree->is_root_leaf, offsets, rd, PetscSqr(tolerance), &indices[p], &distances[p])); 44798480730SJames Wright distances[p] = PetscSqrtReal(distances[p]); 44898480730SJames Wright } 44998480730SJames Wright PetscCall(PetscFree(offsets)); 45098480730SJames Wright PetscCall(PetscLogEventEnd(PetscKDTree_Query, 0, 0, 0, 0)); 45198480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 45298480730SJames Wright } 45398480730SJames Wright 45498480730SJames Wright /*@C 45598480730SJames Wright PetscKDTreeView - view a `PetscKDTree` 45698480730SJames Wright 45798480730SJames Wright Not Collective, No Fortran Support 45898480730SJames Wright 45998480730SJames Wright Input Parameters: 46098480730SJames Wright + tree - tree to view 46198480730SJames Wright - viewer - visualization context 46298480730SJames Wright 46398480730SJames Wright Level: advanced 46498480730SJames Wright 46598480730SJames Wright .seealso: `PetscKDTree`, `PetscKDTreeCreate()`, `PetscViewer` 46698480730SJames Wright @*/ 46798480730SJames Wright PetscErrorCode PetscKDTreeView(PetscKDTree tree, PetscViewer viewer) 46898480730SJames Wright { 46998480730SJames Wright PetscFunctionBeginUser; 47098480730SJames Wright if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 47198480730SJames Wright else PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &viewer)); 47298480730SJames Wright if (tree == NULL) PetscFunctionReturn(PETSC_SUCCESS); 47398480730SJames Wright 47498480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "KDTree:\n")); 47598480730SJames Wright PetscCall(PetscViewerASCIIPushTab(viewer)); // KDTree: 47698480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "Stems:\n")); 47798480730SJames Wright PetscCall(PetscViewerASCIIPushTab(viewer)); // Stems: 47898480730SJames Wright for (PetscCount i = 0; i < tree->num_stems; i++) { 47998480730SJames Wright KDStem stem = tree->stems[i]; 4808279f6d9SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "Stem %" PetscCount_FMT ": Axis=%" PetscInt_FMT " Split=%g Greater_%s=%" PetscCount_FMT " Lesser_Equal_%s=%" PetscCount_FMT "\n", i, (PetscInt)stem.axis, (double)stem.split, PetscBTLookup(&stem.are_handles_leaves, GREATER_BIT) ? "Leaf" : "Stem", 4818279f6d9SJames Wright stem.greater_handle, PetscBTLookup(&stem.are_handles_leaves, LESS_EQUAL_BIT) ? "Leaf" : "Stem", stem.less_equal_handle)); 48298480730SJames Wright } 48398480730SJames Wright PetscCall(PetscViewerASCIIPopTab(viewer)); // Stems: 48498480730SJames Wright 48598480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "Leaves:\n")); 48698480730SJames Wright PetscCall(PetscViewerASCIIPushTab(viewer)); // Leaves: 48798480730SJames Wright for (PetscCount i = 0; i < tree->num_leaves; i++) { 48898480730SJames Wright KDLeaf leaf = tree->leaves[i]; 48998480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "Leaf %" PetscCount_FMT ": Count=%" PetscInt_FMT, i, leaf.count)); 49098480730SJames Wright PetscCall(PetscViewerASCIIPushTab(viewer)); // Coords: 49198480730SJames Wright for (PetscInt j = 0; j < leaf.count; j++) { 49298480730SJames Wright PetscInt tabs; 49398480730SJames Wright PetscCount bucket_index = tree->bucket_indices[leaf.indices_handle + j]; 49498480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 49598480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscCount_FMT ": ", bucket_index)); 49698480730SJames Wright 49798480730SJames Wright PetscCall(PetscViewerASCIIGetTab(viewer, &tabs)); 49898480730SJames Wright PetscCall(PetscViewerASCIISetTab(viewer, 0)); 49998480730SJames Wright if (leaf.coords_handle > -1) { 50098480730SJames Wright for (PetscInt k = 0; k < tree->dim; k++) PetscCall(PetscViewerASCIIPrintf(viewer, "%g ", (double)tree->coords[leaf.coords_handle + leaf.count * k + j])); 50198480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, " (stored at leaf)")); 50298480730SJames Wright } else { 50398480730SJames Wright for (PetscInt k = 0; k < tree->dim; k++) PetscCall(PetscViewerASCIIPrintf(viewer, "%g ", (double)tree->coords[bucket_index * tree->dim + k])); 50498480730SJames Wright } 50598480730SJames Wright PetscCall(PetscViewerASCIISetTab(viewer, tabs)); 50698480730SJames Wright } 50798480730SJames Wright PetscCall(PetscViewerASCIIPopTab(viewer)); // Coords: 50898480730SJames Wright PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 50998480730SJames Wright } 51098480730SJames Wright PetscCall(PetscViewerASCIIPopTab(viewer)); // Leaves: 51198480730SJames Wright PetscCall(PetscViewerASCIIPopTab(viewer)); // KDTree: 51298480730SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 51398480730SJames Wright } 514