Compadre  1.2.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
NeighborSearchTest.cpp
Go to the documentation of this file.
1 #include <vector>
2 #include <map>
3 #include <random>
4 
5 #include <Compadre_Config.h>
6 #include <Compadre_GMLS.hpp>
9 
10 #ifdef COMPADRE_USE_MPI
11 #include <mpi.h>
12 #endif
13 
14 template<typename T1, typename T2>
15 std::pair<T2,T1> swap_pair(const std::pair<T1,T2> &item)
16 {
17  return std::pair<T2,T1>(item.second, item.first);
18 }
19 
20 template<typename T1, typename T2>
21 std::multimap<T2,T1> invert_map(const std::map<T1,T2> &original_map)
22 {
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>);
25  return new_map;
26 }
27 
28 using namespace Compadre;
29 
30 int main (int argc, char* args[]) {
31 
32 // initializes MPI (if available) with command line arguments given
33 #ifdef COMPADRE_USE_MPI
34 MPI_Init(&argc, &args);
35 #endif
36 
37 // initializes Kokkos with command line arguments given
38 auto kp = KokkosParser(argc, args, true);
39 
40 // becomes false if the computed solution not within the failure_threshold of the actual solution
41 bool all_passed = true;
42 
43 {
44  // seed random generator
45  srand(1234);
46  Kokkos::Timer timer;
47 
48  int search_type = 0; // 0 - radius, 1 - knn
49  if (argc >= 5) {
50  int arg5toi = atoi(args[4]);
51  if (arg5toi >= 0) {
52  search_type = arg5toi;
53  }
54  }
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);
59 
60  // check if 4 arguments are given from the command line
61  // multiplier times h spacing for search
62  double h_multiplier = 3.0; // dimension 3 by default
63  if (argc >= 4) {
64  double arg4tof = atof(args[3]);
65  if (arg4tof > 0) {
66  h_multiplier = arg4tof;
67  }
68  }
69 
70  // check if 3 arguments are given from the command line
71  // set the number of target sites where we will reconstruct the target functionals at
72  int number_target_coords = 200; // 200 target sites by default
73  if (argc >= 3) {
74  int arg3toi = atoi(args[2]);
75  if (arg3toi > 0) {
76  number_target_coords = arg3toi;
77  }
78  }
79 
80  // check if 2 arguments are given from the command line
81  // set the number of dimensions for points and the search
82  int dimension = 3; // 3D by default
83  if (argc >= 2) {
84  int arg2toi = atoi(args[1]);
85  if (arg2toi > 0) {
86  dimension = arg2toi;
87  }
88  }
89 
90  //// minimum neighbors for unisolvency is the same as the size of the polynomial basis
91  //const int min_neighbors = Compadre::GMLS::getNP(order, dimension);
92 
93  // approximate spacing of source sites
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; // always odd
96 
97  // number of source coordinate sites that will fill a box of [-1,1]x[-1,1]x[-1,1] with a spacing approximately h
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);
102 
103  // coordinates of source sites
104  Kokkos::View<double**, Kokkos::DefaultHostExecutionSpace> source_coords("source coordinates",
105  number_source_coords, 3);
106 
107  // coordinates of target sites
108  Kokkos::View<double**, Kokkos::DefaultHostExecutionSpace> target_coords("target coordinates", number_target_coords, 3);
109 
110  // fill source coordinates with a uniform grid
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;
119  if (dimension==3) {
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];
123  source_index++;
124  }
125  }
126  if (dimension==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;
130  source_index++;
131  }
132  }
133  if (dimension==1) {
134  source_coords(source_index,0) = this_coord[0];
135  source_coords(source_index,1) = 0;
136  source_coords(source_index,2) = 0;
137  source_index++;
138  }
139  }
140 
141  // fill target coords somewhere inside of [-0.5,0.5]x[-0.5,0.5]x[-0.5,0.5]
142  for(int i=0; i<number_target_coords; i++){
143 
144  // first, we get a uniformly random distributed direction
145  double rand_dir[3] = {0,0,0};
146 
147  for (int j=0; j<dimension; ++j) {
148  // rand_dir[j] is in [-0.5, 0.5]
149  rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
150  }
151 
152  // then we get a uniformly random radius
153  for (int j=0; j<dimension; ++j) {
154  target_coords(i,j) = rand_dir[j];
155  }
156 
157  }
158 
159  if (do_radius_search) {
160  timer.reset();
161  { // 1D compressed row neighbor lists
162  printf("\n1D compressed row neighbor lists:\n");
163  // Point cloud construction for neighbor search
164  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
165  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
166 
167  Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
168  number_target_coords); // first column is # of neighbors
169  Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> number_neighbors_list("number of neighbors list",
170  number_target_coords); // first column is # of neighbors
171 
172  // each target site has a window size
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);
175 
176  Kokkos::parallel_for("random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
177  (0,number_target_coords), [&](const int i) {
178  // value in [-.25, .25]
179  random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
180  epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
181  });
182  Kokkos::fence();
183 
184  // query the point cloud to generate the neighbor lists using a radius search
185 
186  // start with a dry-run, but without enough room to store the results
187  size_t total_num_neighbors = point_cloud_search.generateCRNeighborListsFromRadiusSearch(true /* dry run */,
188  target_coords, neighbor_lists, number_neighbors_list, epsilon);
189  printf("total num neighbors: %lu\n", total_num_neighbors);
190 
191  // resize neighbor lists to be large enough to hold the results
192  Kokkos::resize(neighbor_lists, total_num_neighbors);
193 
194  // search again, now that we know that there is enough room to store the results
195  point_cloud_search.generateCRNeighborListsFromRadiusSearch(false /* dry run */,
196  target_coords, neighbor_lists, number_neighbors_list, epsilon);
197 
198  auto nla(CreateNeighborLists(neighbor_lists, number_neighbors_list));
199 
200  double radius_search_time = timer.seconds();
201  printf("nanoflann search time: %f s\n", radius_search_time);
202 
203  // convert point cloud search to vector of maps
204  timer.reset();
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;
211  }
212  });
213  double point_cloud_convert_time = timer.seconds();
214  printf("point cloud convert time: %f s\n", point_cloud_convert_time);
215  Kokkos::fence();
216 
217  // loop over targets, finding all of their neighbors by brute force
218  timer.reset();
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) {
224  double dist = 0;
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));
227  }
228  dist = std::sqrt(dist);
229 
230  if (dist<(h_multiplier + random_vec(i))*h_spacing) {
231  brute_force_neighbor_list_i[j]=dist;
232  }
233  }
234  });
235  double brute_force_search_time = timer.seconds();
236  printf("brute force search time: %f s\n", brute_force_search_time);
237  Kokkos::fence();
238 
239  timer.reset();
240  // check that all neighbors in brute force list are in point cloud search list
241  int num_passed = 0;
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) {
247  all_found = false;
248  }
249  }
250  if (all_found) t_num_passed++;
251  }, Kokkos::Sum<int>(num_passed));
252  Kokkos::fence();
253  if (num_passed != number_target_coords) {
254  printf("Brute force neighbor not found in point cloud list\n");
255  all_passed = false;
256  }
257 
258  num_passed = 0;
259  // check that all neighbors in point cloud search list are in brute force list
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) {
265  all_found = false;
266  }
267  }
268  if (all_found) t_num_passed++;
269  }, Kokkos::Sum<int>(num_passed));
270  Kokkos::fence();
271  if (num_passed != number_target_coords) {
272  printf("Point cloud neighbor not found in brute force list\n");
273  all_passed = false;
274  }
275 
276  double check_time = timer.seconds();
277  printf("check time: %f s\n", check_time);
278  }
279  timer.reset();
280  { // 2D neighbor lists
281  printf("\n2D neighbor lists:\n");
282  // Point cloud construction for neighbor search
283  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
284  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
285 
286  Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
287  number_target_coords, 2); // first column is # of neighbors
288 
289  // each target site has a window size
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);
292 
293  Kokkos::parallel_for("random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
294  (0,number_target_coords), [&](const int i) {
295  // value in [-.25, .25]
296  random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
297  epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
298  });
299  Kokkos::fence();
300 
301  // query the point cloud to generate the neighbor lists using a radius search
302 
303  // start with a dry-run, but without enough room to store the results
304  auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(true /* dry run */,
305  target_coords, neighbor_lists, epsilon);
306  printf("max num neighbors: %lu\n", max_num_neighbors);
307 
308  // resize neighbor lists to be large enough to hold the results
309  neighbor_lists = Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace>("neighbor lists",
310  number_target_coords, 1+max_num_neighbors); // first column is # of neighbors
311 
312  // search again, now that we know that there is enough room to store the results
313  max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(false /* dry run */,
314  target_coords, neighbor_lists, epsilon);
315 
316  double radius_search_time = timer.seconds();
317  printf("nanoflann search time: %f s\n", radius_search_time);
318 
319  // convert point cloud search to vector of maps
320  timer.reset();
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;
327  }
328  });
329  double point_cloud_convert_time = timer.seconds();
330  printf("point cloud convert time: %f s\n", point_cloud_convert_time);
331  Kokkos::fence();
332 
333  // loop over targets, finding all of their neighbors by brute force
334  timer.reset();
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) {
340  double dist = 0;
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));
343  }
344  dist = std::sqrt(dist);
345 
346  if (dist<(h_multiplier + random_vec(i))*h_spacing) {
347  brute_force_neighbor_list_i[j]=dist;
348  }
349  }
350  });
351  double brute_force_search_time = timer.seconds();
352  printf("brute force search time: %f s\n", brute_force_search_time);
353  Kokkos::fence();
354 
355  timer.reset();
356  // check that all neighbors in brute force list are in point cloud search list
357  int num_passed = 0;
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) {
363  all_found = false;
364  }
365  }
366  if (all_found) t_num_passed++;
367  }, Kokkos::Sum<int>(num_passed));
368  Kokkos::fence();
369  if (num_passed != number_target_coords) {
370  printf("Brute force neighbor not found in point cloud list\n");
371  all_passed = false;
372  }
373 
374  num_passed = 0;
375  // check that all neighbors in point cloud search list are in brute force list
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) {
381  all_found = false;
382  }
383  }
384  if (all_found) t_num_passed++;
385  }, Kokkos::Sum<int>(num_passed));
386  Kokkos::fence();
387  if (num_passed != number_target_coords) {
388  printf("Point cloud neighbor not found in brute force list\n");
389  all_passed = false;
390  }
391 
392  double check_time = timer.seconds();
393  printf("check time: %f s\n", check_time);
394  }
395  timer.reset();
396  { // convert 2D neighbor lists to 1D compressed row neighbor lists
397  printf("\n2D neighbor lists converted to 1D compressed row neighbor lists:\n");
398  // Point cloud construction for neighbor search
399  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
400  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
401 
402  Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
403  number_target_coords, 2); // first column is # of neighbors
404 
405  // each target site has a window size
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);
408 
409  Kokkos::parallel_for("random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
410  (0,number_target_coords), [&](const int i) {
411  // value in [-.25, .25]
412  random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
413  epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
414  });
415  Kokkos::fence();
416 
417  // query the point cloud to generate the neighbor lists using a radius search
418 
419  // start with a dry-run, but without enough room to store the results
420  auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(true /* dry run */,
421  target_coords, neighbor_lists, epsilon);
422  printf("max num neighbors: %lu\n", max_num_neighbors);
423 
424  // resize neighbor lists to be large enough to hold the results
425  neighbor_lists = Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace>("neighbor lists",
426  number_target_coords, 1+max_num_neighbors); // first column is # of neighbors
427 
428  // search again, now that we know that there is enough room to store the results
429  max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(false /* dry run */,
430  target_coords, neighbor_lists, epsilon);
431 
432  auto nla = Convert2DToCompressedRowNeighborLists(neighbor_lists);
433 
434  double radius_search_time = timer.seconds();
435  printf("nanoflann search time: %f s\n", radius_search_time);
436 
437  // convert point cloud search to vector of maps
438  timer.reset();
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;
445  }
446  });
447  double point_cloud_convert_time = timer.seconds();
448  printf("point cloud convert time: %f s\n", point_cloud_convert_time);
449  Kokkos::fence();
450 
451  // loop over targets, finding all of their neighbors by brute force
452  timer.reset();
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) {
458  double dist = 0;
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));
461  }
462  dist = std::sqrt(dist);
463 
464  if (dist<(h_multiplier + random_vec(i))*h_spacing) {
465  brute_force_neighbor_list_i[j]=dist;
466  }
467  }
468  });
469  double brute_force_search_time = timer.seconds();
470  printf("brute force search time: %f s\n", brute_force_search_time);
471  Kokkos::fence();
472 
473  timer.reset();
474  // check that all neighbors in brute force list are in point cloud search list
475  int num_passed = 0;
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) {
481  all_found = false;
482  }
483  }
484  if (all_found) t_num_passed++;
485  }, Kokkos::Sum<int>(num_passed));
486  Kokkos::fence();
487  if (num_passed != number_target_coords) {
488  printf("Brute force neighbor not found in point cloud list\n");
489  all_passed = false;
490  }
491 
492  num_passed = 0;
493  // check that all neighbors in point cloud search list are in brute force list
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) {
499  all_found = false;
500  }
501  }
502  if (all_found) t_num_passed++;
503  }, Kokkos::Sum<int>(num_passed));
504  Kokkos::fence();
505  if (num_passed != number_target_coords) {
506  printf("Point cloud neighbor not found in brute force list\n");
507  all_passed = false;
508  }
509 
510  double check_time = timer.seconds();
511  printf("check time: %f s\n", check_time);
512  }
513  }
514  else if (do_knn_search) {
515 
516  // do knn search
517  timer.reset();
518  { // 1D compressed row neighbor lists
519  printf("\n1D compressed row neighbor lists:\n");
520  // Point cloud construction for neighbor search
521  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
522  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
523 
524  const int min_neighbors = Compadre::GMLS::getNP(3, dimension);
525  Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
526  number_target_coords); // first column is # of neighbors
527  Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> number_neighbors_list("number of neighbors list",
528  number_target_coords); // first column is # of neighbors
529 
530  // each target site has a window size
531  Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
532 
533  // query the point cloud to generate the neighbor lists using a KNN search
534 
535  // start with a dry-run, but without enough room to store the results
536  auto total_num_neighbors = point_cloud_search.generateCRNeighborListsFromKNNSearch(true /* dry-run for sizes */,
537  target_coords, neighbor_lists, number_neighbors_list, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
538  printf("total num neighbors: %lu\n", total_num_neighbors);
539 
540  // resize with room to store results
541  Kokkos::resize(neighbor_lists, total_num_neighbors);
542 
543  // real knn search with space to store
544  point_cloud_search.generateCRNeighborListsFromKNNSearch(false /*not dry run*/,
545  target_coords, neighbor_lists, number_neighbors_list, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
546 
547  auto nla(CreateNeighborLists(neighbor_lists, number_neighbors_list));
548 
549  // convert point cloud search to vector of maps
550  timer.reset();
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) {
556  double dist = 0;
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));
559  }
560  dist = std::sqrt(dist);
561  point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = dist;
562  }
563  });
564  double point_cloud_convert_time = timer.seconds();
565  printf("point cloud convert time: %f s\n", point_cloud_convert_time);
566  Kokkos::fence();
567 
568  double eps = 1e-14;
569  // need to do a sort by dist, then verify epsilon(i) = 1.5 * knn_distance
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);
575  int j = 0;
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++;
580  }
581  break;
582  }
583  j++;
584  }
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");
588  all_passed = false;
589  }
590 
591  // loop over targets, finding knn by brute force
592  // point cloud multiplied kth by 1.5, but we are searching by brute force to ensure that for the distance h,
593  // we find at least as many neighbors as min_neighbors
594  // (not for the 1.5x as many, but rather for 1.0x as many)
595  timer.reset();
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) {
602  double dist = 0;
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));
605  }
606  dist = std::sqrt(dist);
607 
608  // epsilon was set to 1.5 * distance to kth neighbor, but we need distance to kth neighbor
609  if (dist<(epsilon(i)/1.5 + eps)) {
610  brute_force_neighbor_list_i[j]=dist;
611  }
612  }
613  // verify k neighbors found
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");
618  all_passed = false;
619  }
620  double brute_force_search_time = timer.seconds();
621  printf("brute force search time: %f s\n", brute_force_search_time);
622  Kokkos::fence();
623  }
624  timer.reset();
625  { // 2D neighbor lists
626  printf("\n2D neighbor lists:\n");
627  // Point cloud construction for neighbor search
628  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
629  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
630 
631  const int min_neighbors = Compadre::GMLS::getNP(3, dimension);
632  Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
633  number_target_coords, 1+min_neighbors); // first column is # of neighbors
634 
635  // each target site has a window size
636  Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
637 
638  // query the point cloud to generate the neighbor lists using a KNN search
639  // start with a dry-run, but without enough room to store the results
640  auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(true /* dry-run for sizes */,
641  target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
642  printf("max num neighbors: %lu\n", max_num_neighbors);
643 
644  // resize with room to store results
645  Kokkos::resize(neighbor_lists, neighbor_lists.extent(0), 1+max_num_neighbors);
646 
647  // real knn search with space to store
648  max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(false /*not dry run*/,
649  target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
650 
651  // convert point cloud search to vector of maps
652  timer.reset();
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) {
658  double dist = 0;
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));
661  }
662  dist = std::sqrt(dist);
663  point_cloud_neighbor_list_i[neighbor_lists(i,j)] = dist;
664  }
665  });
666  double point_cloud_convert_time = timer.seconds();
667  printf("point cloud convert time: %f s\n", point_cloud_convert_time);
668  Kokkos::fence();
669 
670  double eps = 1e-14;
671  // need to do a sort by dist, then verify epsilon(i) = 1.5 * knn_distance
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);
677  int j = 0;
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++;
682  }
683  break;
684  }
685  j++;
686  }
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");
690  all_passed = false;
691  }
692 
693  // loop over targets, finding knn by brute force
694  // point cloud multiplied kth by 1.5, but we are searching by brute force to ensure that for the distance h,
695  // we find at least as many neighbors as min_neighbors
696  // (not for the 1.5x as many, but rather for 1.0x as many)
697  timer.reset();
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) {
704  double dist = 0;
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));
707  }
708  dist = std::sqrt(dist);
709 
710  // epsilon was set to 1.5 * distance to kth neighbor, but we need distance to kth neighbor
711  if (dist<(epsilon(i)/1.5 + eps)) {
712  brute_force_neighbor_list_i[j]=dist;
713  }
714  }
715  // verify k neighbors found
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");
720  all_passed = false;
721  }
722  double brute_force_search_time = timer.seconds();
723  printf("brute force search time: %f s\n", brute_force_search_time);
724  Kokkos::fence();
725  }
726  timer.reset();
727  { // convert 2D neighbor lists to 1D compressed row neighbor lists
728  printf("\n2D neighbor lists converted to 1D compressed row neighbor lists:\n");
729  // Point cloud construction for neighbor search
730  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
731  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
732 
733  const int min_neighbors = Compadre::GMLS::getNP(3, dimension);
734  Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
735  number_target_coords, 1+min_neighbors); // first column is # of neighbors
736 
737  // each target site has a window size
738  Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
739 
740  // query the point cloud to generate the neighbor lists using a KNN search
741  // start with a dry-run, but without enough room to store the results
742  auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(true /* dry-run for sizes */,
743  target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
744  printf("max num neighbors: %lu\n", max_num_neighbors);
745 
746  // resize with room to store results
747  Kokkos::resize(neighbor_lists, neighbor_lists.extent(0), 1+max_num_neighbors);
748 
749  // real knn search with space to store
750  max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(false /*not dry run*/,
751  target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
752 
753  auto nla = Convert2DToCompressedRowNeighborLists(neighbor_lists);
754 
755  // convert point cloud search to vector of maps
756  timer.reset();
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) {
762  double dist = 0;
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));
765  }
766  dist = std::sqrt(dist);
767  point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = dist;
768  }
769  });
770  double point_cloud_convert_time = timer.seconds();
771  printf("point cloud convert time: %f s\n", point_cloud_convert_time);
772  Kokkos::fence();
773 
774  double eps = 1e-14;
775  // need to do a sort by dist, then verify epsilon(i) = 1.5 * knn_distance
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);
781  int j = 0;
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++;
786  }
787  break;
788  }
789  j++;
790  }
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");
794  all_passed = false;
795  }
796 
797  // loop over targets, finding knn by brute force
798  // point cloud multiplied kth by 1.5, but we are searching by brute force to ensure that for the distance h,
799  // we find at least as many neighbors as min_neighbors
800  // (not for the 1.5x as many, but rather for 1.0x as many)
801  timer.reset();
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) {
808  double dist = 0;
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));
811  }
812  dist = std::sqrt(dist);
813 
814  // epsilon was set to 1.5 * distance to kth neighbor, but we need distance to kth neighbor
815  if (dist<(epsilon(i)/1.5 + eps)) {
816  brute_force_neighbor_list_i[j]=dist;
817  }
818  }
819  // verify k neighbors found
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");
824  all_passed = false;
825  }
826  double brute_force_search_time = timer.seconds();
827  printf("brute force search time: %f s\n", brute_force_search_time);
828  Kokkos::fence();
829  }
830  } else {
831  all_passed = false;
832  printf("do_radius_search and do_knn_search both false. Invalid combination.\n");
833  }
834 }
835 
836 // finalize Kokkos and MPI (if available)
837 kp.finalize();
838 #ifdef COMPADRE_USE_MPI
839 MPI_Finalize();
840 #endif
841 
842 // output to user that test passed or failed
843 if(all_passed) {
844  fprintf(stdout, "Passed test \n");
845  return 0;
846 } else {
847  fprintf(stdout, "Failed test \n");
848  return -1;
849 }
850 
851 } // main
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]
Definition: GMLS_Device.cpp:29
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.