Compadre 1.6.4
Loading...
Searching...
No Matches
GMLS_Manifold.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>
21
22#include "GMLS_Manifold.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
32using namespace Compadre;
33
34//! [Parse Command Line Arguments]
35
36// called from command line
37int main (int argc, char* args[]) {
38
39// initializes MPI (if available) with command line arguments given
40#ifdef COMPADRE_USE_MPI
41MPI_Init(&argc, &args);
42#endif
43
44// initializes Kokkos with command line arguments given
45Kokkos::initialize(argc, args);
46
47// code block to reduce scope for all Kokkos View allocations
48// otherwise, Views may be deallocating when we call Kokkos::finalize() later
49{
50
51 CommandLineProcessor clp(argc, args);
52 auto order = clp.order;
53 auto dimension = clp.dimension;
54 auto number_target_coords = clp.number_target_coords;
55 auto constraint_name = clp.constraint_name;
56 auto solver_name = clp.solver_name;
57 auto problem_name = clp.problem_name;
58 int N_pts_on_sphere = (clp.number_source_coords>=0) ? clp.number_source_coords : 1000;
59
60 // minimum neighbors for unisolvency is the same as the size of the polynomial basis
61 // dimension has one subtracted because it is a D-1 manifold represented in D dimensions
62 const int min_neighbors = Compadre::GMLS::getNP(order, dimension-1);
63
64 //! [Parse Command Line Arguments]
65 Kokkos::Timer timer;
66 Kokkos::Profiling::pushRegion("Setup Point Data");
67 //! [Setting Up The Point Cloud]
68
69
70 // coordinates of source sites, bigger than needed then resized later
71 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device("source coordinates",
72 1.25*N_pts_on_sphere, 3);
73 Kokkos::View<double**>::host_mirror_type source_coords = Kokkos::create_mirror_view(source_coords_device);
74
75 double r = 1.0;
76
77 // set number of source coordinates from what was calculated
78 int number_source_coords;
79
80 {
81 // fill source coordinates from surface of a sphere with quasiuniform points
82 // algorithm described at https://www.cmu.edu/biolphys/deserno/pdf/sphere_equi.pdf
83 int N_count = 0;
84 double a = 4*PI*r*r/N_pts_on_sphere;
85 double d = std::sqrt(a);
86 int M_theta = std::round(PI/d);
87 double d_theta = PI/M_theta;
88 double d_phi = a/d_theta;
89 for (int i=0; i<M_theta; ++i) {
90 double theta = PI*(i + 0.5)/M_theta;
91 int M_phi = std::ceil(2*PI*std::sin(theta)/d_phi);
92 for (int j=0; j<M_phi; ++j) {
93 double phi = 2*PI*j/M_phi;
94 source_coords(N_count, 0) = theta;
95 source_coords(N_count, 1) = phi;
96 N_count++;
97 }
98 }
99
100 for (int i=0; i<N_count; ++i) {
101 double theta = source_coords(i,0);
102 double phi = source_coords(i,1);
103 source_coords(i,0) = r*std::sin(theta)*std::cos(phi);
104 source_coords(i,1) = r*std::sin(theta)*std::sin(phi);
105 source_coords(i,2) = r*cos(theta);
106 //printf("%f %f %f\n", source_coords(i,0), source_coords(i,1), source_coords(i,2));
107 }
108 number_source_coords = N_count;
109 }
110
111 // resize source_coords to the size actually needed
112 Kokkos::resize(source_coords, number_source_coords, 3);
113 Kokkos::resize(source_coords_device, number_source_coords, 3);
114
115 // coordinates of target sites
116 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device ("target coordinates",
117 number_target_coords, 3);
118 Kokkos::View<double**>::host_mirror_type target_coords = Kokkos::create_mirror_view(target_coords_device);
119
120 {
121 bool enough_pts_found = false;
122 // fill target coordinates from surface of a sphere with quasiuniform points
123 // stop when enough points are found
124 int N_count = 0;
125 double a = 4.0*PI*r*r/((double)(1.0*number_target_coords));
126 double d = std::sqrt(a);
127 int M_theta = std::round(PI/d);
128 double d_theta = PI/((double)M_theta);
129 double d_phi = a/d_theta;
130 for (int i=0; i<M_theta; ++i) {
131 double theta = PI*(i + 0.5)/M_theta;
132 int M_phi = std::ceil(2*PI*std::sin(theta)/d_phi);
133 for (int j=0; j<M_phi; ++j) {
134 double phi = 2*PI*j/M_phi;
135 target_coords(N_count, 0) = theta;
136 target_coords(N_count, 1) = phi;
137 N_count++;
138 if (N_count == number_target_coords) {
139 enough_pts_found = true;
140 break;
141 }
142 }
143 if (enough_pts_found) break;
144 }
145
146 for (int i=0; i<N_count; ++i) {
147 double theta = target_coords(i,0);
148 double phi = target_coords(i,1);
149 target_coords(i,0) = r*std::sin(theta)*std::cos(phi);
150 target_coords(i,1) = r*std::sin(theta)*std::sin(phi);
151 target_coords(i,2) = r*cos(theta);
152 //printf("%f %f %f\n", target_coords(i,0), target_coords(i,1), target_coords(i,2));
153 }
154 //for (int i=0; i<number_target_coords; ++i) {
155 // printf("%d: %f %f %f\n", i, target_coords(i,0), target_coords(i,1), target_coords(i,2));
156 //}
157 }
158
159
160 //! [Setting Up The Point Cloud]
161
162 Kokkos::Profiling::popRegion();
163 Kokkos::fence();
164 Kokkos::Profiling::pushRegion("Creating Data");
165
166 //! [Creating The Data]
167
168
169 // source coordinates need copied to device before using to construct sampling data
170 Kokkos::deep_copy(source_coords_device, source_coords);
171
172 // target coordinates copied next, because it is a convenient time to send them to device
173 Kokkos::deep_copy(target_coords_device, target_coords);
174
175 // ensure that source coordinates are sent to device before evaluating sampling data based on them
176 Kokkos::fence();
177
178
179 // need Kokkos View storing true solution (for samples)
180 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device("samples of true solution",
181 source_coords_device.extent(0));
182
183 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> ones_data_device("samples of ones",
184 source_coords_device.extent(0));
185 Kokkos::deep_copy(ones_data_device, 1.0);
186
187 // need Kokkos View storing true vector solution (for samples)
188 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> sampling_vector_data_device("samples of vector true solution",
189 source_coords_device.extent(0), 3);
190
191 Kokkos::parallel_for("Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
192 (0,source_coords.extent(0)), KOKKOS_LAMBDA(const int i) {
193
194 // coordinates of source site i
195 double xval = source_coords_device(i,0);
196 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
197 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
198
199 // data for targets with scalar input
200 sampling_data_device(i) = sphere_harmonic54(xval, yval, zval);
201
202 for (int j=0; j<3; ++j) {
203 double gradient[3] = {0,0,0};
204 gradient_sphereHarmonic54_ambient(gradient, xval, yval, zval);
205 sampling_vector_data_device(i,j) = gradient[j];
206 }
207 });
208
209
210 //! [Creating The Data]
211
212 Kokkos::Profiling::popRegion();
213 Kokkos::Profiling::pushRegion("Neighbor Search");
214
215 //! [Performing Neighbor Search]
216
217 double epsilon_multiplier = 1.9;
218
219 // Point cloud construction for neighbor search
220 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
221 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
222
223 Kokkos::View<int*> neighbor_lists_device("neighbor lists",
224 0); // first column is # of neighbors
225 Kokkos::View<int*>::host_mirror_type neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
226 // number_of_neighbors_list must be the same size as the number of target sites so that it can be populated
227 // with the number of neighbors for each target site.
228 Kokkos::View<int*> number_of_neighbors_list_device("number of neighbor lists",
229 number_target_coords); // first column is # of neighbors
230 Kokkos::View<int*>::host_mirror_type number_of_neighbors_list = Kokkos::create_mirror_view(number_of_neighbors_list_device);
231
232 // each target site has a window size
233 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device("h supports", number_target_coords);
234 Kokkos::View<double*>::host_mirror_type epsilon = Kokkos::create_mirror_view(epsilon_device);
235
236 size_t storage_size = point_cloud_search.generateCRNeighborListsFromKNNSearch(true /*dry run*/, target_coords, neighbor_lists,
237 number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
238
239 // resize neighbor_lists_device so as to be large enough to contain all neighborhoods
240 Kokkos::resize(neighbor_lists_device, storage_size);
241 neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
242
243 // query the point cloud a second time, but this time storing results into neighbor_lists
244 point_cloud_search.generateCRNeighborListsFromKNNSearch(false /*not dry run*/, target_coords, neighbor_lists,
245 number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
246
247 // Diagnostic for neighbors found
248 //int maxn = 0;
249 //int maxi = -1;
250 //for (int i=0; i<number_of_neighbors_list.extent(0); ++i) {
251 // printf("%d: %d\n", i, number_of_neighbors_list[i]);
252 // maxi = (number_of_neighbors_list[i] > maxn) ? i : maxi;
253 // maxn = (number_of_neighbors_list[i] > maxn) ? number_of_neighbors_list[i] : maxn;
254 //}
255 //printf("max at %d: %d", maxi, maxn);
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(number_of_neighbors_list_device, number_of_neighbors_list);
272 Kokkos::deep_copy(epsilon_device, epsilon);
273
274 // initialize an instance of the GMLS class for problems with a scalar basis and
275 // traditional point sampling as the sampling functional
276 GMLS my_GMLS_scalar(order, dimension,
277 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
278 order /*manifold order*/);
279
280 // pass in neighbor lists, source coordinates, target coordinates, and window sizes
281 //
282 // neighbor lists have the format:
283 // dimensions: (# number of target sites) X (# maximum number of neighbors for any given target + 1)
284 // the first column contains the number of neighbors for that rows corresponding target index
285 //
286 // source coordinates have the format:
287 // dimensions: (# number of source sites) X (dimension)
288 // entries in the neighbor lists (integers) correspond to rows of this 2D array
289 //
290 // target coordinates have the format:
291 // dimensions: (# number of target sites) X (dimension)
292 // # of target sites is same as # of rows of neighbor lists
293 //
294 my_GMLS_scalar.setProblemData(neighbor_lists_device, number_of_neighbors_list_device, source_coords_device, target_coords_device, epsilon_device);
295
296 // set a reference outward normal direction, causing a surface orientation after
297 // the GMLS instance computes an approximate tangent bundle
298 // on a sphere, the ambient coordinates are the outward normal direction
299 my_GMLS_scalar.setReferenceOutwardNormalDirection(target_coords_device, true /* use to orient surface */);
300
301 // create a vector of target operations
302 std::vector<TargetOperation> lro_scalar(4);
303 lro_scalar[0] = ScalarPointEvaluation;
304 lro_scalar[1] = LaplacianOfScalarPointEvaluation;
305 lro_scalar[2] = GradientOfScalarPointEvaluation;
306 lro_scalar[3] = GaussianCurvaturePointEvaluation;
307
308 // and then pass them to the GMLS class
309 my_GMLS_scalar.addTargets(lro_scalar);
310
311 // sets the weighting kernel function from WeightingFunctionType for curvature
312 my_GMLS_scalar.setCurvatureWeightingType(WeightingFunctionType::Power);
313
314 // power to use in the weighting kernel function for curvature coefficients
315 my_GMLS_scalar.setCurvatureWeightingParameter(2);
316
317 // sets the weighting kernel function from WeightingFunctionType
318 my_GMLS_scalar.setWeightingType(WeightingFunctionType::Power);
319
320 // power to use in that weighting kernel function
321 my_GMLS_scalar.setWeightingParameter(2);
322
323 // generate the alphas that to be combined with data for each target operation requested in lro
324 my_GMLS_scalar.generateAlphas();
325
326 Kokkos::Profiling::pushRegion("Full Polynomial Basis GMLS Solution");
327 // initialize another instance of the GMLS class for problems with a vector basis on a manifold and point
328 // evaluation of that vector as the sampling functional
329 // VectorTaylorPolynomial indicates that the basis will be a polynomial with as many components as the
330 // dimension of the manifold. This differs from another possibility, which follows this class.
331 GMLS my_GMLS_vector(ReconstructionSpace::VectorTaylorPolynomial, ManifoldVectorPointSample,
332 order, dimension,
333 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
334 order /*manifold order*/);
335
336 my_GMLS_vector.setProblemData(neighbor_lists_device, number_of_neighbors_list_device, source_coords_device, target_coords_device, epsilon_device);
337 std::vector<TargetOperation> lro_vector(2);
338 lro_vector[0] = VectorPointEvaluation;
339 lro_vector[1] = DivergenceOfVectorPointEvaluation;
340 my_GMLS_vector.addTargets(lro_vector);
341 my_GMLS_vector.setCurvatureWeightingType(WeightingFunctionType::Power);
342 my_GMLS_vector.setCurvatureWeightingParameter(2);
343 my_GMLS_vector.setWeightingType(WeightingFunctionType::Power);
344 my_GMLS_vector.setWeightingParameter(2);
345 my_GMLS_vector.generateAlphas();
346 Kokkos::Profiling::popRegion();
347
348 Kokkos::Profiling::pushRegion("Scalar Polynomial Basis Repeated to Form a Vector GMLS Solution");
349 // initialize another instance of the GMLS class for problems with a vector basis on a manifold and point
350 // evaluation of that vector as the sampling functional
351 // VectorOfScalarClonesTaylorPolynomial indicates a scalar polynomial will be solved for, since
352 // each component of the reconstructed vector are independent. This results in solving a smaller system
353 // for each GMLS problem, and is the suggested way to do vector reconstructions when sampling functionals
354 // acting on the basis would not create non-zero offdiagonal blocks.
355 //
356 // As an example, consider a 2D manifold in 3D ambient space. The GMLS problem is posed in the local chart,
357 // meaning that the problem being solved looks like
358 //
359 // [P_0 0 | where P_0 has dimension #number of neighbors for a target X #dimension of a scalar basis
360 // | 0 P_1]
361 //
362 // P_1 is similarly shaped, and for sampling functional that is a point evaluation, P_0 and P_1 are
363 // identical and their degrees of freedom in this system are disjoint, allowing us to solve for the
364 // degrees of freedom for either block independently. Additionally, the will produce the exact
365 // same polynomial coefficients for the point sampling functional, therefore it makes sense to use
366 // VectorOfScalarClonesTaylorPolynomial.
367 //
368 // In the print-out for this program, we include the timings and errors on this and VectorTaylorPolynomial
369 // in order to demonstrate that they produce exactly the same answer, but that one is much more efficient.
370 //
371 GMLS my_GMLS_vector_of_scalar_clones(ReconstructionSpace::VectorOfScalarClonesTaylorPolynomial, ManifoldVectorPointSample,
372 order, dimension,
373 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
374 order /*manifold order*/);
375
376 my_GMLS_vector_of_scalar_clones.setProblemData(neighbor_lists_device, number_of_neighbors_list_device, source_coords_device, target_coords_device, epsilon_device);
377 std::vector<TargetOperation> lro_vector_of_scalar_clones(2);
378 lro_vector_of_scalar_clones[0] = VectorPointEvaluation;
379 lro_vector_of_scalar_clones[1] = DivergenceOfVectorPointEvaluation;
380 my_GMLS_vector_of_scalar_clones.addTargets(lro_vector_of_scalar_clones);
381 my_GMLS_vector_of_scalar_clones.setCurvatureWeightingType(WeightingFunctionType::Power);
382 my_GMLS_vector_of_scalar_clones.setCurvatureWeightingParameter(2);
383 my_GMLS_vector_of_scalar_clones.setWeightingType(WeightingFunctionType::Power);
384 my_GMLS_vector_of_scalar_clones.setWeightingParameter(2);
385 my_GMLS_vector_of_scalar_clones.generateAlphas();
386 Kokkos::Profiling::popRegion();
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
399 // it is important to note that if you expect to use the data as a 1D view, then you should use double*
400 // however, if you know that the target operation will result in a 2D view (vector or matrix output),
401 // then you should template with double** as this is something that can not be infered from the input data
402 // or the target operator at compile time. Additionally, a template argument is required indicating either
403 // Kokkos::HostSpace or Kokkos::DefaultExecutionSpace::memory_space()
404
405 // The Evaluator class takes care of handling input data views as well as the output data views.
406 // It uses information from the GMLS class to determine how many components are in the input
407 // as well as output for any choice of target functionals and then performs the contactions
408 // on the data using the alpha coefficients generated by the GMLS class, all on the device.
409 Evaluator scalar_gmls_evaluator(&my_GMLS_scalar);
410 Evaluator vector_gmls_evaluator(&my_GMLS_vector);
411 Evaluator vector_gmls_evaluator_of_scalar_clones(&my_GMLS_vector_of_scalar_clones);
412
413 auto output_value = scalar_gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
414 (sampling_data_device, ScalarPointEvaluation);
415
416 auto output_laplacian = scalar_gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
417 (sampling_data_device, LaplacianOfScalarPointEvaluation);
418
419 auto output_gradient = scalar_gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
420 (sampling_data_device, GradientOfScalarPointEvaluation);
421
422 auto output_gc = scalar_gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
423 (ones_data_device, GaussianCurvaturePointEvaluation);
424
425 auto output_vector = vector_gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
426 (sampling_vector_data_device, VectorPointEvaluation, ManifoldVectorPointSample);
427
428 auto output_divergence = vector_gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
429 (sampling_vector_data_device, DivergenceOfVectorPointEvaluation, ManifoldVectorPointSample);
430
431 auto output_vector_of_scalar_clones =
432 vector_gmls_evaluator_of_scalar_clones.applyAlphasToDataAllComponentsAllTargetSites<double**,
433 Kokkos::HostSpace>(sampling_vector_data_device, VectorPointEvaluation, ManifoldVectorPointSample);
434
435 auto output_divergence_of_scalar_clones =
436 vector_gmls_evaluator_of_scalar_clones.applyAlphasToDataAllComponentsAllTargetSites<double*,
437 Kokkos::HostSpace>(sampling_vector_data_device, DivergenceOfVectorPointEvaluation, ManifoldVectorPointSample);
438
439
440 // Kokkos::fence(); // let application of alphas to data finish before using results
441 //
442 //// move gradient data to device so that it can be transformed into velocity
443 //auto output_gradient_device_mirror = Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(), output_gradient);
444 //Kokkos::deep_copy(output_gradient_device_mirror, output_gradient);
445 //Kokkos::parallel_for("Create Velocity From Surface Gradient", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
446 // (0,target_coords.extent(0)), KOKKOS_LAMBDA(const int i) {
447 //
448 // // coordinates of target site i
449 // double xval = target_coords_device(i,0);
450 // double yval = (dimension>1) ? target_coords_device(i,1) : 0;
451 // double zval = (dimension>2) ? target_coords_device(i,2) : 0;
452
453 // double gradx = output_gradient_device_mirror(i,0);
454 // double grady = output_gradient_device_mirror(i,1);
455 // double gradz = output_gradient_device_mirror(i,2);
456 //
457 // // overwrites gradient with velocity
458 // output_gradient_device_mirror(i,0) = (grady*zval - yval*gradz);
459 // output_gradient_device_mirror(i,1) = (-(gradx*zval - xval*gradz));
460 // output_gradient_device_mirror(i,2) = (gradx*yval - xval*grady);
461 //
462 //});
463 //Kokkos::deep_copy(output_gradient, output_gradient_device_mirror);
464
465
466 //! [Apply GMLS Alphas To Data]
467
468 Kokkos::fence(); // let application of alphas to data finish before using results
469 Kokkos::Profiling::popRegion();
470 // times the Comparison in Kokkos
471 Kokkos::Profiling::pushRegion("Comparison");
472
473 //! [Check That Solutions Are Correct]
474
475 double tangent_bundle_error = 0;
476 double tangent_bundle_norm = 0;
477 double values_error = 0;
478 double values_norm = 0;
479 double laplacian_error = 0;
480 double laplacian_norm = 0;
481 double gradient_ambient_error = 0;
482 double gradient_ambient_norm = 0;
483 double gc_error = 0;
484 double gc_norm = 0;
485 double vector_ambient_error = 0;
486 double vector_ambient_norm = 0;
487 double divergence_ambient_error = 0;
488 double divergence_ambient_norm = 0;
489 double vector_of_scalar_clones_ambient_error = 0;
490 double vector_of_scalar_clones_ambient_norm = 0;
491 double divergence_of_scalar_clones_ambient_error = 0;
492 double divergence_of_scalar_clones_ambient_norm = 0;
493
494 // loop through the target sites
495 for (int i=0; i<number_target_coords; i++) {
496
497 // load value from output
498 double GMLS_value = output_value(i);
499 double GMLS_gc = output_gc(i);
500
501 // load laplacian from output
502 double GMLS_Laplacian = output_laplacian(i);
503
504 // target site i's coordinate
505 double xval = target_coords(i,0);
506 double yval = (dimension>1) ? target_coords(i,1) : 0;
507 double zval = (dimension>2) ? target_coords(i,2) : 0;
508 double coord[3] = {xval, yval, zval};
509
510 // get tangent vector and see if orthgonal to coordinate (it should be on a sphere)
511 for (unsigned int j=0; j<static_cast<unsigned int>(dimension-1); ++j) {
512 // gcc 7 chokes on int(dimension-1) with -Waggressive-loop-optimizations
513 // so we use unsigned int instead
514 double tangent_inner_prod = 0;
515 for (int k=0; k<dimension; ++k) {
516 tangent_inner_prod += coord[k] * my_GMLS_scalar.getTangentBundle(i, j, k);
517 }
518 tangent_bundle_error += tangent_inner_prod * tangent_inner_prod;
519 }
520 double normal_inner_prod = 0;
521 for (int k=0; k<dimension; ++k) {
522 normal_inner_prod += coord[k] * my_GMLS_scalar.getTangentBundle(i, dimension-1, k);
523 }
524 // inner product could be plus or minus 1 (depends on tangent direction ordering)
525 double normal_error_to_sum = (normal_inner_prod > 0) ? normal_inner_prod - 1 : normal_inner_prod + 1;
526 tangent_bundle_error += normal_error_to_sum * normal_error_to_sum;
527 tangent_bundle_norm += 1;
528
529 // evaluation of various exact solutions
530 double actual_value = sphere_harmonic54(xval, yval, zval);
531 double actual_Laplacian = laplace_beltrami_sphere_harmonic54(xval, yval, zval);
532 double actual_Gradient_ambient[3] = {0,0,0}; // initialized for 3, but only filled up to dimension
533 gradient_sphereHarmonic54_ambient(actual_Gradient_ambient, xval, yval, zval);
534 //velocity_sphereHarmonic54_ambient(actual_Gradient_ambient, xval, yval, zval);
535
536 values_error += (GMLS_value - actual_value)*(GMLS_value - actual_value);
537 values_norm += actual_value*actual_value;
538
539 laplacian_error += (GMLS_Laplacian - actual_Laplacian)*(GMLS_Laplacian - actual_Laplacian);
540 laplacian_norm += actual_Laplacian*actual_Laplacian;
541
542 double actual_gc = 1.0;
543 gc_error += (GMLS_gc - actual_gc)*(GMLS_gc - actual_gc);
544 gc_norm += actual_gc*actual_gc;
545
546 for (int j=0; j<dimension; ++j) {
547 gradient_ambient_error += (output_gradient(i,j) - actual_Gradient_ambient[j])*(output_gradient(i,j) - actual_Gradient_ambient[j]);
548 gradient_ambient_norm += actual_Gradient_ambient[j]*actual_Gradient_ambient[j];
549 }
550
551 for (int j=0; j<dimension; ++j) {
552 vector_ambient_error += (output_vector(i,j) - actual_Gradient_ambient[j])*(output_vector(i,j) - actual_Gradient_ambient[j]);
553 vector_ambient_norm += actual_Gradient_ambient[j]*actual_Gradient_ambient[j];
554 }
555
556 divergence_ambient_error += (output_divergence(i) - actual_Laplacian)*(output_divergence(i) - actual_Laplacian);
557 divergence_ambient_norm += actual_Laplacian*actual_Laplacian;
558
559 for (int j=0; j<dimension; ++j) {
560 vector_of_scalar_clones_ambient_error += (output_vector_of_scalar_clones(i,j) - actual_Gradient_ambient[j])*(output_vector_of_scalar_clones(i,j) - actual_Gradient_ambient[j]);
561 vector_of_scalar_clones_ambient_norm += actual_Gradient_ambient[j]*actual_Gradient_ambient[j];
562 }
563
564 divergence_of_scalar_clones_ambient_error += (output_divergence_of_scalar_clones(i) - actual_Laplacian)*(output_divergence_of_scalar_clones(i) - actual_Laplacian);
565 divergence_of_scalar_clones_ambient_norm += actual_Laplacian*actual_Laplacian;
566
567 }
568
569 tangent_bundle_error /= number_target_coords;
570 tangent_bundle_error = std::sqrt(tangent_bundle_error);
571 tangent_bundle_norm /= number_target_coords;
572 tangent_bundle_norm = std::sqrt(tangent_bundle_norm);
573
574 values_error /= number_target_coords;
575 values_error = std::sqrt(values_error);
576 values_norm /= number_target_coords;
577 values_norm = std::sqrt(values_norm);
578
579 laplacian_error /= number_target_coords;
580 laplacian_error = std::sqrt(laplacian_error);
581 laplacian_norm /= number_target_coords;
582 laplacian_norm = std::sqrt(laplacian_norm);
583
584 gradient_ambient_error /= number_target_coords;
585 gradient_ambient_error = std::sqrt(gradient_ambient_error);
586 gradient_ambient_norm /= number_target_coords;
587 gradient_ambient_norm = std::sqrt(gradient_ambient_norm);
588
589 gc_error /= number_target_coords;
590 gc_error = std::sqrt(gc_error);
591 gc_norm /= number_target_coords;
592 gc_norm = std::sqrt(gc_norm);
593
594 vector_ambient_error /= number_target_coords;
595 vector_ambient_error = std::sqrt(vector_ambient_error);
596 vector_ambient_norm /= number_target_coords;
597 vector_ambient_norm = std::sqrt(vector_ambient_norm);
598
599 divergence_ambient_error /= number_target_coords;
600 divergence_ambient_error = std::sqrt(divergence_ambient_error);
601 divergence_ambient_norm /= number_target_coords;
602 divergence_ambient_norm = std::sqrt(divergence_ambient_norm);
603
604 vector_of_scalar_clones_ambient_error /= number_target_coords;
605 vector_of_scalar_clones_ambient_error = std::sqrt(vector_of_scalar_clones_ambient_error);
606 vector_of_scalar_clones_ambient_norm /= number_target_coords;
607 vector_of_scalar_clones_ambient_norm = std::sqrt(vector_of_scalar_clones_ambient_norm);
608
609 divergence_of_scalar_clones_ambient_error /= number_target_coords;
610 divergence_of_scalar_clones_ambient_error = std::sqrt(divergence_of_scalar_clones_ambient_error);
611 divergence_of_scalar_clones_ambient_norm /= number_target_coords;
612 divergence_of_scalar_clones_ambient_norm = std::sqrt(divergence_of_scalar_clones_ambient_norm);
613
614 printf("Tangent Bundle Error: %g\n", tangent_bundle_error / tangent_bundle_norm);
615 printf("Point Value Error: %g\n", values_error / values_norm);
616 printf("Laplace-Beltrami Error: %g\n", laplacian_error / laplacian_norm);
617 printf("Gaussian Curvature Error: %g\n", gc_error / gc_norm);
618 printf("Surface Gradient (Ambient) Error: %g\n", gradient_ambient_error / gradient_ambient_norm);
619 printf("Surface Vector (VectorBasis) Error: %g\n", vector_ambient_error / vector_ambient_norm);
620 printf("Surface Divergence (VectorBasis) Error: %g\n", divergence_ambient_error / divergence_ambient_norm);
621 printf("Surface Vector (ScalarClones) Error: %g\n",
622 vector_of_scalar_clones_ambient_error / vector_of_scalar_clones_ambient_norm);
623 printf("Surface Divergence (ScalarClones) Error: %g\n",
624 divergence_of_scalar_clones_ambient_error / divergence_of_scalar_clones_ambient_norm);
625 //! [Check That Solutions Are Correct]
626 // popRegion hidden from tutorial
627 // stop timing comparison loop
628 Kokkos::Profiling::popRegion();
629 //! [Finalize Program]
630
631
632} // end of code block to reduce scope, causing Kokkos View de-allocations
633// otherwise, Views may be deallocating when we call Kokkos::finalize() later
634
635// finalize Kokkos and MPI (if available)
636Kokkos::finalize();
637#ifdef COMPADRE_USE_MPI
638MPI_Finalize();
639#endif
640
641return 0;
642
643} // main
644
645
646//! [Finalize Program]
int main(int argc, char *args[])
[Parse Command Line Arguments]
KOKKOS_INLINE_FUNCTION double sphere_harmonic54(double x, double y, double z)
#define PI
KOKKOS_INLINE_FUNCTION void gradient_sphereHarmonic54_ambient(double *gradient, double x, double y, double z)
KOKKOS_INLINE_FUNCTION double laplace_beltrami_sphere_harmonic54(double x, double y, double z)
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 setReferenceOutwardNormalDirection(view_type outward_normal_directions, bool use_to_orient_surface=true)
(OPTIONAL) Sets outward normal direction.
void setCurvatureWeightingParameter(int wp, int index=0)
Parameter for weighting kernel for curvature index = 0 sets p paramater for weighting kernel index = ...
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....
void setCurvatureWeightingType(const std::string &wt)
Type for weighting kernel for curvature.
double getTangentBundle(const int target_index, const int direction, const int component) const
Get component of tangent or normal directions for manifold problems.
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...
constexpr SamplingFunctional ManifoldVectorPointSample
Point evaluations of the entire vector source function (but on a manifold, so it includes a transform...
@ LaplacianOfScalarPointEvaluation
Point evaluation of the laplacian of a scalar (could be on a manifold or not)
@ GradientOfScalarPointEvaluation
Point evaluation of the gradient of a scalar.
@ GaussianCurvaturePointEvaluation
Point evaluation of Gaussian curvature.
@ DivergenceOfVectorPointEvaluation
Point evaluation of the divergence of a vector (results in a scalar)
@ VectorPointEvaluation
Point evaluation of a vector (reconstructs entire vector at once, requiring a ReconstructionSpace hav...
@ ScalarPointEvaluation
Point evaluation of a scalar.