Compadre 1.6.4
Loading...
Searching...
No Matches
GMLS_SmallBatchReuse_Device.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/*
10 *
11 * This examples tests the ability to reuse a GMLS class instance,
12 * changing the target and neighbor list for the new target.
13 *
14 */
15
16#include <iostream>
17#include <string>
18#include <vector>
19#include <map>
20#include <stdlib.h>
21#include <cstdio>
22#include <random>
23
24#include <Compadre_Config.h>
25#include <Compadre_GMLS.hpp>
28
29#include "GMLS_Tutorial.hpp"
31
32#ifdef COMPADRE_USE_MPI
33#include <mpi.h>
34#endif
35
36#include <Kokkos_Timer.hpp>
37#include <Kokkos_Core.hpp>
38
39using namespace Compadre;
40
41//! [Parse Command Line Arguments]
42
43// called from command line
44int main (int argc, char* args[]) {
45
46// initializes MPI (if available) with command line arguments given
47#ifdef COMPADRE_USE_MPI
48MPI_Init(&argc, &args);
49#endif
50
51// initializes Kokkos with command line arguments given
52Kokkos::initialize(argc, args);
53
54// becomes false if the computed solution not within the failure_threshold of the actual solution
55bool all_passed = true;
56
57// code block to reduce scope for all Kokkos View allocations
58// otherwise, Views may be deallocating when we call Kokkos::finalize() later
59{
60
61 CommandLineProcessor clp(argc, args);
62 auto order = clp.order;
63 auto dimension = clp.dimension;
64 auto number_target_coords = clp.number_target_coords;
65 auto constraint_name = clp.constraint_name;
66 auto solver_name = clp.solver_name;
67 auto problem_name = clp.problem_name;
68
69 // the functions we will be seeking to reconstruct are in the span of the basis
70 // of the reconstruction space we choose for GMLS, so the error should be very small
71 const double failure_tolerance = 1e-9;
72
73 // Laplacian is a second order differential operator, which we expect to be slightly less accurate
74 const double laplacian_failure_tolerance = 1e-9;
75
76 // minimum neighbors for unisolvency is the same as the size of the polynomial basis
77 const int min_neighbors = Compadre::GMLS::getNP(order, dimension);
78
79 //! [Parse Command Line Arguments]
80 Kokkos::Timer timer;
81 //! [Setting Up The Point Cloud]
82
83 // approximate spacing of source sites
84 double h_spacing = 0.05;
85 int n_neg1_to_1 = 2*(1/h_spacing) + 1; // always odd
86
87 // number of source coordinate sites that will fill a box of [-1,1]x[-1,1]x[-1,1] with a spacing approximately h
88 const int number_source_coords = std::pow(n_neg1_to_1, dimension);
89
90 // coordinates of source sites
91 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device("source coordinates",
92 number_source_coords, 3);
93 Kokkos::View<double**>::host_mirror_type source_coords = Kokkos::create_mirror_view(source_coords_device);
94
95 // coordinates of target sites
96 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device ("target coordinates", number_target_coords, 3);
97 Kokkos::View<double**>::host_mirror_type target_coords = Kokkos::create_mirror_view(target_coords_device);
98
99
100 // fill source coordinates with a uniform grid
101 int source_index = 0;
102 double this_coord[3] = {0,0,0};
103 for (int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
104 this_coord[0] = i*h_spacing;
105 for (int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
106 this_coord[1] = j*h_spacing;
107 for (int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
108 this_coord[2] = k*h_spacing;
109 if (dimension==3) {
110 source_coords(source_index,0) = this_coord[0];
111 source_coords(source_index,1) = this_coord[1];
112 source_coords(source_index,2) = this_coord[2];
113 source_index++;
114 }
115 }
116 if (dimension==2) {
117 source_coords(source_index,0) = this_coord[0];
118 source_coords(source_index,1) = this_coord[1];
119 source_coords(source_index,2) = 0;
120 source_index++;
121 }
122 }
123 if (dimension==1) {
124 source_coords(source_index,0) = this_coord[0];
125 source_coords(source_index,1) = 0;
126 source_coords(source_index,2) = 0;
127 source_index++;
128 }
129 }
130
131 // fill target coords somewhere inside of [-0.5,0.5]x[-0.5,0.5]x[-0.5,0.5]
132 for(int i=0; i<number_target_coords; i++){
133
134 // first, we get a uniformly random distributed direction
135 double rand_dir[3] = {0,0,0};
136
137 for (int j=0; j<dimension; ++j) {
138 // rand_dir[j] is in [-0.5, 0.5]
139 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
140 }
141
142 // then we get a uniformly random radius
143 for (int j=0; j<dimension; ++j) {
144 target_coords(i,j) = rand_dir[j];
145 }
146
147 }
148
149
150 //! [Setting Up The Point Cloud]
151
152
153 //! [Creating The Data]
154
155
156 // source coordinates need copied to device before using to construct sampling data
157 Kokkos::deep_copy(source_coords_device, source_coords);
158
159 // target coordinates copied next, because it is a convenient time to send them to device
160 Kokkos::deep_copy(target_coords_device, target_coords);
161
162 // need Kokkos View storing true solution
163 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device("samples of true solution",
164 source_coords_device.extent(0));
165
166 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device("samples of true gradient",
167 source_coords_device.extent(0), dimension);
168
169 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
170 ("samples of true solution for divergence test", source_coords_device.extent(0), dimension);
171
172 Kokkos::parallel_for("Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
173 (0,source_coords.extent(0)), KOKKOS_LAMBDA(const int i) {
174
175 // coordinates of source site i
176 double xval = source_coords_device(i,0);
177 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
178 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
179
180 // data for targets with scalar input
181 sampling_data_device(i) = trueSolution(xval, yval, zval, order, dimension);
182
183 // data for targets with vector input (divergence)
184 double true_grad[3] = {0,0,0};
185 trueGradient(true_grad, xval, yval,zval, order, dimension);
186
187 for (int j=0; j<dimension; ++j) {
188 gradient_sampling_data_device(i,j) = true_grad[j];
189
190 // data for target with vector input (curl)
191 divergence_sampling_data_device(i,j) = divergenceTestSamples(xval, yval, zval, j, dimension);
192 }
193
194 });
195
196
197 //! [Creating The Data]
198
199
200 //! [Setting Up The GMLS Object]
201
202 // initialize an instance of the GMLS class
204 order, dimension,
205 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
206 2 /*manifold order*/);
207
208 // create a vector of target operations
209 std::vector<TargetOperation> lro(5);
210 lro[0] = ScalarPointEvaluation;
215
216 // and then pass them to the GMLS class
217 my_GMLS.addTargets(lro);
218
219 // sets the weighting kernel function from WeightingFunctionType
220 my_GMLS.setWeightingType(WeightingFunctionType::Power);
221
222 // power to use in that weighting kernel function
223 my_GMLS.setWeightingParameter(2);
224
225 // set source sites once for all targets
226 my_GMLS.setSourceSites(source_coords_device);
227
228
229 // Point cloud construction for neighbor search
230 // CreatePointCloudSearch constructs an object of type PointCloudSearch, but deduces the templates for you
231 auto point_cloud_search(CreatePointCloudSearch(source_coords, dimension));
232
233 // loop through the target sites
234 for (int i=0; i<number_target_coords; i++) {
235 timer.reset();
236 // pass in neighbor lists, source coordinates, target coordinates, and window sizes
237 //
238 // single neighbor lists have the format:
239 // dimensions: (1 single neighbor list for one target) X (# maximum number of neighbors for any given target + 1)
240 // the first column contains the number of neighbors for that rows corresponding target index
241 //
242 // source coordinates have the format:
243 // dimensions: (# number of source sites) X (dimension)
244 // entries in the neighbor lists (integers) correspond to rows of this 2D array
245 //
246 // single target coordinates have the format:
247 // dimensions: (1 single target site) X (dimension)
248 // # of target sites is same as # of rows of neighbor lists
249 //
250
251
252 // coordinates of single target sites
253 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> single_target_coords_device ("single target coordinates", 1, 3);
254 Kokkos::View<double**>::host_mirror_type single_target_coords = Kokkos::create_mirror_view(single_target_coords_device);
255 for (int j=0; j<3; ++j) {
256 single_target_coords(0,j) = target_coords(i,j);
257 }
258 // target coordinates copied next, because it is a convenient time to send them to device
259 Kokkos::deep_copy(single_target_coords_device, single_target_coords);
260 Kokkos::fence();
261
262 //! [Performing Neighbor Search]
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> single_neighbor_lists_device("neighbor lists",
271 1, estimated_upper_bound_number_neighbors); // first column is # of neighbors
272 Kokkos::View<int**>::host_mirror_type single_neighbor_lists = Kokkos::create_mirror_view(single_neighbor_lists_device);
273
274 // each target site has a window size
275 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> single_epsilon_device("h supports", 1);
276 Kokkos::View<double*>::host_mirror_type single_epsilon = Kokkos::create_mirror_view(single_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*/, single_target_coords,
282 single_neighbor_lists, single_epsilon, min_neighbors, epsilon_multiplier);
283
284 //! [Performing Neighbor Search]
285
286
287 // Copy data back to device (they were filled on the host)
288 // We could have filled Kokkos Views with memory space on the host
289 // and used these instead, and then the copying of data to the device
290 // would be performed in the GMLS class
291 Kokkos::deep_copy(single_neighbor_lists_device, single_neighbor_lists);
292 Kokkos::deep_copy(single_epsilon_device, single_epsilon);
293 Kokkos::fence(); // let call to build neighbor lists complete before copying back to device
294
295 // Create temporary small views to hold just one target's information at a time
296 my_GMLS.setNeighborLists(single_neighbor_lists_device);
297 my_GMLS.setTargetSites(single_target_coords_device);
298 my_GMLS.setWindowSizes(single_epsilon_device);
299
300
301 // generate the alphas that to be combined with data for each target operation requested in lro
302 my_GMLS.generateAlphas(1, true /* keep polynomial coefficients, only needed for a test later in this program */);
303
304
305 //! [Setting Up The GMLS Object]
306
307 double instantiation_time = timer.seconds();
308 std::cout << "Took " << instantiation_time << "s to complete alphas generation." << std::endl;
309 Kokkos::fence(); // let generateAlphas finish up before using alphas
310
311
312 //! [Apply GMLS Alphas To Data]
313
314 // it is important to note that if you expect to use the data as a 1D view, then you should use double*
315 // however, if you know that the target operation will result in a 2D view (vector or matrix output),
316 // then you should template with double** as this is something that can not be infered from the input data
317 // or the target operator at compile time. Additionally, a template argument is required indicating either
318 // Kokkos::HostSpace or Kokkos::DefaultExecutionSpace::memory_space()
319
320 // The Evaluator class takes care of handling input data views as well as the output data views.
321 // It uses information from the GMLS class to determine how many components are in the input
322 // as well as output for any choice of target functionals and then performs the contactions
323 // on the data using the alpha coefficients generated by the GMLS class, all on the device.
324 Evaluator gmls_evaluator(&my_GMLS);
325
326 auto output_value = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
327 (sampling_data_device, ScalarPointEvaluation);
328
329 auto output_laplacian = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
330 (sampling_data_device, LaplacianOfScalarPointEvaluation);
331
332 auto output_gradient = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
333 (sampling_data_device, GradientOfScalarPointEvaluation);
334
335 auto output_divergence = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double*, Kokkos::HostSpace>
336 (gradient_sampling_data_device, DivergenceOfVectorPointEvaluation, VectorPointSample);
337
338 auto output_curl = gmls_evaluator.applyAlphasToDataAllComponentsAllTargetSites<double**, Kokkos::HostSpace>
339 (divergence_sampling_data_device, CurlOfVectorPointEvaluation);
340
341 // retrieves polynomial coefficients instead of remapped field
342 auto scalar_coefficients = gmls_evaluator.applyFullPolynomialCoefficientsBasisToDataAllComponents<double**, Kokkos::HostSpace>
343 (sampling_data_device);
344
345 //! [Apply GMLS Alphas To Data]
346
347 Kokkos::fence(); // let application of alphas to data finish before using results
348 // times the Comparison in Kokkos
349
350 //! [Check That Solutions Are Correct]
351
352 // load value from output
353 double GMLS_value = output_value(0);
354
355 // load laplacian from output
356 double GMLS_Laplacian = output_laplacian(0);
357
358 // load partial x from gradient
359 // this is a test that the scalar_coefficients 2d array returned hold valid entries
360 // scalar_coefficients(i,1)*1./epsilon(i) is equivalent to the target operation acting
361 // on the polynomials applied to the polynomial coefficients
362 double GMLS_GradX = scalar_coefficients(0,1)*1./single_epsilon(0);
363 //output_gradient(i,0);
364
365 // load partial y from gradient
366 double GMLS_GradY = (dimension>1) ? output_gradient(0,1) : 0;
367
368 // load partial z from gradient
369 double GMLS_GradZ = (dimension>2) ? output_gradient(0,2) : 0;
370
371 // load divergence from output
372 double GMLS_Divergence = output_divergence(0);
373
374 // load curl from output
375 double GMLS_CurlX = (dimension>1) ? output_curl(0,0) : 0;
376 double GMLS_CurlY = (dimension>1) ? output_curl(0,1) : 0;
377 double GMLS_CurlZ = (dimension>2) ? output_curl(0,2) : 0;
378
379
380 // target site i's coordinate
381 double xval = target_coords(i,0);
382 double yval = (dimension>1) ? target_coords(i,1) : 0;
383 double zval = (dimension>2) ? target_coords(i,2) : 0;
384
385 // evaluation of various exact solutions
386 double actual_value = trueSolution(xval, yval, zval, order, dimension);
387 double actual_Laplacian = trueLaplacian(xval, yval, zval, order, dimension);
388
389 double actual_Gradient[3] = {0,0,0}; // initialized for 3, but only filled up to dimension
390 trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
391
392 double actual_Divergence;
393 actual_Divergence = trueLaplacian(xval, yval, zval, order, dimension);
394
395 double actual_Curl[3] = {0,0,0}; // initialized for 3, but only filled up to dimension
396 // (and not at all for dimimension = 1)
397 if (dimension>1) {
398 actual_Curl[0] = curlTestSolution(xval, yval, zval, 0, dimension);
399 actual_Curl[1] = curlTestSolution(xval, yval, zval, 1, dimension);
400 if (dimension>2) {
401 actual_Curl[2] = curlTestSolution(xval, yval, zval, 2, dimension);
402 }
403 }
404
405 // check actual function value
406 if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
407 all_passed = false;
408 std::cout << i << " Failed Actual by: " << std::abs(actual_value - GMLS_value) << std::endl;
409 }
410
411 // check Laplacian
412 if(std::abs(actual_Laplacian - GMLS_Laplacian) > laplacian_failure_tolerance) {
413 all_passed = false;
414 std::cout << i <<" Failed Laplacian by: " << std::abs(actual_Laplacian - GMLS_Laplacian) << std::endl;
415 }
416
417 // check gradient
418 if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
419 all_passed = false;
420 std::cout << i << " Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << std::endl;
421 if (dimension>1) {
422 if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
423 all_passed = false;
424 std::cout << i << " Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << std::endl;
425 }
426 }
427 if (dimension>2) {
428 if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
429 all_passed = false;
430 std::cout << i << " Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << std::endl;
431 }
432 }
433 }
434
435 // check divergence
436 if(std::abs(actual_Divergence - GMLS_Divergence) > failure_tolerance) {
437 all_passed = false;
438 std::cout << i << " Failed Divergence by: " << std::abs(actual_Divergence - GMLS_Divergence) << std::endl;
439 }
440
441 // check curl
442 if (order > 2) { // reconstructed solution not in basis unless order greater than 2 used
443 double tmp_diff = 0;
444 if (dimension>1)
445 tmp_diff += std::abs(actual_Curl[0] - GMLS_CurlX) + std::abs(actual_Curl[1] - GMLS_CurlY);
446 if (dimension>2)
447 tmp_diff += std::abs(actual_Curl[2] - GMLS_CurlZ);
448 if(std::abs(tmp_diff) > failure_tolerance) {
449 all_passed = false;
450 std::cout << i << " Failed Curl by: " << std::abs(tmp_diff) << std::endl;
451 }
452 }
453 }
454
455
456 //! [Check That Solutions Are Correct]
457 // stop timing comparison loop
458 //! [Finalize Program]
459
460
461} // end of code block to reduce scope, causing Kokkos View de-allocations
462// otherwise, Views may be deallocating when we call Kokkos::finalize() later
463
464// finalize Kokkos and MPI (if available)
465Kokkos::finalize();
466#ifdef COMPADRE_USE_MPI
467MPI_Finalize();
468#endif
469
470// output to user that test passed or failed
471if(all_passed) {
472 fprintf(stdout, "Passed test \n");
473 return 0;
474} else {
475 fprintf(stdout, "Failed test \n");
476 return -1;
477}
478
479} // main
480
481
482//! [Finalize Program]
int main(int argc, char *args[])
[Parse Command Line Arguments]
KOKKOS_INLINE_FUNCTION double trueSolution(double x, double y, double z, int order, int dimension)
KOKKOS_INLINE_FUNCTION void trueGradient(double *ans, double x, double y, double z, int order, int dimension)
KOKKOS_INLINE_FUNCTION double 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 double trueLaplacian(double x, double y, double z, int order, int dimension)
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
Kokkos::View< output_data_type, output_array_layout, output_memory_space > applyAlphasToDataAllComponentsAllTargetSites(view_type_input_data sampling_data, TargetOperation lro, const SamplingFunctional sro_in=PointSample, bool scalar_as_vector_if_needed=true, const int evaluation_site_local_index=0) const
Transformation of data under GMLS (allocates memory for output)
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)
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...
@ VectorOfScalarClonesTaylorPolynomial
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.
@ 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.
@ CurlOfVectorPointEvaluation
Point evaluation of the curl of a vector (results in a vector)
@ DivergenceOfVectorPointEvaluation
Point evaluation of the divergence of a vector (results in a scalar)
@ ScalarPointEvaluation
Point evaluation of a scalar.