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