198480730SJames Wright static const char help[] = "Test KDTree\n\n";
298480730SJames Wright
398480730SJames Wright #include <petsc.h>
498480730SJames Wright
Distance(PetscInt dim,const PetscReal * PETSC_RESTRICT x,const PetscReal * PETSC_RESTRICT y)598480730SJames Wright static inline PetscReal Distance(PetscInt dim, const PetscReal *PETSC_RESTRICT x, const PetscReal *PETSC_RESTRICT y)
698480730SJames Wright {
798480730SJames Wright PetscReal dist = 0;
898480730SJames Wright for (PetscInt j = 0; j < dim; j++) dist += PetscSqr(x[j] - y[j]);
998480730SJames Wright return PetscSqrtReal(dist);
1098480730SJames Wright }
1198480730SJames Wright
main(int argc,char ** argv)1298480730SJames Wright int main(int argc, char **argv)
1398480730SJames Wright {
1498480730SJames Wright MPI_Comm comm;
1598480730SJames Wright PetscInt num_coords, dim, num_rand_points = 0, bucket_size = PETSC_DECIDE;
1698480730SJames Wright PetscRandom random;
1798480730SJames Wright PetscReal *coords;
1898480730SJames Wright PetscInt coords_size, num_points_queried = 0, num_trees_built = 0, loops = 1;
1998480730SJames 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;
2098480730SJames Wright PetscCopyMode copy_mode = PETSC_OWN_POINTER;
2198480730SJames Wright PetscKDTree tree;
2298480730SJames Wright
2398480730SJames Wright PetscFunctionBeginUser;
2498480730SJames Wright PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2598480730SJames Wright PetscCall(PetscLogDefaultBegin()); // So we can use PetscLogEventGetPerfInfo without -log_view
2698480730SJames Wright comm = PETSC_COMM_WORLD;
2798480730SJames Wright
2898480730SJames Wright PetscOptionsBegin(PETSC_COMM_WORLD, "", "Test Options", "none");
2998480730SJames Wright PetscCall(PetscOptionsInt("-num_coords", "Number of coordinates", "", PETSC_FALSE, &num_coords, NULL));
3098480730SJames Wright PetscCall(PetscOptionsInt("-dim", "Dimension of the coordinates", "", PETSC_FALSE, &dim, NULL));
3198480730SJames Wright PetscCall(PetscOptionsInt("-bucket_size", "Size of the bucket at each leaf", "", bucket_size, &bucket_size, NULL));
3298480730SJames Wright PetscCall(PetscOptionsBoundedInt("-loops", "Number of times to loop through KDTree creation and querying", "", loops, &loops, NULL, 1));
3398480730SJames Wright PetscCall(PetscOptionsEnum("-copy_mode", "Copy mode to use with KDTree", "PetscKDTreeCreate", PetscCopyModes, (PetscEnum)copy_mode, (PetscEnum *)©_mode, NULL));
3498480730SJames Wright PetscCall(PetscOptionsBool("-view_tree", "View the KDTree", "", view_tree, &view_tree, NULL));
3598480730SJames Wright PetscCall(PetscOptionsBool("-view_performance", "View the performance speed of the KDTree", "", view_performance, &view_performance, NULL));
3698480730SJames Wright PetscCall(PetscOptionsBool("-test_tree_points", "Test querying tree points against itself", "", test_tree_points, &test_tree_points, NULL));
3798480730SJames Wright PetscCall(PetscOptionsBool("-test_rand_points", "Test querying random points via brute force", "", test_rand_points, &test_rand_points, NULL));
3898480730SJames Wright PetscCall(PetscOptionsBool("-query_rand_points", "Query querying random points", "", query_rand_points, &query_rand_points, NULL));
3998480730SJames 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));
4098480730SJames Wright PetscOptionsEnd();
4198480730SJames Wright
4298480730SJames Wright coords_size = num_coords * dim;
4398480730SJames Wright PetscCall(PetscMalloc1(coords_size, &coords));
4498480730SJames Wright PetscCall(PetscRandomCreate(comm, &random));
4598480730SJames Wright
4698480730SJames Wright for (PetscInt loop_count = 0; loop_count < loops; loop_count++) {
4798480730SJames Wright PetscCall(PetscRandomGetValuesReal(random, coords_size, coords));
4898480730SJames Wright
4998480730SJames Wright PetscCall(PetscKDTreeCreate(num_coords, dim, coords, copy_mode, bucket_size, &tree));
5098480730SJames Wright num_trees_built++;
5198480730SJames Wright if (view_tree) PetscCall(PetscKDTreeView(tree, NULL));
5298480730SJames Wright
5398480730SJames Wright if (test_tree_points) { // Test querying the current coordinates
5498480730SJames Wright PetscCount *indices;
5598480730SJames Wright PetscReal *distances;
5698480730SJames Wright num_points_queried += num_coords;
5798480730SJames Wright
5898480730SJames Wright PetscCall(PetscMalloc2(num_coords, &indices, num_coords, &distances));
5998480730SJames Wright PetscCall(PetscKDTreeQueryPointsNearestNeighbor(tree, num_coords, coords, PETSC_MACHINE_EPSILON * 1e3, indices, distances));
6098480730SJames Wright for (PetscInt i = 0; i < num_coords; i++) {
6198480730SJames 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]));
6298480730SJames Wright }
6398480730SJames Wright PetscCall(PetscFree2(indices, distances));
6498480730SJames Wright }
6598480730SJames Wright
6698480730SJames Wright if (num_rand_points > 0) {
6798480730SJames Wright PetscCount *indices;
6898480730SJames Wright PetscReal *distances;
6998480730SJames Wright PetscReal *rand_points;
7098480730SJames Wright PetscInt rand_queries_size = num_rand_points * dim;
7198480730SJames Wright
7298480730SJames Wright num_points_queried += num_rand_points;
7398480730SJames Wright PetscCall(PetscMalloc3(rand_queries_size, &rand_points, num_rand_points, &indices, num_rand_points, &distances));
7498480730SJames Wright PetscCall(PetscRandomGetValuesReal(random, rand_queries_size, rand_points));
7598480730SJames Wright PetscCall(PetscKDTreeQueryPointsNearestNeighbor(tree, num_rand_points, rand_points, 0., indices, distances));
7698480730SJames Wright
7798480730SJames Wright if (test_rand_points) {
7898480730SJames Wright for (PetscInt i = 0; i < num_rand_points; i++) {
7998480730SJames Wright PetscReal *rand_point = &rand_points[dim * i], nearest_distance = PETSC_MAX_REAL;
8098480730SJames Wright PetscInt index = -1;
8198480730SJames Wright for (PetscInt j = 0; j < num_coords; j++) {
8298480730SJames Wright PetscReal dist = Distance(dim, rand_point, &coords[dim * j]);
8398480730SJames Wright if (dist < nearest_distance) {
8498480730SJames Wright nearest_distance = dist;
8598480730SJames Wright index = j;
8698480730SJames Wright }
8798480730SJames Wright }
8898480730SJames Wright if (indices[i] != index)
8998480730SJames 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));
9098480730SJames Wright }
9198480730SJames Wright }
9298480730SJames Wright PetscCall(PetscFree3(rand_points, indices, distances));
9398480730SJames Wright }
9498480730SJames Wright }
9598480730SJames Wright
9698480730SJames Wright if (view_performance) {
9798480730SJames Wright PetscLogEvent kdtree_build_log, kdtree_query_log;
9898480730SJames Wright PetscEventPerfInfo build_perf_info, query_perf_info;
9998480730SJames Wright
10098480730SJames Wright PetscCall(PetscLogEventGetId("KDTreeBuild", &kdtree_build_log));
10198480730SJames Wright PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, kdtree_build_log, &build_perf_info));
10298480730SJames Wright PetscCall(PetscLogEventGetId("KDTreeQuery", &kdtree_query_log));
10398480730SJames Wright PetscCall(PetscLogEventGetPerfInfo(PETSC_DETERMINE, kdtree_query_log, &query_perf_info));
10498480730SJames Wright PetscCall(PetscPrintf(comm, "KDTreeBuild %.6e for %" PetscInt_FMT " trees built\n", build_perf_info.time, num_trees_built));
10598480730SJames Wright PetscCall(PetscPrintf(comm, "\tTime per tree: %.6e\n", build_perf_info.time / num_trees_built));
10698480730SJames Wright PetscCall(PetscPrintf(comm, "KDTreeQuery %.6e for %" PetscCount_FMT " queries\n", query_perf_info.time, (PetscCount)num_points_queried));
10798480730SJames Wright PetscCall(PetscPrintf(comm, "\tTime per query: %.6e\n", query_perf_info.time / num_points_queried));
10898480730SJames Wright }
10998480730SJames Wright
11098480730SJames Wright if (copy_mode != PETSC_OWN_POINTER) PetscCall(PetscFree(coords));
11198480730SJames Wright PetscCall(PetscKDTreeDestroy(&tree));
11298480730SJames Wright PetscCall(PetscRandomDestroy(&random));
11398480730SJames Wright PetscCall(PetscFinalize());
11498480730SJames Wright return 0;
11598480730SJames Wright }
11698480730SJames Wright
11798480730SJames Wright /*TEST
11898480730SJames Wright testset:
11998480730SJames Wright suffix: kdtree
12098480730SJames 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
12198480730SJames Wright test:
12298480730SJames Wright suffix: 1D
12398480730SJames Wright args: -dim 1
12498480730SJames Wright test:
12598480730SJames Wright suffix: 2D
12698480730SJames Wright args: -dim 2
12798480730SJames Wright test:
12898480730SJames Wright suffix: 3D
12998480730SJames Wright args: -dim 3
13098480730SJames Wright test:
13198480730SJames Wright suffix: 3D_small
13298480730SJames Wright args: -dim 3 -num_coords 13
13398480730SJames Wright
13498480730SJames Wright testset:
13598480730SJames Wright suffix: kdtree_3D_large
13698480730SJames Wright args: -dim 3 -num_coords 350 -test_tree_points -test_rand_points -num_rand_points 300 -view_performance false -kdtree_debug
137*3886731fSPierre Jolivet output_file: output/empty.out
13898480730SJames Wright test:
13998480730SJames Wright test:
14098480730SJames Wright suffix: copy
14198480730SJames Wright args: -copy_mode copy_values
14298480730SJames Wright test:
14398480730SJames Wright suffix: use
14498480730SJames Wright args: -copy_mode use_pointer
14598480730SJames Wright TEST*/
146