1*98480730SJames Wright static const char help[] = "Test KDTree\n\n"; 2*98480730SJames Wright 3*98480730SJames Wright #include <petsc.h> 4*98480730SJames Wright 5*98480730SJames Wright static inline PetscReal Distance(PetscInt dim, const PetscReal *PETSC_RESTRICT x, const PetscReal *PETSC_RESTRICT y) 6*98480730SJames Wright { 7*98480730SJames Wright PetscReal dist = 0; 8*98480730SJames Wright for (PetscInt j = 0; j < dim; j++) dist += PetscSqr(x[j] - y[j]); 9*98480730SJames Wright return PetscSqrtReal(dist); 10*98480730SJames Wright } 11*98480730SJames Wright 12*98480730SJames Wright int main(int argc, char **argv) 13*98480730SJames Wright { 14*98480730SJames Wright MPI_Comm comm; 15*98480730SJames Wright PetscInt num_coords, dim, num_rand_points = 0, bucket_size = PETSC_DECIDE; 16*98480730SJames Wright PetscRandom random; 17*98480730SJames Wright PetscReal *coords; 18*98480730SJames Wright PetscInt coords_size, num_points_queried = 0, num_trees_built = 0, loops = 1; 19*98480730SJames Wright PetscBool view_tree = PETSC_FALSE, view_performance = PETSC_TRUE, test_tree_points = PETSC_FALSE, test_rand_points = PETSC_FALSE, query_rand_points = PETSC_FALSE; 20*98480730SJames Wright PetscCopyMode copy_mode = PETSC_OWN_POINTER; 21*98480730SJames Wright PetscKDTree tree; 22*98480730SJames Wright 23*98480730SJames Wright PetscFunctionBeginUser; 24*98480730SJames Wright PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 25*98480730SJames Wright PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogEventGetPerfInfo without -log_view 26*98480730SJames Wright comm = PETSC_COMM_WORLD; 27*98480730SJames Wright 28*98480730SJames Wright PetscOptionsBegin(PETSC_COMM_WORLD, "", "Test Options", "none"); 29*98480730SJames Wright PetscCall(PetscOptionsInt("-num_coords", "Number of coordinates", "", PETSC_FALSE, &num_coords, NULL)); 30*98480730SJames Wright PetscCall(PetscOptionsInt("-dim", "Dimension of the coordinates", "", PETSC_FALSE, &dim, NULL)); 31*98480730SJames Wright PetscCall(PetscOptionsInt("-bucket_size", "Size of the bucket at each leaf", "", bucket_size, &bucket_size, NULL)); 32*98480730SJames Wright PetscCall(PetscOptionsBoundedInt("-loops", "Number of times to loop through KDTree creation and querying", "", loops, &loops, NULL, 1)); 33*98480730SJames Wright PetscCall(PetscOptionsEnum("-copy_mode", "Copy mode to use with KDTree", "PetscKDTreeCreate", PetscCopyModes, (PetscEnum)copy_mode, (PetscEnum *)©_mode, NULL)); 34*98480730SJames Wright PetscCall(PetscOptionsBool("-view_tree", "View the KDTree", "", view_tree, &view_tree, NULL)); 35*98480730SJames Wright PetscCall(PetscOptionsBool("-view_performance", "View the performance speed of the KDTree", "", view_performance, &view_performance, NULL)); 36*98480730SJames Wright PetscCall(PetscOptionsBool("-test_tree_points", "Test querying tree points against itself", "", test_tree_points, &test_tree_points, NULL)); 37*98480730SJames Wright PetscCall(PetscOptionsBool("-test_rand_points", "Test querying random points via brute force", "", test_rand_points, &test_rand_points, NULL)); 38*98480730SJames Wright PetscCall(PetscOptionsBool("-query_rand_points", "Query querying random points", "", query_rand_points, &query_rand_points, NULL)); 39*98480730SJames Wright if (test_rand_points || query_rand_points) PetscCall(PetscOptionsInt("-num_rand_points", "Number of random points to test with", "", num_rand_points, &num_rand_points, NULL)); 40*98480730SJames Wright PetscOptionsEnd(); 41*98480730SJames Wright 42*98480730SJames Wright coords_size = num_coords * dim; 43*98480730SJames Wright PetscCall(PetscMalloc1(coords_size, &coords)); 44*98480730SJames Wright PetscCall(PetscRandomCreate(comm, &random)); 45*98480730SJames Wright 46*98480730SJames Wright for (PetscInt loop_count = 0; loop_count < loops; loop_count++) { 47*98480730SJames Wright PetscCall(PetscRandomGetValuesReal(random, coords_size, coords)); 48*98480730SJames Wright 49*98480730SJames Wright PetscCall(PetscKDTreeCreate(num_coords, dim, coords, copy_mode, bucket_size, &tree)); 50*98480730SJames Wright num_trees_built++; 51*98480730SJames Wright if (view_tree) PetscCall(PetscKDTreeView(tree, NULL)); 52*98480730SJames Wright 53*98480730SJames Wright if (test_tree_points) { // Test querying the current coordinates 54*98480730SJames Wright PetscCount *indices; 55*98480730SJames Wright PetscReal *distances; 56*98480730SJames Wright num_points_queried += num_coords; 57*98480730SJames Wright 58*98480730SJames Wright PetscCall(PetscMalloc2(num_coords, &indices, num_coords, &distances)); 59*98480730SJames Wright PetscCall(PetscKDTreeQueryPointsNearestNeighbor(tree, num_coords, coords, PETSC_MACHINE_EPSILON * 1e3, indices, distances)); 60*98480730SJames Wright for (PetscInt i = 0; i < num_coords; i++) { 61*98480730SJames Wright if (indices[i] != i) PetscCall(PetscPrintf(comm, "Query failed for coord %" PetscInt_FMT ", query returned index %" PetscCount_FMT " with distance %g\n", i, indices[i], (double)distances[i])); 62*98480730SJames Wright } 63*98480730SJames Wright PetscCall(PetscFree2(indices, distances)); 64*98480730SJames Wright } 65*98480730SJames Wright 66*98480730SJames Wright if (num_rand_points > 0) { 67*98480730SJames Wright PetscCount *indices; 68*98480730SJames Wright PetscReal *distances; 69*98480730SJames Wright PetscReal *rand_points; 70*98480730SJames Wright PetscInt rand_queries_size = num_rand_points * dim; 71*98480730SJames Wright 72*98480730SJames Wright num_points_queried += num_rand_points; 73*98480730SJames Wright PetscCall(PetscMalloc3(rand_queries_size, &rand_points, num_rand_points, &indices, num_rand_points, &distances)); 74*98480730SJames Wright PetscCall(PetscRandomGetValuesReal(random, rand_queries_size, rand_points)); 75*98480730SJames Wright PetscCall(PetscKDTreeQueryPointsNearestNeighbor(tree, num_rand_points, rand_points, 0., indices, distances)); 76*98480730SJames Wright 77*98480730SJames Wright if (test_rand_points) { 78*98480730SJames Wright for (PetscInt i = 0; i < num_rand_points; i++) { 79*98480730SJames Wright PetscReal *rand_point = &rand_points[dim * i], nearest_distance = PETSC_MAX_REAL; 80*98480730SJames Wright PetscInt index = -1; 81*98480730SJames Wright for (PetscInt j = 0; j < num_coords; j++) { 82*98480730SJames Wright PetscReal dist = Distance(dim, rand_point, &coords[dim * j]); 83*98480730SJames Wright if (dist < nearest_distance) { 84*98480730SJames Wright nearest_distance = dist; 85*98480730SJames Wright index = j; 86*98480730SJames Wright } 87*98480730SJames Wright } 88*98480730SJames Wright if (indices[i] != index) 89*98480730SJames Wright PetscCall(PetscPrintf(comm, "Query failed for random point %" PetscInt_FMT ". Query returned index %" PetscCount_FMT " with distance %g, but coordinate %" PetscInt_FMT " with distance %g is closer\n", i, indices[i], (double)distances[i], index, (double)nearest_distance)); 90*98480730SJames Wright } 91*98480730SJames Wright } 92*98480730SJames Wright PetscCall(PetscFree3(rand_points, indices, distances)); 93*98480730SJames Wright } 94*98480730SJames Wright } 95*98480730SJames Wright 96*98480730SJames Wright if (view_performance) { 97*98480730SJames Wright PetscLogEvent kdtree_build_log, kdtree_query_log; 98*98480730SJames Wright PetscEventPerfInfo build_perf_info, query_perf_info; 99*98480730SJames Wright 100*98480730SJames Wright PetscCall(PetscLogEventGetId("KDTreeBuild", &kdtree_build_log)); 101*98480730SJames Wright PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, kdtree_build_log, &build_perf_info)); 102*98480730SJames Wright PetscCall(PetscLogEventGetId("KDTreeQuery", &kdtree_query_log)); 103*98480730SJames Wright PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, kdtree_query_log, &query_perf_info)); 104*98480730SJames Wright PetscCall(PetscPrintf(comm, "KDTreeBuild %.6e for %" PetscInt_FMT " trees built\n", build_perf_info.time, num_trees_built)); 105*98480730SJames Wright PetscCall(PetscPrintf(comm, "\tTime per tree: %.6e\n", build_perf_info.time / num_trees_built)); 106*98480730SJames Wright PetscCall(PetscPrintf(comm, "KDTreeQuery %.6e for %" PetscCount_FMT " queries\n", query_perf_info.time, (PetscCount)num_points_queried)); 107*98480730SJames Wright PetscCall(PetscPrintf(comm, "\tTime per query: %.6e\n", query_perf_info.time / num_points_queried)); 108*98480730SJames Wright } 109*98480730SJames Wright 110*98480730SJames Wright if (copy_mode != PETSC_OWN_POINTER) PetscCall(PetscFree(coords)); 111*98480730SJames Wright PetscCall(PetscKDTreeDestroy(&tree)); 112*98480730SJames Wright PetscCall(PetscRandomDestroy(&random)); 113*98480730SJames Wright PetscCall(PetscFinalize()); 114*98480730SJames Wright return 0; 115*98480730SJames Wright } 116*98480730SJames Wright 117*98480730SJames Wright /*TEST 118*98480730SJames Wright testset: 119*98480730SJames Wright suffix: kdtree 120*98480730SJames Wright args: -num_coords 35 -test_tree_points -test_rand_points -num_rand_points 300 -bucket_size 13 -view_performance false -view_tree true -kdtree_debug 121*98480730SJames Wright test: 122*98480730SJames Wright suffix: 1D 123*98480730SJames Wright args: -dim 1 124*98480730SJames Wright test: 125*98480730SJames Wright suffix: 2D 126*98480730SJames Wright args: -dim 2 127*98480730SJames Wright test: 128*98480730SJames Wright suffix: 3D 129*98480730SJames Wright args: -dim 3 130*98480730SJames Wright test: 131*98480730SJames Wright suffix: 3D_small 132*98480730SJames Wright args: -dim 3 -num_coords 13 133*98480730SJames Wright 134*98480730SJames Wright testset: 135*98480730SJames Wright suffix: kdtree_3D_large 136*98480730SJames Wright args: -dim 3 -num_coords 350 -test_tree_points -test_rand_points -num_rand_points 300 -view_performance false -kdtree_debug 137*98480730SJames Wright test: 138*98480730SJames Wright test: 139*98480730SJames Wright suffix: copy 140*98480730SJames Wright args: -copy_mode copy_values 141*98480730SJames Wright test: 142*98480730SJames Wright suffix: use 143*98480730SJames Wright args: -copy_mode use_pointer 144*98480730SJames Wright TEST*/ 145