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