Compadre  1.6.0
GMLS_NeumannGradScalar.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 
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  // tangent bundle for each target sites
92  Kokkos::View<double***, Kokkos::DefaultExecutionSpace> tangent_bundles_device ("tangent bundles", number_target_coords, dimension, dimension);
93  Kokkos::View<double***>::HostMirror tangent_bundles = Kokkos::create_mirror_view(tangent_bundles_device);
94 
95  // fill source coordinates with a uniform grid
96  int source_index = 0;
97  double this_coord[3] = {0,0,0};
98  for (int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
99  this_coord[0] = i*h_spacing;
100  for (int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
101  this_coord[1] = j*h_spacing;
102  for (int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
103  this_coord[2] = k*h_spacing;
104  if (dimension==3) {
105  source_coords(source_index,0) = this_coord[0];
106  source_coords(source_index,1) = this_coord[1];
107  source_coords(source_index,2) = this_coord[2];
108  source_index++;
109  }
110  }
111  if (dimension==2) {
112  source_coords(source_index,0) = this_coord[0];
113  source_coords(source_index,1) = this_coord[1];
114  source_coords(source_index,2) = 0;
115  source_index++;
116  }
117  }
118  if (dimension==1) {
119  source_coords(source_index,0) = this_coord[0];
120  source_coords(source_index,1) = 0;
121  source_coords(source_index,2) = 0;
122  source_index++;
123  }
124  }
125 
126  // fill target coords somewhere inside of [-0.5,0.5]x[-0.5,0.5]x[-0.5,0.5]
127  for(int i=0; i<number_target_coords; i++){
128 
129  // first, we get a uniformly random distributed direction
130  double rand_dir[3] = {0,0,0};
131 
132  for (int j=0; j<dimension; ++j) {
133  // rand_dir[j] is in [-0.5, 0.5]
134  rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
135  }
136 
137  // then we get a uniformly random radius
138  for (int j=0; j<dimension; ++j) {
139  target_coords(i,j) = rand_dir[j];
140  }
141  // target_coords(i, 2) = 1.0;
142 
143  // Set tangent bundles
144  if (dimension == 2) {
145  tangent_bundles(i, 0, 0) = 0.0;
146  tangent_bundles(i, 0, 1) = 0.0;
147  tangent_bundles(i, 1, 0) = 1.0/(sqrt(2.0));
148  tangent_bundles(i, 1, 1) = 1.0/(sqrt(2.0));
149  } else if (dimension == 3) {
150  tangent_bundles(i, 0, 0) = 0.0;
151  tangent_bundles(i, 0, 1) = 0.0;
152  tangent_bundles(i, 0, 2) = 0.0;
153  tangent_bundles(i, 1, 0) = 0.0;
154  tangent_bundles(i, 1, 1) = 0.0;
155  tangent_bundles(i, 1, 2) = 0.0;
156  tangent_bundles(i, 2, 0) = 1.0/(sqrt(3.0));
157  tangent_bundles(i, 2, 1) = 1.0/(sqrt(3.0));
158  tangent_bundles(i, 2, 2) = 1.0/(sqrt(3.0));
159  }
160  }
161 
162  //! [Setting Up The Point Cloud]
163 
164  Kokkos::Profiling::popRegion();
165  Kokkos::Profiling::pushRegion("Creating Data");
166 
167  //! [Creating The Data]
168 
169 
170  // source coordinates need copied to device before using to construct sampling data
171  Kokkos::deep_copy(source_coords_device, source_coords);
172 
173  // target coordinates copied next, because it is a convenient time to send them to device
174  Kokkos::deep_copy(target_coords_device, target_coords);
175 
176  // tangent bundles copied next, because it is a convenient time to send them to device
177  Kokkos::deep_copy(tangent_bundles_device, tangent_bundles);
178 
179  // need Kokkos View storing true solution
180  Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device("samples of true solution",
181  source_coords_device.extent(0));
182 
183  Kokkos::parallel_for("Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
184  (0,source_coords.extent(0)), KOKKOS_LAMBDA(const int i) {
185 
186  // coordinates of source site i
187  double xval = source_coords_device(i,0);
188  double yval = (dimension>1) ? source_coords_device(i,1) : 0;
189  double zval = (dimension>2) ? source_coords_device(i,2) : 0;
190 
191  // data for targets with scalar input
192  sampling_data_device(i) = trueSolution(xval, yval, zval, order, dimension);
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  // each row is a neighbor list for a target site, with the first column of each row containing
208  // the number of neighbors for that rows corresponding target site
209  double epsilon_multiplier = 1.8;
210  int estimated_upper_bound_number_neighbors =
211  point_cloud_search.getEstimatedNumberNeighborsUpperBound(min_neighbors, dimension, epsilon_multiplier);
212 
213  Kokkos::View<int**, Kokkos::DefaultExecutionSpace> neighbor_lists_device("neighbor lists",
214  number_target_coords, estimated_upper_bound_number_neighbors); // first column is # of neighbors
215  Kokkos::View<int**>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
216 
217  // each target site has a window size
218  Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device("h supports", number_target_coords);
219  Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
220 
221  // query the point cloud to generate the neighbor lists using a kdtree to produce the n nearest neighbor
222  // to each target site, adding (epsilon_multiplier-1)*100% to whatever the distance away the further neighbor used is from
223  // each target to the view for epsilon
224  point_cloud_search.generate2DNeighborListsFromKNNSearch(false /*not dry run*/, target_coords, neighbor_lists,
225  epsilon, min_neighbors, epsilon_multiplier);
226 
227  //! [Performing Neighbor Search]
228 
229  Kokkos::Profiling::popRegion();
230  Kokkos::fence(); // let call to build neighbor lists complete before copying back to device
231  timer.reset();
232 
233  //! [Setting Up The GMLS Object]
234 
235 
236  // Copy data back to device (they were filled on the host)
237  // We could have filled Kokkos Views with memory space on the host
238  // and used these instead, and then the copying of data to the device
239  // would be performed in the GMLS class
240  Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
241  Kokkos::deep_copy(epsilon_device, epsilon);
242 
243  // initialize an instance of the GMLS class
245  PointSample,
246  order, dimension,
247  solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
248  0 /*manifold order*/);
249 
250  // pass in neighbor lists, source coordinates, target coordinates, and window sizes
251  //
252  // neighbor lists have the format:
253  // dimensions: (# number of target sites) X (# maximum number of neighbors for any given target + 1)
254  // the first column contains the number of neighbors for that rows corresponding target index
255  //
256  // source coordinates have the format:
257  // dimensions: (# number of source sites) X (dimension)
258  // entries in the neighbor lists (integers) correspond to rows of this 2D array
259  //
260  // target coordinates have the format:
261  // dimensions: (# number of target sites) X (dimension)
262  // # of target sites is same as # of rows of neighbor lists
263  //
264  my_GMLS.setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
265  my_GMLS.setTangentBundle(tangent_bundles_device);
266 
267  // create a vector of target operations
268  TargetOperation lro;
270 
271  // and then pass them to the GMLS class
272  my_GMLS.addTargets(lro);
273 
274  // sets the weighting kernel function from WeightingFunctionType
276 
277  // power to use in that weighting kernel function
278  my_GMLS.setWeightingParameter(2);
279 
280  // generate the alphas that to be combined with data for each target operation requested in lro
281  my_GMLS.generateAlphas(number_of_batches);
282 
283  //! [Setting Up The GMLS Object]
284 
285  double instantiation_time = timer.seconds();
286  std::cout << "Took " << instantiation_time << "s to complete normal vectors generation." << std::endl;
287  Kokkos::fence(); // let generateNormalVectors finish up before using alphas
288  Kokkos::Profiling::pushRegion("Apply Alphas to Data");
289 
290  //! [Apply GMLS Alphas To Data]
291 
292  // it is important to note that if you expect to use the data as a 1D view, then you should use double*
293  // however, if you know that the target operation will result in a 2D view (vector or matrix output),
294  // then you should template with double** as this is something that can not be infered from the input data
295  // or the target operator at compile time. Additionally, a template argument is required indicating either
296  // Kokkos::HostSpace or Kokkos::DefaultExecutionSpace::memory_space()
297 
298  // The Evaluator class takes care of handling input data views as well as the output data views.
299  // It uses information from the GMLS class to determine how many components are in the input
300  // as well as output for any choice of target functionals and then performs the contactions
301  // on the data using the alpha coefficients generated by the GMLS class, all on the device.
302  Evaluator gmls_evaluator(&my_GMLS);
303 
304  auto output_value = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
305  (sampling_data_device, LaplacianOfScalarPointEvaluation);
306 
307  Kokkos::fence(); // let application of alphas to data finish before using results
308  Kokkos::Profiling::popRegion();
309  // times the Comparison in Kokkos
310  Kokkos::Profiling::pushRegion("Comparison");
311 
312  //! [Check That Solutions Are Correct]
313 
314  // loop through the target sites
315  for (int i=0; i<number_target_coords; i++) {
316  // target site i's coordinate
317  double xval = target_coords(i,0);
318  double yval = (dimension>1) ? target_coords(i,1) : 0;
319  double zval = (dimension>2) ? target_coords(i,2) : 0;
320 
321  // 0th entry is # of neighbors, which is the index beyond the last neighbor
322  int num_neigh_i = neighbor_lists(i, 0);
323  double b_i = my_GMLS.getSolutionSetHost()->getAlpha0TensorTo0Tensor(lro, i, num_neigh_i);
324 
325  // load value from output
326  double GMLS_value = output_value(i);
327 
328  // obtain the real Laplacian
329  double actual_Laplacian = trueLaplacian(xval, yval, zval, order, dimension);
330 
331  // calculate value of g to reconstruct the computed Laplacian
332  double actual_Gradient[3] = {0,0,0}; // initialized for 3, but only filled up to dimension
333  trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
334  double g = (dimension == 3) ? (1.0/sqrt(3.0))*(actual_Gradient[0] + actual_Gradient[1] + actual_Gradient[2])
335  : (1.0/sqrt(2.0))*(actual_Gradient[0] + actual_Gradient[1]);
336  double adjusted_value = GMLS_value + b_i*g;
337 
338  // check actual function value
339  if(GMLS_value!=GMLS_value || std::abs(actual_Laplacian - adjusted_value) > failure_tolerance) {
340  all_passed = false;
341  std::cout << i << " Failed Actual by: " << std::abs(actual_Laplacian - adjusted_value) << std::endl;
342  }
343  }
344 
345  //! [Check That Solutions Are Correct]
346  // popRegion hidden from tutorial
347  // stop timing comparison loop
348  Kokkos::Profiling::popRegion();
349  //! [Finalize Program]
350 } // end of code block to reduce scope, causing Kokkos View de-allocations
351 // otherwise, Views may be deallocating when we call Kokkos finalize() later
352 
353 // finalize Kokkos and MPI (if available)
354 Kokkos::finalize();
355 #ifdef COMPADRE_USE_MPI
356 MPI_Finalize();
357 #endif
358 
359 // output to user that test passed or failed
360 if(all_passed) {
361  fprintf(stdout, "Passed test \n");
362  return 0;
363 } else {
364  fprintf(stdout, "Failed test \n");
365  return -1;
366 }
367 
368 } // main
369 
370 
371 //! [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 trueLaplacian(double x, double y, double z, int order, 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 setTangentBundle(view_type tangent_directions)
(OPTIONAL) Sets orthonormal tangent directions for reconstruction on a manifold.
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...
decltype(_h_ss) * getSolutionSetHost(bool alpha_validity_check=true)
Get solution set on host.
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...
@ ScalarTaylorPolynomial
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e....
TargetOperation
Available target functionals.
@ LaplacianOfScalarPointEvaluation
Point evaluation of the laplacian of a scalar (could be on a manifold or not)