Compadre  1.5.5
GMLS_Device.cpp
Go to the documentation of this file.
1 #include <iostream>
2 #include <string>
3 #include <vector>
4 #include <map>
5 #include <stdlib.h>
6 #include <cstdio>
7 #include <random>
8 
9 #include <Compadre_Config.h>
10 #include <Compadre_GMLS.hpp>
11 #include <Compadre_Evaluator.hpp>
13 
14 #include "GMLS_Tutorial.hpp"
15 #include "CommandLineProcessor.hpp"
16 
17 #ifdef COMPADRE_USE_MPI
18 #include <mpi.h>
19 #endif
20 
21 #include <Kokkos_Timer.hpp>
22 #include <Kokkos_Core.hpp>
23 
24 using namespace Compadre;
25 
26 //! [Parse Command Line Arguments]
27 
28 // called from command line
29 int main (int argc, char* args[]) {
30 
31 // initializes MPI (if available) with command line arguments given
32 #ifdef COMPADRE_USE_MPI
33 MPI_Init(&argc, &args);
34 #endif
35 
36 // initializes Kokkos with command line arguments given
37 Kokkos::initialize(argc, args);
38 
39 // becomes false if the computed solution not within the failure_threshold of the actual solution
40 bool all_passed = true;
41 
42 // code block to reduce scope for all Kokkos View allocations
43 // otherwise, Views may be deallocating when we call Kokkos finalize() later
44 {
45 
46  CommandLineProcessor clp(argc, args);
47  auto order = clp.order;
48  auto dimension = clp.dimension;
49  auto number_target_coords = clp.number_target_coords;
50  auto constraint_name = clp.constraint_name;
51  auto solver_name = clp.solver_name;
52  auto problem_name = clp.problem_name;
53  auto number_of_batches = clp.number_of_batches;
54  bool keep_coefficients = number_of_batches==1;
55 
56  // the functions we will be seeking to reconstruct are in the span of the basis
57  // of the reconstruction space we choose for GMLS, so the error should be very small
58  const double failure_tolerance = 1e-9;
59 
60  // Laplacian is a second order differential operator, which we expect to be slightly less accurate
61  const double laplacian_failure_tolerance = 1e-9;
62 
63  // minimum neighbors for unisolvency is the same as the size of the polynomial basis
64  const int min_neighbors = Compadre::GMLS::getNP(order, dimension);
65 
66  //! [Parse Command Line Arguments]
67  Kokkos::Timer timer;
68  Kokkos::Profiling::pushRegion("Setup Point Data");
69  //! [Setting Up The Point Cloud]
70 
71  // approximate spacing of source sites
72  double h_spacing = 0.05;
73  int n_neg1_to_1 = 2*(1/h_spacing) + 1; // always odd
74 
75  // number of source coordinate sites that will fill a box of [-1,1]x[-1,1]x[-1,1] with a spacing approximately h
76  const int number_source_coords = std::pow(n_neg1_to_1, dimension);
77 
78  // coordinates of source sites
79  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device("source coordinates",
80  number_source_coords, 3);
81  Kokkos::View<double**>::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
82 
83  // coordinates of target sites
84  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device ("target coordinates", number_target_coords, 3);
85  Kokkos::View<double**>::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
86 
87 
88  // fill source coordinates with a uniform grid
89  int source_index = 0;
90  double this_coord[3] = {0,0,0};
91  for (int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
92  this_coord[0] = i*h_spacing;
93  for (int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
94  this_coord[1] = j*h_spacing;
95  for (int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
96  this_coord[2] = k*h_spacing;
97  if (dimension==3) {
98  source_coords(source_index,0) = this_coord[0];
99  source_coords(source_index,1) = this_coord[1];
100  source_coords(source_index,2) = this_coord[2];
101  source_index++;
102  }
103  }
104  if (dimension==2) {
105  source_coords(source_index,0) = this_coord[0];
106  source_coords(source_index,1) = this_coord[1];
107  source_coords(source_index,2) = 0;
108  source_index++;
109  }
110  }
111  if (dimension==1) {
112  source_coords(source_index,0) = this_coord[0];
113  source_coords(source_index,1) = 0;
114  source_coords(source_index,2) = 0;
115  source_index++;
116  }
117  }
118 
119  // fill target coords somewhere inside of [-0.5,0.5]x[-0.5,0.5]x[-0.5,0.5]
120  for(int i=0; i<number_target_coords; i++){
121 
122  // first, we get a uniformly random distributed direction
123  double rand_dir[3] = {0,0,0};
124 
125  for (int j=0; j<dimension; ++j) {
126  // rand_dir[j] is in [-0.5, 0.5]
127  rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
128  }
129 
130  // then we get a uniformly random radius
131  for (int j=0; j<dimension; ++j) {
132  target_coords(i,j) = rand_dir[j];
133  }
134 
135  }
136 
137 
138  //! [Setting Up The Point Cloud]
139 
140  Kokkos::Profiling::popRegion();
141  Kokkos::Profiling::pushRegion("Creating Data");
142 
143  //! [Creating The Data]
144 
145 
146  // source coordinates need copied to device before using to construct sampling data
147  Kokkos::deep_copy(source_coords_device, source_coords);
148 
149  // target coordinates copied next, because it is a convenient time to send them to device
150  Kokkos::deep_copy(target_coords_device, target_coords);
151 
152  // need Kokkos View storing true solution
153  Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device("samples of true solution",
154  source_coords_device.extent(0));
155 
156  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device("samples of true gradient",
157  source_coords_device.extent(0), dimension);
158 
159  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
160  ("samples of true solution for divergence test", source_coords_device.extent(0), dimension);
161 
162  Kokkos::parallel_for("Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
163  (0,source_coords.extent(0)), KOKKOS_LAMBDA(const int i) {
164 
165  // coordinates of source site i
166  double xval = source_coords_device(i,0);
167  double yval = (dimension>1) ? source_coords_device(i,1) : 0;
168  double zval = (dimension>2) ? source_coords_device(i,2) : 0;
169 
170  // data for targets with scalar input
171  sampling_data_device(i) = trueSolution(xval, yval, zval, order, dimension);
172 
173  // data for targets with vector input (divergence)
174  double true_grad[3] = {0,0,0};
175  trueGradient(true_grad, xval, yval,zval, order, dimension);
176 
177  for (int j=0; j<dimension; ++j) {
178  gradient_sampling_data_device(i,j) = true_grad[j];
179 
180  // data for target with vector input (curl)
181  divergence_sampling_data_device(i,j) = divergenceTestSamples(xval, yval, zval, j, dimension);
182  }
183 
184  });
185 
186 
187  //! [Creating The Data]
188 
189  Kokkos::Profiling::popRegion();
190  Kokkos::Profiling::pushRegion("Neighbor Search");
191 
192  //! [Performing Neighbor Search]
193 
194 
195  // Point cloud construction for neighbor search
196  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
197  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
198 
199  double epsilon_multiplier = 1.4;
200 
201  // neighbor_lists_device will contain all neighbor lists (for each target site) in a compressed row format
202  // Initially, we do a dry-run to calculate neighborhood sizes before actually storing the result. This is
203  // why we can start with a neighbor_lists_device size of 0.
204  Kokkos::View<int*> neighbor_lists_device("neighbor lists",
205  0); // first column is # of neighbors
206  Kokkos::View<int*>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
207  // number_of_neighbors_list must be the same size as the number of target sites so that it can be populated
208  // with the number of neighbors for each target site.
209  Kokkos::View<int*> number_of_neighbors_list_device("number of neighbor lists",
210  number_target_coords); // first column is # of neighbors
211  Kokkos::View<int*>::HostMirror number_of_neighbors_list = Kokkos::create_mirror_view(number_of_neighbors_list_device);
212 
213  // each target site has a window size
214  Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device("h supports", number_target_coords);
215  Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
216 
217  // query the point cloud to generate the neighbor lists using a kdtree to produce the n nearest neighbor
218  // to each target site, adding (epsilon_multiplier-1)*100% to whatever the distance away the further neighbor used is from
219  // each target to the view for epsilon
220  //
221  // This dry run populates number_of_neighbors_list with neighborhood sizes
222  size_t storage_size = point_cloud_search.generateCRNeighborListsFromKNNSearch(true /*dry run*/, target_coords, neighbor_lists,
223  number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
224 
225  // resize neighbor_lists_device so as to be large enough to contain all neighborhoods
226  Kokkos::resize(neighbor_lists_device, storage_size);
227  neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
228 
229  // query the point cloud a second time, but this time storing results into neighbor_lists
230  point_cloud_search.generateCRNeighborListsFromKNNSearch(false /*not dry run*/, target_coords, neighbor_lists,
231  number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
232 
233  //! [Performing Neighbor Search]
234 
235  Kokkos::Profiling::popRegion();
236  Kokkos::fence(); // let call to build neighbor lists complete before copying back to device
237  timer.reset();
238 
239  //! [Setting Up The GMLS Object]
240 
241 
242  // Copy data back to device (they were filled on the host)
243  // We could have filled Kokkos Views with memory space on the host
244  // and used these instead, and then the copying of data to the device
245  // would be performed in the GMLS class
246  Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
247  Kokkos::deep_copy(number_of_neighbors_list_device, number_of_neighbors_list);
248  Kokkos::deep_copy(epsilon_device, epsilon);
249 
250  // initialize an instance of the GMLS class
252  order, dimension,
253  solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
254  2 /*manifold order*/);
255 
256  // pass in neighbor lists, number of neighbor lists, source coordinates, target coordinates, and window sizes
257  //
258  // neighbor lists has a compressed row format and is a 1D view:
259  // dimensions: ? (determined by neighbor search, but total size of the sum of the number of neighbors over all target sites)
260  //
261  // number of neighbors list is a 1D view:
262  // dimensions: (# number of target sites)
263  // each entry contains the number of neighbors for a target site
264  //
265  // source coordinates have the format:
266  // dimensions: (# number of source sites) X (dimension)
267  // entries in the neighbor lists (integers) correspond to rows of this 2D array
268  //
269  // target coordinates have the format:
270  // dimensions: (# number of target sites) X (dimension)
271  // # of target sites is same as # of rows of neighbor lists
272  //
273  my_GMLS.setProblemData(neighbor_lists_device, number_of_neighbors_list_device, source_coords_device, target_coords_device, epsilon_device);
274 
275  // create a vector of target operations
276  std::vector<TargetOperation> lro(5);
277  lro[0] = ScalarPointEvaluation;
282 
283  // and then pass them to the GMLS class
284  my_GMLS.addTargets(lro);
285 
286  // sets the weighting kernel function from WeightingFunctionType
287  my_GMLS.setWeightingType(WeightingFunctionType::Power);
288 
289  // power to use in that weighting kernel function
290  my_GMLS.setWeightingParameter(2);
291 
292  // generate the alphas that to be combined with data for each target operation requested in lro
293  my_GMLS.generateAlphas(number_of_batches, keep_coefficients /* keep polynomial coefficients, only needed for a test later in this program */);
294 
295 
296  //! [Setting Up The GMLS Object]
297 
298  double instantiation_time = timer.seconds();
299  std::cout << "Took " << instantiation_time << "s to complete alphas generation." << std::endl;
300  Kokkos::fence(); // let generateAlphas finish up before using alphas
301  Kokkos::Profiling::pushRegion("Apply Alphas to Data");
302 
303  //! [Apply GMLS Alphas To Data]
304 
305  // it is important to note that if you expect to use the data as a 1D view, then you should use double*
306  // however, if you know that the target operation will result in a 2D view (vector or matrix output),
307  // then you should template with double** as this is something that can not be infered from the input data
308  // or the target operator at compile time. Additionally, a template argument is required indicating either
309  // Kokkos::HostSpace or Kokkos::DefaultExecutionSpace::memory_space()
310 
311  // The Evaluator class takes care of handling input data views as well as the output data views.
312  // It uses information from the GMLS class to determine how many components are in the input
313  // as well as output for any choice of target functionals and then performs the contactions
314  // on the data using the alpha coefficients generated by the GMLS class, all on the device.
315  Evaluator gmls_evaluator(&my_GMLS);
316 
317  auto output_value = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
318  (sampling_data_device, ScalarPointEvaluation);
319 
320  auto output_laplacian = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
321  (sampling_data_device, LaplacianOfScalarPointEvaluation);
322 
323  auto output_gradient = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
324  (sampling_data_device, GradientOfScalarPointEvaluation);
325 
326  auto output_divergence = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
327  (gradient_sampling_data_device, DivergenceOfVectorPointEvaluation, VectorPointSample);
328 
329  auto output_curl = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
330  (divergence_sampling_data_device, CurlOfVectorPointEvaluation);
331 
332  // retrieves polynomial coefficients instead of remapped field
333  decltype(output_curl) scalar_coefficients;
334  if (number_of_batches==1)
335  scalar_coefficients =
336  gmls_evaluator.applyFullPolynomialCoefficientsBasisToDataAllComponents<double**, Kokkos::HostSpace>
337  (sampling_data_device);
338 
339  //! [Apply GMLS Alphas To Data]
340 
341  Kokkos::fence(); // let application of alphas to data finish before using results
342  Kokkos::Profiling::popRegion();
343  // times the Comparison in Kokkos
344  Kokkos::Profiling::pushRegion("Comparison");
345 
346  //! [Check That Solutions Are Correct]
347 
348 
349  // loop through the target sites
350  for (int i=0; i<number_target_coords; i++) {
351 
352  // load value from output
353  double GMLS_value = output_value(i);
354 
355  // load laplacian from output
356  double GMLS_Laplacian = output_laplacian(i);
357 
358  // load partial x from gradient
359  // this is a test that the scalar_coefficients 2d array returned hold valid entries
360  // scalar_coefficients(i,1)*1./epsilon(i) is equivalent to the target operation acting
361  // on the polynomials applied to the polynomial coefficients
362  double GMLS_GradX = (number_of_batches==1) ? scalar_coefficients(i,1)*1./epsilon(i) : output_gradient(i,0);
363 
364  // load partial y from gradient
365  double GMLS_GradY = (dimension>1) ? output_gradient(i,1) : 0;
366 
367  // load partial z from gradient
368  double GMLS_GradZ = (dimension>2) ? output_gradient(i,2) : 0;
369 
370  // load divergence from output
371  double GMLS_Divergence = output_divergence(i);
372 
373  // load curl from output
374  double GMLS_CurlX = (dimension>1) ? output_curl(i,0) : 0;
375  double GMLS_CurlY = (dimension>1) ? output_curl(i,1) : 0;
376  double GMLS_CurlZ = (dimension>2) ? output_curl(i,2) : 0;
377 
378 
379  // target site i's coordinate
380  double xval = target_coords(i,0);
381  double yval = (dimension>1) ? target_coords(i,1) : 0;
382  double zval = (dimension>2) ? target_coords(i,2) : 0;
383 
384  // evaluation of various exact solutions
385  double actual_value = trueSolution(xval, yval, zval, order, dimension);
386  double actual_Laplacian = trueLaplacian(xval, yval, zval, order, dimension);
387 
388  double actual_Gradient[3] = {0,0,0}; // initialized for 3, but only filled up to dimension
389  trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
390 
391  double actual_Divergence;
392  actual_Divergence = trueLaplacian(xval, yval, zval, order, dimension);
393 
394  double actual_Curl[3] = {0,0,0}; // initialized for 3, but only filled up to dimension
395  // (and not at all for dimimension = 1)
396  if (dimension>1) {
397  actual_Curl[0] = curlTestSolution(xval, yval, zval, 0, dimension);
398  actual_Curl[1] = curlTestSolution(xval, yval, zval, 1, dimension);
399  if (dimension>2) {
400  actual_Curl[2] = curlTestSolution(xval, yval, zval, 2, dimension);
401  }
402  }
403 
404  // check actual function value
405  if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
406  all_passed = false;
407  std::cout << i << " Failed Actual by: " << std::abs(actual_value - GMLS_value) << std::endl;
408  }
409 
410  // check Laplacian
411  if(std::abs(actual_Laplacian - GMLS_Laplacian) > laplacian_failure_tolerance) {
412  all_passed = false;
413  std::cout << i <<" Failed Laplacian by: " << std::abs(actual_Laplacian - GMLS_Laplacian) << std::endl;
414  }
415 
416  // check gradient
417  if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
418  all_passed = false;
419  std::cout << i << " Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << std::endl;
420  if (dimension>1) {
421  if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
422  all_passed = false;
423  std::cout << i << " Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << std::endl;
424  }
425  }
426  if (dimension>2) {
427  if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
428  all_passed = false;
429  std::cout << i << " Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << std::endl;
430  }
431  }
432  }
433 
434  // check divergence
435  if(std::abs(actual_Divergence - GMLS_Divergence) > failure_tolerance) {
436  all_passed = false;
437  std::cout << i << " Failed Divergence by: " << std::abs(actual_Divergence - GMLS_Divergence) << std::endl;
438  }
439 
440  // check curl
441  if (order > 2) { // reconstructed solution not in basis unless order greater than 2 used
442  double tmp_diff = 0;
443  if (dimension>1)
444  tmp_diff += std::abs(actual_Curl[0] - GMLS_CurlX) + std::abs(actual_Curl[1] - GMLS_CurlY);
445  if (dimension>2)
446  tmp_diff += std::abs(actual_Curl[2] - GMLS_CurlZ);
447  if(std::abs(tmp_diff) > failure_tolerance) {
448  all_passed = false;
449  std::cout << i << " Failed Curl by: " << std::abs(tmp_diff) << std::endl;
450  }
451  }
452  }
453 
454 
455  //! [Check That Solutions Are Correct]
456  // popRegion hidden from tutorial
457  // stop timing comparison loop
458  Kokkos::Profiling::popRegion();
459  //! [Finalize Program]
460 
461 
462 } // end of code block to reduce scope, causing Kokkos View de-allocations
463 // otherwise, Views may be deallocating when we call Kokkos finalize() later
464 
465 // finalize Kokkos and MPI (if available)
466 Kokkos::finalize();
467 #ifdef COMPADRE_USE_MPI
468 MPI_Finalize();
469 #endif
470 
471 // output to user that test passed or failed
472 if(all_passed) {
473  fprintf(stdout, "Passed test \n");
474  return 0;
475 } else {
476  fprintf(stdout, "Failed test \n");
477  return -1;
478 }
479 
480 } // main
481 
482 
483 //! [Finalize Program]
Point evaluation of the curl of a vector (results in a vector)
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
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...
int main(int argc, char *args[])
[Parse Command Line Arguments]
Definition: GMLS_Device.cpp:29
Point evaluation of a scalar.
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...
KOKKOS_INLINE_FUNCTION void trueGradient(double *ans, double x, double y, double z, int order, int dimension)
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
Kokkos::View< output_data_type, output_array_layout, output_memory_space > applyAlphasToDataAllComponentsAllTargetSites(view_type_input_data sampling_data, TargetOperation lro, const SamplingFunctional sro_in=PointSample, bool scalar_as_vector_if_needed=true, const int evaluation_site_local_index=0) const
Transformation of data under GMLS (allocates memory for output)
Kokkos::View< output_data_type, output_array_layout, output_memory_space > applyFullPolynomialCoefficientsBasisToDataAllComponents(view_type_input_data sampling_data, bool scalar_as_vector_if_needed=true) const
Generation of polynomial reconstruction coefficients by applying to data in GMLS (allocates memory fo...
Point evaluation of the laplacian of a scalar (could be on a manifold or not)
Point evaluation of the divergence of a vector (results in a scalar)
Point evaluation of the gradient of a scalar.
KOKKOS_INLINE_FUNCTION double trueLaplacian(double x, double y, double z, int order, int dimension)
Generalized Moving Least Squares (GMLS)
KOKKOS_INLINE_FUNCTION double curlTestSolution(double x, double y, double z, int component, int dimension)
KOKKOS_INLINE_FUNCTION double divergenceTestSamples(double x, double y, double z, int component, int dimension)
void setProblemData(view_type_1 neighbor_lists, view_type_2 source_coordinates, view_type_3 target_coordinates, view_type_4 epsilons)
Sets basic problem data (neighbor lists, source coordinates, and target coordinates) ...
KOKKOS_INLINE_FUNCTION double trueSolution(double x, double y, double z, int order, int dimension)
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.