Compadre 1.6.4
Loading...
Searching...
No Matches
NeighborSearchTest.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Compadre: COMpatible PArticle Discretization and REmap Toolkit
4//
5// Copyright 2018 NTESS and the Compadre contributors.
6// SPDX-License-Identifier: BSD-2-Clause
7// *****************************************************************************
8// @HEADER
9#include <vector>
10#include <map>
11#include <random>
12
13#include <Compadre_Config.h>
14#include <Compadre_GMLS.hpp>
16
17#ifdef COMPADRE_USE_MPI
18#include <mpi.h>
19#endif
20
21template<typename T1, typename T2>
22std::pair<T2,T1> swap_pair(const std::pair<T1,T2> &item)
23{
24 return std::pair<T2,T1>(item.second, item.first);
25}
26
27template<typename T1, typename T2>
28std::multimap<T2,T1> invert_map(const std::map<T1,T2> &original_map)
29{
30 std::multimap<T2,T1> new_map;
31 std::transform(original_map.begin(), original_map.end(), std::inserter(new_map, new_map.begin()), swap_pair<T1,T2>);
32 return new_map;
33}
34
35using namespace Compadre;
36
37int main (int argc, char* args[]) {
38
39// initializes MPI (if available) with command line arguments given
40#ifdef COMPADRE_USE_MPI
41MPI_Init(&argc, &args);
42#endif
43
44// initializes Kokkos with command line arguments given
45Kokkos::initialize(argc, args);
46
47// becomes false if the computed solution not within the failure_threshold of the actual solution
48bool all_passed = true;
49
50{
51 // seed random generator
52 srand(1234);
53 Kokkos::Timer timer;
54
55 int search_type = 0; // 0 - radius, 1 - knn
56 if (argc >= 5) {
57 int arg5toi = atoi(args[4]);
58 if (arg5toi >= 0) {
59 search_type = arg5toi;
60 }
61 }
62 bool do_radius_search = (search_type==0);
63 bool do_knn_search = (search_type==1);
64 printf("do_radius_search: %d\n", do_radius_search);
65 printf("do_knn_search: %d\n", do_knn_search);
66
67 // check if 4 arguments are given from the command line
68 // multiplier times h spacing for search
69 double h_multiplier = 3.0; // dimension 3 by default
70 if (argc >= 4) {
71 double arg4tof = atof(args[3]);
72 if (arg4tof > 0) {
73 h_multiplier = arg4tof;
74 }
75 }
76
77 // check if 3 arguments are given from the command line
78 // set the number of target sites where we will reconstruct the target functionals at
79 int number_target_coords = 200; // 200 target sites by default
80 if (argc >= 3) {
81 int arg3toi = atoi(args[2]);
82 if (arg3toi > 0) {
83 number_target_coords = arg3toi;
84 }
85 }
86
87 // check if 2 arguments are given from the command line
88 // set the number of dimensions for points and the search
89 int dimension = 3; // 3D by default
90 if (argc >= 2) {
91 int arg2toi = atoi(args[1]);
92 if (arg2toi > 0) {
93 dimension = arg2toi;
94 }
95 }
96
97 //// minimum neighbors for unisolvency is the same as the size of the polynomial basis
98 //const int min_neighbors = Compadre::GMLS::getNP(order, dimension);
99
100 // approximate spacing of source sites
101 double h_spacing = 1./std::pow(number_target_coords, 0.50 + 0.2*(dimension==2));
102 int n_neg1_to_1 = 2*(1/h_spacing) + 1; // always odd
103
104 // number of source coordinate sites that will fill a box of [-1,1]x[-1,1]x[-1,1] with a spacing approximately h
105 int number_source_coords = (n_neg1_to_1+1)*(n_neg1_to_1+1);
106 if (dimension==3) number_source_coords *= (n_neg1_to_1+1);
107 printf("target coords: %d\n", number_target_coords);
108 printf("source coords: %d\n", number_source_coords);
109
110 // coordinates of source sites
111 Kokkos::View<double**, Kokkos::DefaultHostExecutionSpace> source_coords("source coordinates",
112 number_source_coords, 3);
113
114 // coordinates of target sites
115 Kokkos::View<double**, Kokkos::DefaultHostExecutionSpace> target_coords("target coordinates", number_target_coords, 3);
116
117 // fill source coordinates with a uniform grid
118 int source_index = 0;
119 double this_coord[3] = {0,0,0};
120 for (int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
121 this_coord[0] = i*h_spacing;
122 for (int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
123 this_coord[1] = j*h_spacing;
124 for (int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
125 this_coord[2] = k*h_spacing;
126 if (dimension==3) {
127 source_coords(source_index,0) = this_coord[0];
128 source_coords(source_index,1) = this_coord[1];
129 source_coords(source_index,2) = this_coord[2];
130 source_index++;
131 }
132 }
133 if (dimension==2) {
134 source_coords(source_index,0) = this_coord[0];
135 source_coords(source_index,1) = this_coord[1];
136 source_coords(source_index,2) = 0;
137 source_index++;
138 }
139 }
140 if (dimension==1) {
141 source_coords(source_index,0) = this_coord[0];
142 source_coords(source_index,1) = 0;
143 source_coords(source_index,2) = 0;
144 source_index++;
145 }
146 }
147
148 // fill target coords somewhere inside of [-0.5,0.5]x[-0.5,0.5]x[-0.5,0.5]
149 for(int i=0; i<number_target_coords; i++){
150
151 // first, we get a uniformly random distributed direction
152 double rand_dir[3] = {0,0,0};
153
154 for (int j=0; j<dimension; ++j) {
155 // rand_dir[j] is in [-0.5, 0.5]
156 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
157 }
158
159 // then we get a uniformly random radius
160 for (int j=0; j<dimension; ++j) {
161 target_coords(i,j) = rand_dir[j];
162 }
163
164 }
165
166 if (do_radius_search) {
167 timer.reset();
168 { // 1D compressed row neighbor lists
169 printf("\n1D compressed row neighbor lists:\n");
170 // Point cloud construction for neighbor search
171 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
172 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
173
174 Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
175 number_target_coords); // first column is # of neighbors
176 Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> number_neighbors_list("number of neighbors list",
177 number_target_coords); // first column is # of neighbors
178
179 // each target site has a window size
180 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
181 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> random_vec("random vector", number_target_coords);
182
183 Kokkos::parallel_for("random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
184 (0,number_target_coords), [&](const int i) {
185 // value in [-.25, .25]
186 random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
187 epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
188 });
189 Kokkos::fence();
190
191 // query the point cloud to generate the neighbor lists using a radius search
192
193 // start with a dry-run, but without enough room to store the results
194 size_t total_num_neighbors = point_cloud_search.generateCRNeighborListsFromRadiusSearch(true /* dry run */,
195 target_coords, neighbor_lists, number_neighbors_list, epsilon);
196 printf("total num neighbors: %lu\n", total_num_neighbors);
197
198 // resize neighbor lists to be large enough to hold the results
199 Kokkos::resize(neighbor_lists, total_num_neighbors);
200
201 // search again, now that we know that there is enough room to store the results
202 point_cloud_search.generateCRNeighborListsFromRadiusSearch(false /* dry run */,
203 target_coords, neighbor_lists, number_neighbors_list, epsilon);
204
205 auto nla(CreateNeighborLists(neighbor_lists, number_neighbors_list));
206
207 double radius_search_time = timer.seconds();
208 printf("nanoflann search time: %f s\n", radius_search_time);
209
210 // convert point cloud search to vector of maps
211 timer.reset();
212 std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
213 Kokkos::parallel_for("point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
214 (0,number_target_coords), [&](const int i) {
215 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
216 for (int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
217 point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = 1.0;
218 }
219 });
220 double point_cloud_convert_time = timer.seconds();
221 printf("point cloud convert time: %f s\n", point_cloud_convert_time);
222 Kokkos::fence();
223
224 // loop over targets, finding all of their neighbors by brute force
225 timer.reset();
226 std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
227 Kokkos::parallel_for("brute force radius search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
228 (0,number_target_coords), [&](const int i) {
229 auto &brute_force_neighbor_list_i = const_cast<std::map<int, double>& >(brute_force_neighbor_list[i]);
230 for (int j=0; j<number_source_coords; ++j) {
231 double dist = 0;
232 for (int k=0; k<dimension; ++k) {
233 dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
234 }
235 dist = std::sqrt(dist);
236
237 if (dist<(h_multiplier + random_vec(i))*h_spacing) {
238 brute_force_neighbor_list_i[j]=dist;
239 }
240 }
241 });
242 double brute_force_search_time = timer.seconds();
243 printf("brute force search time: %f s\n", brute_force_search_time);
244 Kokkos::fence();
245
246 timer.reset();
247 // check that all neighbors in brute force list are in point cloud search list
248 int num_passed = 0;
249 Kokkos::parallel_reduce("check brute force search in point cloud search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
250 (0,number_target_coords), [&](const int i, int& t_num_passed) {
251 bool all_found = true;
252 for (auto it=brute_force_neighbor_list[i].begin(); it!=brute_force_neighbor_list[i].end(); ++it) {
253 if (point_cloud_neighbor_list[i].count(it->first)!=1) {
254 all_found = false;
255 }
256 }
257 if (all_found) t_num_passed++;
258 }, Kokkos::Sum<int>(num_passed));
259 Kokkos::fence();
260 if (num_passed != number_target_coords) {
261 printf("Brute force neighbor not found in point cloud list\n");
262 all_passed = false;
263 }
264
265 num_passed = 0;
266 // check that all neighbors in point cloud search list are in brute force list
267 Kokkos::parallel_reduce("original in brute force search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
268 (0,number_target_coords), KOKKOS_LAMBDA(const int i, int& t_num_passed) {
269 bool all_found = true;
270 for (int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
271 if (brute_force_neighbor_list[i].count(nla.getNeighborHost(i,j))!=1) {
272 all_found = false;
273 }
274 }
275 if (all_found) t_num_passed++;
276 }, Kokkos::Sum<int>(num_passed));
277 Kokkos::fence();
278 if (num_passed != number_target_coords) {
279 printf("Point cloud neighbor not found in brute force list\n");
280 all_passed = false;
281 }
282
283 double check_time = timer.seconds();
284 printf("check time: %f s\n", check_time);
285 }
286 timer.reset();
287 { // 2D neighbor lists
288 printf("\n2D neighbor lists:\n");
289 // Point cloud construction for neighbor search
290 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
291 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
292
293 Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
294 number_target_coords, 2); // first column is # of neighbors
295
296 // each target site has a window size
297 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
298 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> random_vec("random vector", number_target_coords);
299
300 Kokkos::parallel_for("random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
301 (0,number_target_coords), [&](const int i) {
302 // value in [-.25, .25]
303 random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
304 epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
305 });
306 Kokkos::fence();
307
308 // query the point cloud to generate the neighbor lists using a radius search
309
310 // start with a dry-run, but without enough room to store the results
311 auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(true /* dry run */,
312 target_coords, neighbor_lists, epsilon);
313 printf("max num neighbors: %lu\n", max_num_neighbors);
314
315 // resize neighbor lists to be large enough to hold the results
316 neighbor_lists = Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace>("neighbor lists",
317 number_target_coords, 1+max_num_neighbors); // first column is # of neighbors
318
319 // search again, now that we know that there is enough room to store the results
320 max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(false /* dry run */,
321 target_coords, neighbor_lists, epsilon);
322
323 double radius_search_time = timer.seconds();
324 printf("nanoflann search time: %f s\n", radius_search_time);
325
326 // convert point cloud search to vector of maps
327 timer.reset();
328 std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
329 Kokkos::parallel_for("point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
330 (0,number_target_coords), [&](const int i) {
331 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
332 for (int j=1; j<=neighbor_lists(i,0); ++j) {
333 point_cloud_neighbor_list_i[neighbor_lists(i,j)] = 1.0;
334 }
335 });
336 double point_cloud_convert_time = timer.seconds();
337 printf("point cloud convert time: %f s\n", point_cloud_convert_time);
338 Kokkos::fence();
339
340 // loop over targets, finding all of their neighbors by brute force
341 timer.reset();
342 std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
343 Kokkos::parallel_for("brute force radius search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
344 (0,number_target_coords), [&](const int i) {
345 auto &brute_force_neighbor_list_i = const_cast<std::map<int, double>& >(brute_force_neighbor_list[i]);
346 for (int j=0; j<number_source_coords; ++j) {
347 double dist = 0;
348 for (int k=0; k<dimension; ++k) {
349 dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
350 }
351 dist = std::sqrt(dist);
352
353 if (dist<(h_multiplier + random_vec(i))*h_spacing) {
354 brute_force_neighbor_list_i[j]=dist;
355 }
356 }
357 });
358 double brute_force_search_time = timer.seconds();
359 printf("brute force search time: %f s\n", brute_force_search_time);
360 Kokkos::fence();
361
362 timer.reset();
363 // check that all neighbors in brute force list are in point cloud search list
364 int num_passed = 0;
365 Kokkos::parallel_reduce("check brute force search in point cloud search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
366 (0,number_target_coords), [&](const int i, int& t_num_passed) {
367 bool all_found = true;
368 for (auto it=brute_force_neighbor_list[i].begin(); it!=brute_force_neighbor_list[i].end(); ++it) {
369 if (point_cloud_neighbor_list[i].count(it->first)!=1) {
370 all_found = false;
371 }
372 }
373 if (all_found) t_num_passed++;
374 }, Kokkos::Sum<int>(num_passed));
375 Kokkos::fence();
376 if (num_passed != number_target_coords) {
377 printf("Brute force neighbor not found in point cloud list\n");
378 all_passed = false;
379 }
380
381 num_passed = 0;
382 // check that all neighbors in point cloud search list are in brute force list
383 Kokkos::parallel_reduce("original in brute force search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
384 (0,number_target_coords), KOKKOS_LAMBDA(const int i, int& t_num_passed) {
385 bool all_found = true;
386 for (int j=1; j<=neighbor_lists(i,0); ++j) {
387 if (brute_force_neighbor_list[i].count(neighbor_lists(i,j))!=1) {
388 all_found = false;
389 }
390 }
391 if (all_found) t_num_passed++;
392 }, Kokkos::Sum<int>(num_passed));
393 Kokkos::fence();
394 if (num_passed != number_target_coords) {
395 printf("Point cloud neighbor not found in brute force list\n");
396 all_passed = false;
397 }
398
399 double check_time = timer.seconds();
400 printf("check time: %f s\n", check_time);
401 }
402 timer.reset();
403 { // convert 2D neighbor lists to 1D compressed row neighbor lists
404 printf("\n2D neighbor lists converted to 1D compressed row neighbor lists:\n");
405 // Point cloud construction for neighbor search
406 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
407 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
408
409 Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
410 number_target_coords, 2); // first column is # of neighbors
411
412 // each target site has a window size
413 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
414 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> random_vec("random vector", number_target_coords);
415
416 Kokkos::parallel_for("random perturbation of search size", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
417 (0,number_target_coords), [&](const int i) {
418 // value in [-.25, .25]
419 random_vec(i) = 0.5 * ((double)rand() / (double) RAND_MAX) - 0.5;
420 epsilon(i) = (h_multiplier + random_vec(i)) * h_spacing;
421 });
422 Kokkos::fence();
423
424 // query the point cloud to generate the neighbor lists using a radius search
425
426 // start with a dry-run, but without enough room to store the results
427 auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(true /* dry run */,
428 target_coords, neighbor_lists, epsilon);
429 printf("max num neighbors: %lu\n", max_num_neighbors);
430
431 // resize neighbor lists to be large enough to hold the results
432 neighbor_lists = Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace>("neighbor lists",
433 number_target_coords, 1+max_num_neighbors); // first column is # of neighbors
434
435 // search again, now that we know that there is enough room to store the results
436 max_num_neighbors = point_cloud_search.generate2DNeighborListsFromRadiusSearch(false /* dry run */,
437 target_coords, neighbor_lists, epsilon);
438
439 auto nla = Convert2DToCompressedRowNeighborLists(neighbor_lists);
440
441 double radius_search_time = timer.seconds();
442 printf("nanoflann search time: %f s\n", radius_search_time);
443
444 // convert point cloud search to vector of maps
445 timer.reset();
446 std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
447 Kokkos::parallel_for("point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
448 (0,number_target_coords), [&](const int i) {
449 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
450 for (int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
451 point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = 1.0;
452 }
453 });
454 double point_cloud_convert_time = timer.seconds();
455 printf("point cloud convert time: %f s\n", point_cloud_convert_time);
456 Kokkos::fence();
457
458 // loop over targets, finding all of their neighbors by brute force
459 timer.reset();
460 std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
461 Kokkos::parallel_for("brute force radius search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
462 (0,number_target_coords), [&](const int i) {
463 auto &brute_force_neighbor_list_i = const_cast<std::map<int, double>& >(brute_force_neighbor_list[i]);
464 for (int j=0; j<number_source_coords; ++j) {
465 double dist = 0;
466 for (int k=0; k<dimension; ++k) {
467 dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
468 }
469 dist = std::sqrt(dist);
470
471 if (dist<(h_multiplier + random_vec(i))*h_spacing) {
472 brute_force_neighbor_list_i[j]=dist;
473 }
474 }
475 });
476 double brute_force_search_time = timer.seconds();
477 printf("brute force search time: %f s\n", brute_force_search_time);
478 Kokkos::fence();
479
480 timer.reset();
481 // check that all neighbors in brute force list are in point cloud search list
482 int num_passed = 0;
483 Kokkos::parallel_reduce("check brute force search in point cloud search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
484 (0,number_target_coords), [&](const int i, int& t_num_passed) {
485 bool all_found = true;
486 for (auto it=brute_force_neighbor_list[i].begin(); it!=brute_force_neighbor_list[i].end(); ++it) {
487 if (point_cloud_neighbor_list[i].count(it->first)!=1) {
488 all_found = false;
489 }
490 }
491 if (all_found) t_num_passed++;
492 }, Kokkos::Sum<int>(num_passed));
493 Kokkos::fence();
494 if (num_passed != number_target_coords) {
495 printf("Brute force neighbor not found in point cloud list\n");
496 all_passed = false;
497 }
498
499 num_passed = 0;
500 // check that all neighbors in point cloud search list are in brute force list
501 Kokkos::parallel_reduce("original in brute force search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
502 (0,number_target_coords), KOKKOS_LAMBDA(const int i, int& t_num_passed) {
503 bool all_found = true;
504 for (int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
505 if (brute_force_neighbor_list[i].count(nla.getNeighborHost(i,j))!=1) {
506 all_found = false;
507 }
508 }
509 if (all_found) t_num_passed++;
510 }, Kokkos::Sum<int>(num_passed));
511 Kokkos::fence();
512 if (num_passed != number_target_coords) {
513 printf("Point cloud neighbor not found in brute force list\n");
514 all_passed = false;
515 }
516
517 double check_time = timer.seconds();
518 printf("check time: %f s\n", check_time);
519 }
520 }
521 else if (do_knn_search) {
522
523 // do knn search
524 timer.reset();
525 { // 1D compressed row neighbor lists
526 printf("\n1D compressed row neighbor lists:\n");
527 // Point cloud construction for neighbor search
528 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
529 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
530
531 const int min_neighbors = Compadre::GMLS::getNP(3, dimension);
532 Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
533 number_target_coords); // first column is # of neighbors
534 Kokkos::View<int*, Kokkos::DefaultHostExecutionSpace> number_neighbors_list("number of neighbors list",
535 number_target_coords); // first column is # of neighbors
536
537 // each target site has a window size
538 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
539
540 // query the point cloud to generate the neighbor lists using a KNN search
541
542 // start with a dry-run, but without enough room to store the results
543 auto total_num_neighbors = point_cloud_search.generateCRNeighborListsFromKNNSearch(true /* dry-run for sizes */,
544 target_coords, neighbor_lists, number_neighbors_list, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
545 printf("total num neighbors: %lu\n", total_num_neighbors);
546
547 // resize with room to store results
548 Kokkos::resize(neighbor_lists, total_num_neighbors);
549
550 // real knn search with space to store
551 point_cloud_search.generateCRNeighborListsFromKNNSearch(false /*not dry run*/,
552 target_coords, neighbor_lists, number_neighbors_list, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
553
554 auto nla(CreateNeighborLists(neighbor_lists, number_neighbors_list));
555
556 // convert point cloud search to vector of maps
557 timer.reset();
558 std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
559 Kokkos::parallel_for("point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
560 (0,number_target_coords), [&](const int i) {
561 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
562 for (int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
563 double dist = 0;
564 for (int k=0; k<dimension; ++k) {
565 dist += (target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k))*(target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k));
566 }
567 dist = std::sqrt(dist);
568 point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = dist;
569 }
570 });
571 double point_cloud_convert_time = timer.seconds();
572 printf("point cloud convert time: %f s\n", point_cloud_convert_time);
573 Kokkos::fence();
574
575 double eps = 1e-14;
576 // need to do a sort by dist, then verify epsilon(i) = 1.5 * knn_distance
577 int epsilon_diff_from_knn_dist_multiplied = 0;
578 Kokkos::parallel_reduce("check kth neighbor in point cloud", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
579 (0,number_target_coords), [&](const int i, int& t_epsilon_diff_from_knn_dist_multiplied) {
580 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
581 std::multimap<double, int> inverted_map = invert_map(point_cloud_neighbor_list_i);
582 int j = 0;
583 for (auto it=inverted_map.begin(); it!=inverted_map.end(); ++it) {
584 if (j==min_neighbors-1) {
585 if (std::abs(epsilon(i)/1.5 - it->first) > eps) {
586 t_epsilon_diff_from_knn_dist_multiplied++;
587 }
588 break;
589 }
590 j++;
591 }
592 }, Kokkos::Sum<int>(epsilon_diff_from_knn_dist_multiplied));
593 if (epsilon_diff_from_knn_dist_multiplied > 0) {
594 printf("1.5*kth_neighbor_distance differs from calculation in epsilon vector.\n");
595 all_passed = false;
596 }
597
598 // loop over targets, finding knn by brute force
599 // point cloud multiplied kth by 1.5, but we are searching by brute force to ensure that for the distance h,
600 // we find at least as many neighbors as min_neighbors
601 // (not for the 1.5x as many, but rather for 1.0x as many)
602 timer.reset();
603 int total_insufficient_neighbors = 0;
604 std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
605 Kokkos::parallel_reduce("brute force knn search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
606 (0,number_target_coords), [&](const int i, int& t_total_insufficient_neighbors) {
607 auto &brute_force_neighbor_list_i = const_cast<std::map<int, double>& >(brute_force_neighbor_list[i]);
608 for (int j=0; j<number_source_coords; ++j) {
609 double dist = 0;
610 for (int k=0; k<dimension; ++k) {
611 dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
612 }
613 dist = std::sqrt(dist);
614
615 // epsilon was set to 1.5 * distance to kth neighbor, but we need distance to kth neighbor
616 if (dist<(epsilon(i)/1.5 + eps)) {
617 brute_force_neighbor_list_i[j]=dist;
618 }
619 }
620 // verify k neighbors found
621 if (brute_force_neighbor_list_i.size() < (size_t)min_neighbors) t_total_insufficient_neighbors++;
622 }, Kokkos::Sum<int>(total_insufficient_neighbors));
623 if (total_insufficient_neighbors > 0) {
624 printf("less neighbors found using kth_neighbor_distance+eps than using knn search.\n");
625 all_passed = false;
626 }
627 double brute_force_search_time = timer.seconds();
628 printf("brute force search time: %f s\n", brute_force_search_time);
629 Kokkos::fence();
630 }
631 timer.reset();
632 { // 2D neighbor lists
633 printf("\n2D neighbor lists:\n");
634 // Point cloud construction for neighbor search
635 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
636 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
637
638 const int min_neighbors = Compadre::GMLS::getNP(3, dimension);
639 Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
640 number_target_coords, 1+min_neighbors); // first column is # of neighbors
641
642 // each target site has a window size
643 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
644
645 // query the point cloud to generate the neighbor lists using a KNN search
646 // start with a dry-run, but without enough room to store the results
647 auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(true /* dry-run for sizes */,
648 target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
649 printf("max num neighbors: %lu\n", max_num_neighbors);
650
651 // resize with room to store results
652 Kokkos::resize(neighbor_lists, neighbor_lists.extent(0), 1+max_num_neighbors);
653
654 // real knn search with space to store
655 max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(false /*not dry run*/,
656 target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
657
658 // convert point cloud search to vector of maps
659 timer.reset();
660 std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
661 Kokkos::parallel_for("point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
662 (0,number_target_coords), [&](const int i) {
663 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
664 for (int j=1; j<=neighbor_lists(i,0); ++j) {
665 double dist = 0;
666 for (int k=0; k<dimension; ++k) {
667 dist += (target_coords(i,k)-source_coords(neighbor_lists(i,j),k))*(target_coords(i,k)-source_coords(neighbor_lists(i,j),k));
668 }
669 dist = std::sqrt(dist);
670 point_cloud_neighbor_list_i[neighbor_lists(i,j)] = dist;
671 }
672 });
673 double point_cloud_convert_time = timer.seconds();
674 printf("point cloud convert time: %f s\n", point_cloud_convert_time);
675 Kokkos::fence();
676
677 double eps = 1e-14;
678 // need to do a sort by dist, then verify epsilon(i) = 1.5 * knn_distance
679 int epsilon_diff_from_knn_dist_multiplied = 0;
680 Kokkos::parallel_reduce("check kth neighbor in point cloud", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
681 (0,number_target_coords), [&](const int i, int& t_epsilon_diff_from_knn_dist_multiplied) {
682 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
683 std::multimap<double, int> inverted_map = invert_map(point_cloud_neighbor_list_i);
684 int j = 0;
685 for (auto it=inverted_map.begin(); it!=inverted_map.end(); ++it) {
686 if (j==min_neighbors-1) {
687 if (std::abs(epsilon(i)/1.5 - it->first) > eps) {
688 t_epsilon_diff_from_knn_dist_multiplied++;
689 }
690 break;
691 }
692 j++;
693 }
694 }, Kokkos::Sum<int>(epsilon_diff_from_knn_dist_multiplied));
695 if (epsilon_diff_from_knn_dist_multiplied > 0) {
696 printf("1.5*kth_neighbor_distance differs from calculation in epsilon vector.\n");
697 all_passed = false;
698 }
699
700 // loop over targets, finding knn by brute force
701 // point cloud multiplied kth by 1.5, but we are searching by brute force to ensure that for the distance h,
702 // we find at least as many neighbors as min_neighbors
703 // (not for the 1.5x as many, but rather for 1.0x as many)
704 timer.reset();
705 int total_insufficient_neighbors = 0;
706 std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
707 Kokkos::parallel_reduce("brute force knn search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
708 (0,number_target_coords), [&](const int i, int& t_total_insufficient_neighbors) {
709 auto &brute_force_neighbor_list_i = const_cast<std::map<int, double>& >(brute_force_neighbor_list[i]);
710 for (int j=0; j<number_source_coords; ++j) {
711 double dist = 0;
712 for (int k=0; k<dimension; ++k) {
713 dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
714 }
715 dist = std::sqrt(dist);
716
717 // epsilon was set to 1.5 * distance to kth neighbor, but we need distance to kth neighbor
718 if (dist<(epsilon(i)/1.5 + eps)) {
719 brute_force_neighbor_list_i[j]=dist;
720 }
721 }
722 // verify k neighbors found
723 if (brute_force_neighbor_list_i.size() < (size_t)min_neighbors) t_total_insufficient_neighbors++;
724 }, Kokkos::Sum<int>(total_insufficient_neighbors));
725 if (total_insufficient_neighbors > 0) {
726 printf("less neighbors found using kth_neighbor_distance+eps than using knn search.\n");
727 all_passed = false;
728 }
729 double brute_force_search_time = timer.seconds();
730 printf("brute force search time: %f s\n", brute_force_search_time);
731 Kokkos::fence();
732 }
733 timer.reset();
734 { // convert 2D neighbor lists to 1D compressed row neighbor lists
735 printf("\n2D neighbor lists converted to 1D compressed row neighbor lists:\n");
736 // Point cloud construction for neighbor search
737 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
738 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
739
740 const int min_neighbors = Compadre::GMLS::getNP(3, dimension);
741 Kokkos::View<int**, Kokkos::DefaultHostExecutionSpace> neighbor_lists("neighbor lists",
742 number_target_coords, 1+min_neighbors); // first column is # of neighbors
743
744 // each target site has a window size
745 Kokkos::View<double*, Kokkos::DefaultHostExecutionSpace> epsilon("h supports", number_target_coords);
746
747 // query the point cloud to generate the neighbor lists using a KNN search
748 // start with a dry-run, but without enough room to store the results
749 auto max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(true /* dry-run for sizes */,
750 target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
751 printf("max num neighbors: %lu\n", max_num_neighbors);
752
753 // resize with room to store results
754 Kokkos::resize(neighbor_lists, neighbor_lists.extent(0), 1+max_num_neighbors);
755
756 // real knn search with space to store
757 max_num_neighbors = point_cloud_search.generate2DNeighborListsFromKNNSearch(false /*not dry run*/,
758 target_coords, neighbor_lists, epsilon, min_neighbors, 1.5 /* cutoff_multiplier */);
759
760 auto nla = Convert2DToCompressedRowNeighborLists(neighbor_lists);
761
762 // convert point cloud search to vector of maps
763 timer.reset();
764 std::vector<std::map<int, double> > point_cloud_neighbor_list(number_target_coords, std::map<int, double>());
765 Kokkos::parallel_for("point search conversion", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
766 (0,number_target_coords), [&](const int i) {
767 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
768 for (int j=0; j<nla.getNumberOfNeighborsHost(i); ++j) {
769 double dist = 0;
770 for (int k=0; k<dimension; ++k) {
771 dist += (target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k))*(target_coords(i,k)-source_coords(nla.getNeighborHost(i,j),k));
772 }
773 dist = std::sqrt(dist);
774 point_cloud_neighbor_list_i[nla.getNeighborHost(i,j)] = dist;
775 }
776 });
777 double point_cloud_convert_time = timer.seconds();
778 printf("point cloud convert time: %f s\n", point_cloud_convert_time);
779 Kokkos::fence();
780
781 double eps = 1e-14;
782 // need to do a sort by dist, then verify epsilon(i) = 1.5 * knn_distance
783 int epsilon_diff_from_knn_dist_multiplied = 0;
784 Kokkos::parallel_reduce("check kth neighbor in point cloud", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
785 (0,number_target_coords), [&](const int i, int& t_epsilon_diff_from_knn_dist_multiplied) {
786 auto &point_cloud_neighbor_list_i = const_cast<std::map<int, double>& >(point_cloud_neighbor_list[i]);
787 std::multimap<double, int> inverted_map = invert_map(point_cloud_neighbor_list_i);
788 int j = 0;
789 for (auto it=inverted_map.begin(); it!=inverted_map.end(); ++it) {
790 if (j==min_neighbors-1) {
791 if (std::abs(epsilon(i)/1.5 - it->first) > eps) {
792 t_epsilon_diff_from_knn_dist_multiplied++;
793 }
794 break;
795 }
796 j++;
797 }
798 }, Kokkos::Sum<int>(epsilon_diff_from_knn_dist_multiplied));
799 if (epsilon_diff_from_knn_dist_multiplied > 0) {
800 printf("1.5*kth_neighbor_distance differs from calculation in epsilon vector.\n");
801 all_passed = false;
802 }
803
804 // loop over targets, finding knn by brute force
805 // point cloud multiplied kth by 1.5, but we are searching by brute force to ensure that for the distance h,
806 // we find at least as many neighbors as min_neighbors
807 // (not for the 1.5x as many, but rather for 1.0x as many)
808 timer.reset();
809 int total_insufficient_neighbors = 0;
810 std::vector<std::map<int, double> > brute_force_neighbor_list(number_target_coords, std::map<int, double>());
811 Kokkos::parallel_reduce("brute force knn search", Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>
812 (0,number_target_coords), [&](const int i, int& t_total_insufficient_neighbors) {
813 auto &brute_force_neighbor_list_i = const_cast<std::map<int, double>& >(brute_force_neighbor_list[i]);
814 for (int j=0; j<number_source_coords; ++j) {
815 double dist = 0;
816 for (int k=0; k<dimension; ++k) {
817 dist += (target_coords(i,k)-source_coords(j,k))*(target_coords(i,k)-source_coords(j,k));
818 }
819 dist = std::sqrt(dist);
820
821 // epsilon was set to 1.5 * distance to kth neighbor, but we need distance to kth neighbor
822 if (dist<(epsilon(i)/1.5 + eps)) {
823 brute_force_neighbor_list_i[j]=dist;
824 }
825 }
826 // verify k neighbors found
827 if (brute_force_neighbor_list_i.size() < (size_t)min_neighbors) t_total_insufficient_neighbors++;
828 }, Kokkos::Sum<int>(total_insufficient_neighbors));
829 if (total_insufficient_neighbors > 0) {
830 printf("less neighbors found using kth_neighbor_distance+eps than using knn search.\n");
831 all_passed = false;
832 }
833 double brute_force_search_time = timer.seconds();
834 printf("brute force search time: %f s\n", brute_force_search_time);
835 Kokkos::fence();
836 }
837 } else {
838 all_passed = false;
839 printf("do_radius_search and do_knn_search both false. Invalid combination.\n");
840 }
841}
842
843// finalize Kokkos and MPI (if available)
844Kokkos::finalize();
845#ifdef COMPADRE_USE_MPI
846MPI_Finalize();
847#endif
848
849// output to user that test passed or failed
850if(all_passed) {
851 fprintf(stdout, "Passed test \n");
852 return 0;
853} else {
854 fprintf(stdout, "Failed test \n");
855 return -1;
856}
857
858} // main
std::pair< T2, T1 > swap_pair(const std::pair< T1, T2 > &item)
std::multimap< T2, T1 > invert_map(const std::map< T1, T2 > &original_map)
int main(int argc, char *args[])
Manifold GMLS Example.
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....
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...
NeighborLists< view_type > CreateNeighborLists(view_type number_of_neighbors_list)
CreateNeighborLists allows for the construction of an object of type NeighborLists with template dedu...
NeighborLists< view_type_1d > Convert2DToCompressedRowNeighborLists(view_type_2d neighbor_lists)
Converts 2D neighbor lists to compressed row neighbor lists.