Compadre  1.6.0
GMLS_Multiple_Evaluation_Sites.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 
51 // code block to reduce scope for all Kokkos View allocations
52 // otherwise, Views may be deallocating when we call Kokkos::finalize() later
53 {
54 
55  CommandLineProcessor clp(argc, args);
56  auto order = clp.order;
57  auto dimension = clp.dimension;
58  auto number_target_coords = clp.number_target_coords;
59  auto constraint_name = clp.constraint_name;
60  auto solver_name = clp.solver_name;
61  auto problem_name = clp.problem_name;
62 
63  // the functions we will be seeking to reconstruct are in the span of the basis
64  // of the reconstruction space we choose for GMLS, so the error should be very small
65  const double failure_tolerance = 1e-9;
66 
67  // minimum neighbors for unisolvency is the same as the size of the polynomial basis
68  const int min_neighbors = Compadre::GMLS::getNP(order, dimension);
69 
70  //! [Parse Command Line Arguments]
71  Kokkos::Timer timer;
72  Kokkos::Profiling::pushRegion("Setup Point Data");
73  //! [Setting Up The Point Cloud]
74 
75  // approximate spacing of source sites
76  double h_spacing = 0.05;
77  int n_neg1_to_1 = 2*(1/h_spacing) + 1; // always odd
78 
79  // number of source coordinate sites that will fill a box of [-1,1]x[-1,1]x[-1,1] with a spacing approximately h
80  const int number_source_coords = std::pow(n_neg1_to_1, dimension);
81 
82  // coordinates of source sites
83  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device("source coordinates",
84  number_source_coords, 3);
85  Kokkos::View<double**>::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
86 
87  // coordinates of target sites
88  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device ("target coordinates", number_target_coords, 3);
89  Kokkos::View<double**>::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
90 
91  // coordinates of additional evaluation sites
92  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> additional_target_coords_device ("additional target coordinates", 2*number_target_coords /* multiple evaluation sites for each target index */, 3);
93  Kokkos::View<double**>::HostMirror additional_target_coords = Kokkos::create_mirror_view(additional_target_coords_device);
94 
95  // additional target site indices
96  Kokkos::View<int**, Kokkos::DefaultExecutionSpace> additional_target_indices_device ("additional target indices", number_target_coords, 4 /* # of extra evaluation sites plus index for each */);
97  Kokkos::View<int**>::HostMirror additional_target_indices = Kokkos::create_mirror_view(additional_target_indices_device);
98 
99 
100  // fill source coordinates with a uniform grid
101  int source_index = 0;
102  double this_coord[3] = {0,0,0};
103  for (int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
104  this_coord[0] = i*h_spacing;
105  for (int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
106  this_coord[1] = j*h_spacing;
107  for (int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
108  this_coord[2] = k*h_spacing;
109  if (dimension==3) {
110  source_coords(source_index,0) = this_coord[0];
111  source_coords(source_index,1) = this_coord[1];
112  source_coords(source_index,2) = this_coord[2];
113  source_index++;
114  }
115  }
116  if (dimension==2) {
117  source_coords(source_index,0) = this_coord[0];
118  source_coords(source_index,1) = this_coord[1];
119  source_coords(source_index,2) = 0;
120  source_index++;
121  }
122  }
123  if (dimension==1) {
124  source_coords(source_index,0) = this_coord[0];
125  source_coords(source_index,1) = 0;
126  source_coords(source_index,2) = 0;
127  source_index++;
128  }
129  }
130 
131  // fill target coords somewhere inside of [-0.5,0.5]x[-0.5,0.5]x[-0.5,0.5]
132  for(int i=0; i<number_target_coords; i++){
133 
134  // first, we get a uniformly random distributed direction
135  double rand_dir[3] = {0,0,0};
136 
137  for (int j=0; j<dimension; ++j) {
138  // rand_dir[j] is in [-0.5, 0.5]
139  rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
140  }
141 
142  // then we get a uniformly random radius
143  for (int j=0; j<dimension; ++j) {
144  target_coords(i,j) = rand_dir[j];
145  }
146 
147  }
148 
149  // generate coordinates to test multiple site evaluations
150  // strategy is to have a variable number of evaluation sites per target site
151  // so as to fully test the multi-site evaluation
152  int extra_evaluation_coordinates_count = 0;
153  for(int i=0; i<number_target_coords; i++){
154 
155  // set list of indices for extra evaluations
156  additional_target_indices(i,0) = (i%3)+1;
157 
158  // evaluation sites are same as target plus some perturbation
159  for (int k=0; k<(i%3+1); ++k) {
160  for (int j=0; j<dimension; ++j) {
161  additional_target_coords(extra_evaluation_coordinates_count,j) = target_coords(i,j) + (j==0)*1e-3 + (j==1)*1e-2 + (j==1)*(-1e-1);
162  }
163  additional_target_indices(i,k+1) = extra_evaluation_coordinates_count;
164  extra_evaluation_coordinates_count++;
165  }
166  }
167 
168 
169  //! [Setting Up The Point Cloud]
170 
171  Kokkos::Profiling::popRegion();
172  Kokkos::Profiling::pushRegion("Creating Data");
173 
174  //! [Creating The Data]
175 
176 
177  // source coordinates need copied to device before using to construct sampling data
178  Kokkos::deep_copy(source_coords_device, source_coords);
179 
180  // target coordinates copied next, because it is a convenient time to send them to device
181  Kokkos::deep_copy(target_coords_device, target_coords);
182 
183  // additional evaluation coordinates copied next, because it is a convenient time to send them to device
184  Kokkos::deep_copy(additional_target_coords_device, additional_target_coords);
185 
186  // additional evaluation indices copied next, because it is a convenient time to send them to device
187  Kokkos::deep_copy(additional_target_indices_device, additional_target_indices);
188 
189  // need Kokkos View storing true solution
190  Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device("samples of true solution",
191  source_coords_device.extent(0));
192 
193  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device("samples of true gradient",
194  source_coords_device.extent(0), dimension);
195 
196  Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
197  ("samples of true solution for divergence test", source_coords_device.extent(0), dimension);
198 
199  Kokkos::parallel_for("Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
200  (0,source_coords.extent(0)), KOKKOS_LAMBDA(const int i) {
201 
202  // coordinates of source site i
203  double xval = source_coords_device(i,0);
204  double yval = (dimension>1) ? source_coords_device(i,1) : 0;
205  double zval = (dimension>2) ? source_coords_device(i,2) : 0;
206 
207  // data for targets with scalar input
208  sampling_data_device(i) = trueSolution(xval, yval, zval, order, dimension);
209 
210  // data for targets with vector input (divergence)
211  double true_grad[3] = {0,0,0};
212  trueGradient(true_grad, xval, yval,zval, order, dimension);
213 
214  for (int j=0; j<dimension; ++j) {
215  gradient_sampling_data_device(i,j) = true_grad[j];
216 
217  // data for target with vector input (curl)
218  divergence_sampling_data_device(i,j) = divergenceTestSamples(xval, yval, zval, j, dimension);
219  }
220 
221  });
222 
223 
224  //! [Creating The Data]
225 
226  Kokkos::Profiling::popRegion();
227  Kokkos::Profiling::pushRegion("Neighbor Search");
228 
229  //! [Performing Neighbor Search]
230 
231 
232  // Point cloud construction for neighbor search
233  // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
234  auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
235 
236  // each row is a neighbor list for a target site, with the first column of each row containing
237  // the number of neighbors for that rows corresponding target site
238  double epsilon_multiplier = 1.5;
239  int estimated_upper_bound_number_neighbors =
240  point_cloud_search.getEstimatedNumberNeighborsUpperBound(min_neighbors, dimension, epsilon_multiplier);
241 
242  Kokkos::View<int**, Kokkos::DefaultExecutionSpace> neighbor_lists_device("neighbor lists",
243  number_target_coords, estimated_upper_bound_number_neighbors); // first column is # of neighbors
244  Kokkos::View<int**>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
245 
246  // each target site has a window size
247  Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device("h supports", number_target_coords);
248  Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
249 
250  // query the point cloud to generate the neighbor lists using a kdtree to produce the n nearest neighbor
251  // to each target site, adding (epsilon_multiplier-1)*100% to whatever the distance away the further neighbor used is from
252  // each target to the view for epsilon
253  point_cloud_search.generate2DNeighborListsFromKNNSearch(false /*not dry run*/, target_coords, neighbor_lists,
254  epsilon, min_neighbors, epsilon_multiplier);
255 
256 
257  //! [Performing Neighbor Search]
258 
259  Kokkos::Profiling::popRegion();
260  Kokkos::fence(); // let call to build neighbor lists complete before copying back to device
261  timer.reset();
262 
263  //! [Setting Up The GMLS Object]
264 
265 
266  // Copy data back to device (they were filled on the host)
267  // We could have filled Kokkos Views with memory space on the host
268  // and used these instead, and then the copying of data to the device
269  // would be performed in the GMLS class
270  Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
271  Kokkos::deep_copy(epsilon_device, epsilon);
272 
273  // initialize an instance of the GMLS class
275  order, dimension,
276  solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
277  2 /*manifold order*/);
278 
279  // pass in neighbor lists, source coordinates, target coordinates, and window sizes
280  //
281  // neighbor lists have the format:
282  // dimensions: (# number of target sites) X (# maximum number of neighbors for any given target + 1)
283  // the first column contains the number of neighbors for that rows corresponding target index
284  //
285  // source coordinates have the format:
286  // dimensions: (# number of source sites) X (dimension)
287  // entries in the neighbor lists (integers) correspond to rows of this 2D array
288  //
289  // target coordinates have the format:
290  // dimensions: (# number of target sites) X (dimension)
291  // # of target sites is same as # of rows of neighbor lists
292  //
293  my_GMLS.setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
294 
295  // set up additional sites to evaluate target operators
296  my_GMLS.setAdditionalEvaluationSitesData(additional_target_indices_device, additional_target_coords_device);
297 
298  // create a vector of target operations
299  std::vector<TargetOperation> lro(2);
300  lro[0] = ScalarPointEvaluation;
302 
303  // and then pass them to the GMLS class
304  my_GMLS.addTargets(lro);
305 
306  // sets the weighting kernel function from WeightingFunctionType
308 
309  // power to use in that weighting kernel function
310  my_GMLS.setWeightingParameter(2);
311 
312  // generate the alphas that to be combined with data for each target operation requested in lro
313  my_GMLS.generateAlphas();
314 
315 
316  //! [Setting Up The GMLS Object]
317 
318  double instantiation_time = timer.seconds();
319  std::cout << "Took " << instantiation_time << "s to complete alphas generation." << std::endl;
320  Kokkos::fence(); // let generateAlphas finish up before using alphas
321  Kokkos::Profiling::pushRegion("Apply Alphas to Data");
322 
323  //! [Apply GMLS Alphas To Data]
324 
325  // it is important to note that if you expect to use the data as a 1D view, then you should use double*
326  // however, if you know that the target operation will result in a 2D view (vector or matrix output),
327  // then you should template with double** as this is something that can not be infered from the input data
328  // or the target operator at compile time. Additionally, a template argument is required indicating either
329  // Kokkos::HostSpace or Kokkos::DefaultExecutionSpace::memory_space()
330 
331  // The Evaluator class takes care of handling input data views as well as the output data views.
332  // It uses information from the GMLS class to determine how many components are in the input
333  // as well as output for any choice of target functionals and then performs the contactions
334  // on the data using the alpha coefficients generated by the GMLS class, all on the device.
335  Evaluator gmls_evaluator(&my_GMLS);
336 
337  auto output_value1 = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
338  (sampling_data_device, ScalarPointEvaluation, PointSample,
339  true /*scalar_as_vector_if_needed*/, 1 /*evaluation site index*/);
340 
341  auto output_gradient1 = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
342  (sampling_data_device, GradientOfScalarPointEvaluation, PointSample,
343  true /*scalar_as_vector_if_needed*/, 1 /*evaluation site index*/);
344 
345  auto output_value2 = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
346  (sampling_data_device, ScalarPointEvaluation, PointSample,
347  true /*scalar_as_vector_if_needed*/, 2 /*evaluation site index*/);
348 
349  auto output_gradient2 = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
350  (sampling_data_device, GradientOfScalarPointEvaluation, PointSample,
351  true /*scalar_as_vector_if_needed*/, 2 /*evaluation site index*/);
352 
353  auto output_value3 = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
354  (sampling_data_device, ScalarPointEvaluation, PointSample,
355  true /*scalar_as_vector_if_needed*/, 3 /*evaluation site index*/);
356 
357  auto output_gradient3 = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
358  (sampling_data_device, GradientOfScalarPointEvaluation, PointSample,
359  true /*scalar_as_vector_if_needed*/, 3 /*evaluation site index*/);
360 
361  //! [Apply GMLS Alphas To Data]
362 
363  Kokkos::fence(); // let application of alphas to data finish before using results
364  Kokkos::Profiling::popRegion();
365  // times the Comparison in Kokkos
366  Kokkos::Profiling::pushRegion("Comparison");
367 
368  //! [Check That Solutions Are Correct]
369 
370  // load value from output
371  double GMLS_value;
372  // load partial x from gradient
373  double GMLS_GradX;
374  // load partial y from gradient
375  double GMLS_GradY;
376  // load partial z from gradient
377  double GMLS_GradZ;
378 
379  // loop through the target sites
380  extra_evaluation_coordinates_count = 0;
381  for (int i=0; i<number_target_coords; i++) {
382 
383  for (int k=0; k<(i%3)+1; ++k) {
384  if (k==0) {
385  // load value from output
386  GMLS_value = output_value1(i);
387  // load partial x from gradient
388  GMLS_GradX = output_gradient1(i,0);
389  // load partial y from gradient
390  GMLS_GradY = (dimension>1) ? output_gradient1(i,1) : 0;
391  // load partial z from gradient
392  GMLS_GradZ = (dimension>2) ? output_gradient1(i,2) : 0;
393  } else if (k==1) {
394  // load value from output
395  GMLS_value = output_value2(i);
396  // load partial x from gradient
397  GMLS_GradX = output_gradient2(i,0);
398  // load partial y from gradient
399  GMLS_GradY = (dimension>1) ? output_gradient2(i,1) : 0;
400  // load partial z from gradient
401  GMLS_GradZ = (dimension>2) ? output_gradient2(i,2) : 0;
402  } else if (k==2) {
403  // load value from output
404  GMLS_value = output_value3(i);
405  // load partial x from gradient
406  GMLS_GradX = output_gradient3(i,0);
407  // load partial y from gradient
408  GMLS_GradY = (dimension>1) ? output_gradient3(i,1) : 0;
409  // load partial z from gradient
410  GMLS_GradZ = (dimension>2) ? output_gradient3(i,2) : 0;
411  }
412 
413  // target site i's coordinate
414  double xval = additional_target_coords(extra_evaluation_coordinates_count,0);
415  double yval = additional_target_coords(extra_evaluation_coordinates_count,1);
416  double zval = additional_target_coords(extra_evaluation_coordinates_count,2);
417 
418  // evaluation of various exact solutions
419  double actual_value = trueSolution(xval, yval, zval, order, dimension);
420 
421  double actual_Gradient[3] = {0,0,0}; // initialized for 3, but only filled up to dimension
422  trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
423 
424  // check actual function value
425  if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
426  all_passed = false;
427  std::cout << i << " Failed Actual by: " << std::abs(actual_value - GMLS_value) << " for evaluation site: " << k << std::endl;
428  }
429 
430  // check gradient
431  if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
432  all_passed = false;
433  std::cout << i << " Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << " for evaluation site: " << k << std::endl;
434  if (dimension>1) {
435  if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
436  all_passed = false;
437  std::cout << i << " Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << " for evaluation site: " << k << std::endl;
438  }
439  }
440  if (dimension>2) {
441  if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
442  all_passed = false;
443  std::cout << i << " Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << " for evaluation site: " << k << std::endl;
444  }
445  }
446  }
447  extra_evaluation_coordinates_count++;
448  }
449  }
450 
451 
452  //! [Check That Solutions Are Correct]
453  // popRegion hidden from tutorial
454  // stop timing comparison loop
455  Kokkos::Profiling::popRegion();
456  //! [Finalize Program]
457 
458 
459 } // end of code block to reduce scope, causing Kokkos View de-allocations
460 // otherwise, Views may be deallocating when we call Kokkos::finalize() later
461 
462 // finalize Kokkos and MPI (if available)
463 Kokkos::finalize();
464 #ifdef COMPADRE_USE_MPI
465 MPI_Finalize();
466 #endif
467 
468 // output to user that test passed or failed
469 if(all_passed) {
470  fprintf(stdout, "Passed test \n");
471  return 0;
472 } else {
473  fprintf(stdout, "Failed test \n");
474  return -1;
475 }
476 
477 } // main
478 
479 
480 //! [Finalize Program]
int main(int argc, char *args[])
[Parse Command Line Arguments]
KOKKOS_INLINE_FUNCTION double trueSolution(double x, double y, double z, int order, int dimension)
KOKKOS_INLINE_FUNCTION void trueGradient(double *ans, double x, double y, double z, int order, int dimension)
KOKKOS_INLINE_FUNCTION double divergenceTestSamples(double x, double y, double z, int component, int dimension)
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
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)
Generalized Moving Least Squares (GMLS)
void addTargets(TargetOperation lro)
Adds a target to the vector of target functional to be applied to the reconstruction.
void setAdditionalEvaluationSitesData(view_type_1 additional_evaluation_indices, view_type_2 additional_evaluation_coordinates)
(OPTIONAL) Sets additional evaluation sites for each target site
void setWeightingParameter(int wp, int index=0)
Parameter for weighting kernel for GMLS problem index = 0 sets p paramater for weighting kernel index...
void generateAlphas(const int number_of_batches=1, const bool keep_coefficients=false, const bool clear_cache=true)
Meant to calculate target operations and apply the evaluations to the previously constructed polynomi...
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)
void setWeightingType(const std::string &wt)
Type for weighting kernel for GMLS problem.
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....
constexpr SamplingFunctional PointSample
Available sampling functionals.
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...
@ VectorOfScalarClonesTaylorPolynomial
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.
@ GradientOfScalarPointEvaluation
Point evaluation of the gradient of a scalar.
@ ScalarPointEvaluation
Point evaluation of a scalar.