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