5 #include <Compadre_Config.h> 
   10 #ifdef COMPADRE_USE_MPI 
   14 template<
typename T1, 
typename T2>
 
   15 std::pair<T2,T1> 
swap_pair(
const std::pair<T1,T2> &item)
 
   17     return std::pair<T2,T1>(item.second, item.first);
 
   20 template<
typename T1, 
typename T2>
 
   21 std::multimap<T2,T1> 
invert_map(
const std::map<T1,T2> &original_map)
 
   23     std::multimap<T2,T1> new_map;
 
   24     std::transform(original_map.begin(), original_map.end(), std::inserter(new_map, new_map.begin()), swap_pair<T1,T2>);
 
   28 using namespace Compadre;
 
   30 int main (
int argc, 
char* args[]) {
 
   33 #ifdef COMPADRE_USE_MPI 
   34 MPI_Init(&argc, &args);
 
   41 bool all_passed = 
true;
 
   50         int arg5toi = atoi(args[4]);
 
   52             search_type = arg5toi;
 
   55     bool do_radius_search = (search_type==0);
 
   56     bool do_knn_search = (search_type==1);
 
   57     printf(
"do_radius_search: %d\n", do_radius_search);
 
   58     printf(
"do_knn_search: %d\n", do_knn_search);
 
   62     double h_multiplier = 3.0; 
 
   64         double arg4tof = atof(args[3]);
 
   66             h_multiplier = arg4tof;
 
   72     int number_target_coords = 200; 
 
   74         int arg3toi = atoi(args[2]);
 
   76             number_target_coords = arg3toi;
 
   84         int arg2toi = atoi(args[1]);
 
   94     double h_spacing = 1./std::pow(number_target_coords, 0.50 + 0.2*(dimension==2));
 
   95     int n_neg1_to_1 = 2*(1/h_spacing) + 1; 
 
   98     int number_source_coords = (n_neg1_to_1+1)*(n_neg1_to_1+1);
 
   99     if (dimension==3) number_source_coords *= (n_neg1_to_1+1);
 
  100     printf(
"target coords: %d\n", number_target_coords);
 
  101     printf(
"source coords: %d\n", number_source_coords);
 
  104     Kokkos::View<double**, Kokkos::DefaultHostExecutionSpace> source_coords(
"source coordinates", 
 
  105             number_source_coords, 3);
 
  108     Kokkos::View<double**, Kokkos::DefaultHostExecutionSpace> target_coords(
"target coordinates", number_target_coords, 3);
 
  111     int source_index = 0;
 
  112     double this_coord[3] = {0,0,0};
 
  113     for (
int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
 
  114         this_coord[0] = i*h_spacing;
 
  115         for (
int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
 
  116             this_coord[1] = j*h_spacing;
 
  117             for (
int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
 
  118                 this_coord[2] = k*h_spacing;
 
  120                     source_coords(source_index,0) = this_coord[0]; 
 
  121                     source_coords(source_index,1) = this_coord[1]; 
 
  122                     source_coords(source_index,2) = this_coord[2]; 
 
  127                 source_coords(source_index,0) = this_coord[0]; 
 
  128                 source_coords(source_index,1) = this_coord[1]; 
 
  129                 source_coords(source_index,2) = 0;
 
  134             source_coords(source_index,0) = this_coord[0]; 
 
  135             source_coords(source_index,1) = 0;
 
  136             source_coords(source_index,2) = 0;
 
  142     for(
int i=0; i<number_target_coords; i++){
 
  145         double rand_dir[3] = {0,0,0};
 
  147         for (
int j=0; j<dimension; ++j) {
 
  149             rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
 
  153         for (
int j=0; j<dimension; ++j) {
 
  154             target_coords(i,j) = rand_dir[j];
 
  159     if (do_radius_search) {
 
  162             printf(
"\n1D compressed row neighbor lists:\n");
 
  167             Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> neighbor_lists(
"neighbor lists", 
 
  168                     number_target_coords); 
 
  169             Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> number_neighbors_list(
"number of neighbors list", 
 
  170                     number_target_coords); 
 
  173             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon(
"h supports", number_target_coords);
 
  174             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> random_vec(
"random vector", number_target_coords);
 
  176             Kokkos::parallel_for(
"random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  177                     (0,number_target_coords), [&](
const int i) {
 
  179                 random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
 
  180                 epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
 
  187             size_t total_num_neighbors = point_cloud_search.generateCRNeighborListsFromRadiusSearch(
true ,
 
  188                 target_coords, neighbor_lists, number_neighbors_list, epsilon);
 
  189             printf(
"total num neighbors: %lu\n", total_num_neighbors);
 
  192             Kokkos::resize(neighbor_lists, total_num_neighbors);
 
  195             point_cloud_search.generateCRNeighborListsFromRadiusSearch(
false ,
 
  196                 target_coords, neighbor_lists, number_neighbors_list, epsilon);
 
  200             double radius_search_time = timer.seconds();
 
  201             printf(
"nanoflann search time: %f s\n", radius_search_time);
 
  205             std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
 
  206             Kokkos::parallel_for(
"point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  207                     (0,number_target_coords), [&](
const int i) {
 
  208                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  209                 for (
int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
 
  210                     point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = 1.0;
 
  213             double point_cloud_convert_time = timer.seconds();
 
  214             printf(
"point cloud convert time: %f s\n", point_cloud_convert_time);
 
  219             std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
 
  220             Kokkos::parallel_for(
"brute force radius search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  221                     (0,number_target_coords), [&](
const int i) {
 
  222                 auto &brute_force_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(brute_force_neighbor_list[i]);
 
  223                 for (
int j=0; j<number_source_coords; ++j) {
 
  225                     for (
int k=0; k<dimension; ++k) {
 
  226                         dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
 
  228                     dist = std::sqrt(dist);
 
  230                     if (dist<(h_multiplier + random_vec(i))*h_spacing) {
 
  231                         brute_force_neighbor_list_i[j]=dist;
 
  235             double brute_force_search_time = timer.seconds();
 
  236             printf(
"brute force search time: %f s\n", brute_force_search_time);
 
  242             Kokkos::parallel_reduce(
"check brute force search in point cloud search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  243                     (0,number_target_coords), [&](
const int i, 
int& t_num_passed) {
 
  244                 bool all_found = 
true;
 
  245                 for (
auto it=brute_force_neighbor_list[i].begin(); it!=brute_force_neighbor_list[i].end(); ++it) {
 
  246                     if (point_cloud_neighbor_list[i].count(it->first)!=1) {
 
  250                 if (all_found) t_num_passed++;
 
  251             }, Kokkos::Sum<int>(num_passed));
 
  253             if (num_passed != number_target_coords) {
 
  254                 printf(
"Brute force neighbor not found in point cloud list\n");
 
  260             Kokkos::parallel_reduce(
"original in brute force search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  261                     (0,number_target_coords), KOKKOS_LAMBDA(
const int i, 
int& t_num_passed) {
 
  262                 bool all_found = 
true;
 
  263                 for (
int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
 
  264                     if (brute_force_neighbor_list[i].count(nla.getNeighborHost(i,j))!=1) {
 
  268                 if (all_found) t_num_passed++;
 
  269             }, Kokkos::Sum<int>(num_passed));
 
  271             if (num_passed != number_target_coords) {
 
  272                 printf(
"Point cloud neighbor not found in brute force list\n");
 
  276             double check_time = timer.seconds();
 
  277             printf(
"check time: %f s\n", check_time);
 
  281             printf(
"\n2D neighbor lists:\n");
 
  286             Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists(
"neighbor lists", 
 
  287                     number_target_coords, 2); 
 
  290             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon(
"h supports", number_target_coords);
 
  291             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> random_vec(
"random vector", number_target_coords);
 
  293             Kokkos::parallel_for(
"random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  294                     (0,number_target_coords), [&](
const int i) {
 
  296                 random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
 
  297                 epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
 
  304             auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(
true ,
 
  305                 target_coords, neighbor_lists, epsilon);
 
  306             printf(
"max num neighbors: %lu\n", max_num_neighbors);
 
  309             neighbor_lists = Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace>(
"neighbor lists", 
 
  310                 number_target_coords, 1+max_num_neighbors); 
 
  313             max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(
false ,
 
  314                 target_coords, neighbor_lists, epsilon);
 
  316             double radius_search_time = timer.seconds();
 
  317             printf(
"nanoflann search time: %f s\n", radius_search_time);
 
  321             std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
 
  322             Kokkos::parallel_for(
"point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  323                     (0,number_target_coords), [&](
const int i) {
 
  324                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  325                 for (
int j=1; j<=neighbor_lists(i,0); ++j) {
 
  326                     point_cloud_neighbor_list_i[neighbor_lists(i,j)] = 1.0;
 
  329             double point_cloud_convert_time = timer.seconds();
 
  330             printf(
"point cloud convert time: %f s\n", point_cloud_convert_time);
 
  335             std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
 
  336             Kokkos::parallel_for(
"brute force radius search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  337                     (0,number_target_coords), [&](
const int i) {
 
  338                 auto &brute_force_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(brute_force_neighbor_list[i]);
 
  339                 for (
int j=0; j<number_source_coords; ++j) {
 
  341                     for (
int k=0; k<dimension; ++k) {
 
  342                         dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
 
  344                     dist = std::sqrt(dist);
 
  346                     if (dist<(h_multiplier + random_vec(i))*h_spacing) {
 
  347                         brute_force_neighbor_list_i[j]=dist;
 
  351             double brute_force_search_time = timer.seconds();
 
  352             printf(
"brute force search time: %f s\n", brute_force_search_time);
 
  358             Kokkos::parallel_reduce(
"check brute force search in point cloud search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  359                     (0,number_target_coords), [&](
const int i, 
int& t_num_passed) {
 
  360                 bool all_found = 
true;
 
  361                 for (
auto it=brute_force_neighbor_list[i].begin(); it!=brute_force_neighbor_list[i].end(); ++it) {
 
  362                     if (point_cloud_neighbor_list[i].count(it->first)!=1) {
 
  366                 if (all_found) t_num_passed++;
 
  367             }, Kokkos::Sum<int>(num_passed));
 
  369             if (num_passed != number_target_coords) {
 
  370                 printf(
"Brute force neighbor not found in point cloud list\n");
 
  376             Kokkos::parallel_reduce(
"original in brute force search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  377                     (0,number_target_coords), KOKKOS_LAMBDA(
const int i, 
int& t_num_passed) {
 
  378                 bool all_found = 
true;
 
  379                 for (
int j=1; j<=neighbor_lists(i,0); ++j) {
 
  380                     if (brute_force_neighbor_list[i].count(neighbor_lists(i,j))!=1) {
 
  384                 if (all_found) t_num_passed++;
 
  385             }, Kokkos::Sum<int>(num_passed));
 
  387             if (num_passed != number_target_coords) {
 
  388                 printf(
"Point cloud neighbor not found in brute force list\n");
 
  392             double check_time = timer.seconds();
 
  393             printf(
"check time: %f s\n", check_time);
 
  397             printf(
"\n2D neighbor lists converted to 1D compressed row neighbor lists:\n");
 
  402             Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists(
"neighbor lists", 
 
  403                     number_target_coords, 2); 
 
  406             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon(
"h supports", number_target_coords);
 
  407             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> random_vec(
"random vector", number_target_coords);
 
  409             Kokkos::parallel_for(
"random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  410                     (0,number_target_coords), [&](
const int i) {
 
  412                 random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
 
  413                 epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
 
  420             auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(
true ,
 
  421                 target_coords, neighbor_lists, epsilon);
 
  422             printf(
"max num neighbors: %lu\n", max_num_neighbors);
 
  425             neighbor_lists = Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace>(
"neighbor lists", 
 
  426                 number_target_coords, 1+max_num_neighbors); 
 
  429             max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(
false ,
 
  430                 target_coords, neighbor_lists, epsilon);
 
  434             double radius_search_time = timer.seconds();
 
  435             printf(
"nanoflann search time: %f s\n", radius_search_time);
 
  439             std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
 
  440             Kokkos::parallel_for(
"point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  441                     (0,number_target_coords), [&](
const int i) {
 
  442                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  443                 for (
int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
 
  444                     point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = 1.0;
 
  447             double point_cloud_convert_time = timer.seconds();
 
  448             printf(
"point cloud convert time: %f s\n", point_cloud_convert_time);
 
  453             std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
 
  454             Kokkos::parallel_for(
"brute force radius search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  455                     (0,number_target_coords), [&](
const int i) {
 
  456                 auto &brute_force_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(brute_force_neighbor_list[i]);
 
  457                 for (
int j=0; j<number_source_coords; ++j) {
 
  459                     for (
int k=0; k<dimension; ++k) {
 
  460                         dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
 
  462                     dist = std::sqrt(dist);
 
  464                     if (dist<(h_multiplier + random_vec(i))*h_spacing) {
 
  465                         brute_force_neighbor_list_i[j]=dist;
 
  469             double brute_force_search_time = timer.seconds();
 
  470             printf(
"brute force search time: %f s\n", brute_force_search_time);
 
  476             Kokkos::parallel_reduce(
"check brute force search in point cloud search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  477                     (0,number_target_coords), [&](
const int i, 
int& t_num_passed) {
 
  478                 bool all_found = 
true;
 
  479                 for (
auto it=brute_force_neighbor_list[i].begin(); it!=brute_force_neighbor_list[i].end(); ++it) {
 
  480                     if (point_cloud_neighbor_list[i].count(it->first)!=1) {
 
  484                 if (all_found) t_num_passed++;
 
  485             }, Kokkos::Sum<int>(num_passed));
 
  487             if (num_passed != number_target_coords) {
 
  488                 printf(
"Brute force neighbor not found in point cloud list\n");
 
  494             Kokkos::parallel_reduce(
"original in brute force search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  495                     (0,number_target_coords), KOKKOS_LAMBDA(
const int i, 
int& t_num_passed) {
 
  496                 bool all_found = 
true;
 
  497                 for (
int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
 
  498                     if (brute_force_neighbor_list[i].count(nla.getNeighborHost(i,j))!=1) {
 
  502                 if (all_found) t_num_passed++;
 
  503             }, Kokkos::Sum<int>(num_passed));
 
  505             if (num_passed != number_target_coords) {
 
  506                 printf(
"Point cloud neighbor not found in brute force list\n");
 
  510             double check_time = timer.seconds();
 
  511             printf(
"check time: %f s\n", check_time);
 
  514     else if (do_knn_search) {
 
  519             printf(
"\n1D compressed row neighbor lists:\n");
 
  525             Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> neighbor_lists(
"neighbor lists", 
 
  526                     number_target_coords); 
 
  527             Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> number_neighbors_list(
"number of neighbors list", 
 
  528                     number_target_coords); 
 
  531             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon(
"h supports", number_target_coords);
 
  536             auto total_num_neighbors = point_cloud_search.generateCRNeighborListsFromKNNSearch(
true ,
 
  537                     target_coords, neighbor_lists, number_neighbors_list, epsilon, min_neighbors, 1.5 );
 
  538             printf(
"total num neighbors: %lu\n", total_num_neighbors);
 
  541             Kokkos::resize(neighbor_lists, total_num_neighbors);
 
  544             point_cloud_search.generateCRNeighborListsFromKNNSearch(
false , 
 
  545                     target_coords, neighbor_lists, number_neighbors_list, epsilon, min_neighbors, 1.5 );
 
  551             std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
 
  552             Kokkos::parallel_for(
"point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  553                     (0,number_target_coords), [&](
const int i) {
 
  554                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  555                 for (
int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
 
  557                     for (
int k=0; k<dimension; ++k) {
 
  558                         dist += (target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k))*(target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k));
 
  560                     dist = std::sqrt(dist);
 
  561                     point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = dist;
 
  564             double point_cloud_convert_time = timer.seconds();
 
  565             printf(
"point cloud convert time: %f s\n", point_cloud_convert_time);
 
  570             int epsilon_diff_from_knn_dist_multiplied = 0;
 
  571             Kokkos::parallel_reduce(
"check kth neighbor in point cloud", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  572                     (0,number_target_coords), [&](
const int i, 
int& t_epsilon_diff_from_knn_dist_multiplied) {
 
  573                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  574                 std::multimap<double, int> inverted_map = 
invert_map(point_cloud_neighbor_list_i);
 
  576                 for (
auto it=inverted_map.begin(); it!=inverted_map.end(); ++it) {
 
  577                     if (j==min_neighbors-1) {
 
  578                         if (std::abs(epsilon(i)/1.5 - it->first) > eps) {
 
  579                             t_epsilon_diff_from_knn_dist_multiplied++;
 
  585             }, Kokkos::Sum<int>(epsilon_diff_from_knn_dist_multiplied));
 
  586             if (epsilon_diff_from_knn_dist_multiplied > 0) {
 
  587                 printf(
"1.5*kth_neighbor_distance differs from calculation in epsilon vector.\n");
 
  596             int total_insufficient_neighbors = 0;
 
  597             std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
 
  598             Kokkos::parallel_reduce(
"brute force knn search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  599                     (0,number_target_coords), [&](
const int i, 
int& t_total_insufficient_neighbors) {
 
  600                 auto &brute_force_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(brute_force_neighbor_list[i]);
 
  601                 for (
int j=0; j<number_source_coords; ++j) {
 
  603                     for (
int k=0; k<dimension; ++k) {
 
  604                         dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
 
  606                     dist = std::sqrt(dist);
 
  609                     if (dist<(epsilon(i)/1.5 + eps)) {
 
  610                         brute_force_neighbor_list_i[j]=dist;
 
  614                 if (brute_force_neighbor_list_i.size() < (size_t)min_neighbors) t_total_insufficient_neighbors++;
 
  615             }, Kokkos::Sum<int>(total_insufficient_neighbors));
 
  616             if (total_insufficient_neighbors > 0) {
 
  617                 printf(
"less neighbors found using kth_neighbor_distance+eps than using knn search.\n");
 
  620             double brute_force_search_time = timer.seconds();
 
  621             printf(
"brute force search time: %f s\n", brute_force_search_time);
 
  626             printf(
"\n2D neighbor lists:\n");
 
  632             Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists(
"neighbor lists", 
 
  633                     number_target_coords, 1+min_neighbors); 
 
  636             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon(
"h supports", number_target_coords);
 
  640             auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(
true ,
 
  641                     target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 );
 
  642             printf(
"max num neighbors: %lu\n", max_num_neighbors);
 
  645             Kokkos::resize(neighbor_lists, neighbor_lists.extent(0), 1+max_num_neighbors);
 
  648             max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , 
 
  649                     target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 );
 
  653             std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
 
  654             Kokkos::parallel_for(
"point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  655                     (0,number_target_coords), [&](
const int i) {
 
  656                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  657                 for (
int j=1; j<=neighbor_lists(i,0); ++j) {
 
  659                     for (
int k=0; k<dimension; ++k) {
 
  660                         dist += (target_coords(i,k)-source_coords(neighbor_lists(i,j),k))*(target_coords(i,k)-source_coords(neighbor_lists(i,j),k));
 
  662                     dist = std::sqrt(dist);
 
  663                     point_cloud_neighbor_list_i[neighbor_lists(i,j)] = dist;
 
  666             double point_cloud_convert_time = timer.seconds();
 
  667             printf(
"point cloud convert time: %f s\n", point_cloud_convert_time);
 
  672             int epsilon_diff_from_knn_dist_multiplied = 0;
 
  673             Kokkos::parallel_reduce(
"check kth neighbor in point cloud", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  674                     (0,number_target_coords), [&](
const int i, 
int& t_epsilon_diff_from_knn_dist_multiplied) {
 
  675                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  676                 std::multimap<double, int> inverted_map = 
invert_map(point_cloud_neighbor_list_i);
 
  678                 for (
auto it=inverted_map.begin(); it!=inverted_map.end(); ++it) {
 
  679                     if (j==min_neighbors-1) {
 
  680                         if (std::abs(epsilon(i)/1.5 - it->first) > eps) {
 
  681                             t_epsilon_diff_from_knn_dist_multiplied++;
 
  687             }, Kokkos::Sum<int>(epsilon_diff_from_knn_dist_multiplied));
 
  688             if (epsilon_diff_from_knn_dist_multiplied > 0) {
 
  689                 printf(
"1.5*kth_neighbor_distance differs from calculation in epsilon vector.\n");
 
  698             int total_insufficient_neighbors = 0;
 
  699             std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
 
  700             Kokkos::parallel_reduce(
"brute force knn search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  701                     (0,number_target_coords), [&](
const int i, 
int& t_total_insufficient_neighbors) {
 
  702                 auto &brute_force_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(brute_force_neighbor_list[i]);
 
  703                 for (
int j=0; j<number_source_coords; ++j) {
 
  705                     for (
int k=0; k<dimension; ++k) {
 
  706                         dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
 
  708                     dist = std::sqrt(dist);
 
  711                     if (dist<(epsilon(i)/1.5 + eps)) {
 
  712                         brute_force_neighbor_list_i[j]=dist;
 
  716                 if (brute_force_neighbor_list_i.size() < (size_t)min_neighbors) t_total_insufficient_neighbors++;
 
  717             }, Kokkos::Sum<int>(total_insufficient_neighbors));
 
  718             if (total_insufficient_neighbors > 0) {
 
  719                 printf(
"less neighbors found using kth_neighbor_distance+eps than using knn search.\n");
 
  722             double brute_force_search_time = timer.seconds();
 
  723             printf(
"brute force search time: %f s\n", brute_force_search_time);
 
  728             printf(
"\n2D neighbor lists converted to 1D compressed row neighbor lists:\n");
 
  734             Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists(
"neighbor lists", 
 
  735                     number_target_coords, 1+min_neighbors); 
 
  738             Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon(
"h supports", number_target_coords);
 
  742             auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(
true ,
 
  743                     target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 );
 
  744             printf(
"max num neighbors: %lu\n", max_num_neighbors);
 
  747             Kokkos::resize(neighbor_lists, neighbor_lists.extent(0), 1+max_num_neighbors);
 
  750             max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , 
 
  751                     target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 );
 
  757             std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
 
  758             Kokkos::parallel_for(
"point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  759                     (0,number_target_coords), [&](
const int i) {
 
  760                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  761                 for (
int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
 
  763                     for (
int k=0; k<dimension; ++k) {
 
  764                         dist += (target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k))*(target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k));
 
  766                     dist = std::sqrt(dist);
 
  767                     point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = dist;
 
  770             double point_cloud_convert_time = timer.seconds();
 
  771             printf(
"point cloud convert time: %f s\n", point_cloud_convert_time);
 
  776             int epsilon_diff_from_knn_dist_multiplied = 0;
 
  777             Kokkos::parallel_reduce(
"check kth neighbor in point cloud", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  778                     (0,number_target_coords), [&](
const int i, 
int& t_epsilon_diff_from_knn_dist_multiplied) {
 
  779                 auto &point_cloud_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(point_cloud_neighbor_list[i]);
 
  780                 std::multimap<double, int> inverted_map = 
invert_map(point_cloud_neighbor_list_i);
 
  782                 for (
auto it=inverted_map.begin(); it!=inverted_map.end(); ++it) {
 
  783                     if (j==min_neighbors-1) {
 
  784                         if (std::abs(epsilon(i)/1.5 - it->first) > eps) {
 
  785                             t_epsilon_diff_from_knn_dist_multiplied++;
 
  791             }, Kokkos::Sum<int>(epsilon_diff_from_knn_dist_multiplied));
 
  792             if (epsilon_diff_from_knn_dist_multiplied > 0) {
 
  793                 printf(
"1.5*kth_neighbor_distance differs from calculation in epsilon vector.\n");
 
  802             int total_insufficient_neighbors = 0;
 
  803             std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
 
  804             Kokkos::parallel_reduce(
"brute force knn search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
 
  805                     (0,number_target_coords), [&](
const int i, 
int& t_total_insufficient_neighbors) {
 
  806                 auto &brute_force_neighbor_list_i = 
const_cast<std::map<int, double>& 
>(brute_force_neighbor_list[i]);
 
  807                 for (
int j=0; j<number_source_coords; ++j) {
 
  809                     for (
int k=0; k<dimension; ++k) {
 
  810                         dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
 
  812                     dist = std::sqrt(dist);
 
  815                     if (dist<(epsilon(i)/1.5 + eps)) {
 
  816                         brute_force_neighbor_list_i[j]=dist;
 
  820                 if (brute_force_neighbor_list_i.size() < (size_t)min_neighbors) t_total_insufficient_neighbors++;
 
  821             }, Kokkos::Sum<int>(total_insufficient_neighbors));
 
  822             if (total_insufficient_neighbors > 0) {
 
  823                 printf(
"less neighbors found using kth_neighbor_distance+eps than using knn search.\n");
 
  826             double brute_force_search_time = timer.seconds();
 
  827             printf(
"brute force search time: %f s\n", brute_force_search_time);
 
  832         printf(
"do_radius_search and do_knn_search both false. Invalid combination.\n");
 
  838 #ifdef COMPADRE_USE_MPI 
  844     fprintf(stdout, 
"Passed test \n");
 
  847     fprintf(stdout, 
"Failed test \n");
 
Class handling Kokkos command line arguments and returning parameters. 
 
PointCloudSearch< view_type > CreatePointCloudSearch(view_type src_view, const local_index_type dimensions=-1, const local_index_type max_leaf=-1)
CreatePointCloudSearch allows for the construction of an object of type PointCloudSearch with templat...
 
std::pair< T2, T1 > swap_pair(const std::pair< T1, T2 > &item)
 
NeighborLists< view_type > CreateNeighborLists(view_type neighbor_lists, view_type number_of_neighbors_list)
CreateNeighborLists allows for the construction of an object of type NeighborLists with template dedu...
 
std::multimap< T2, T1 > invert_map(const std::map< T1, T2 > &original_map)
 
int main(int argc, char *args[])
[Parse Command Line Arguments] 
 
static KOKKOS_INLINE_FUNCTION int getNP(const int m, const int dimension=3, const ReconstructionSpace r_space=ReconstructionSpace::ScalarTaylorPolynomial)
Returns size of the basis for a given polynomial order and dimension General to dimension 1...
 
NeighborLists< view_type_1d > Convert2DToCompressedRowNeighborLists(view_type_2d neighbor_lists)
Converts 2D neighbor lists to compressed row neighbor lists.