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