9 #include <Compadre_Config.h>
17 #ifdef COMPADRE_USE_MPI
21 #include <Kokkos_Timer.hpp>
22 #include <Kokkos_Core.hpp>
24 using namespace Compadre;
29 int main (
int argc,
char* args[]) {
32 #ifdef COMPADRE_USE_MPI
33 MPI_Init(&argc, &args);
40 bool all_passed =
true;
47 int number_of_batches = 1;
49 int arg8toi = atoi(args[7]);
51 number_of_batches = arg8toi;
54 bool keep_coefficients = number_of_batches==1;
60 int constraint_type = 0;
62 int arg7toi = atoi(args[6]);
64 constraint_type = arg7toi;
74 int arg6toi = atoi(args[5]);
76 problem_type = arg6toi;
87 int arg5toi = atoi(args[4]);
89 solver_type = arg5toi;
97 int arg4toi = atoi(args[3]);
105 int number_target_coords = 200;
107 int arg3toi = atoi(args[2]);
109 number_target_coords = arg3toi;
117 int arg2toi = atoi(args[1]);
125 const double failure_tolerance = 1e-9;
128 const double laplacian_failure_tolerance = 1e-9;
135 Kokkos::Profiling::pushRegion(
"Setup Point Data");
139 double h_spacing = 0.05;
140 int n_neg1_to_1 = 2*(1/h_spacing) + 1;
143 const int number_source_coords = std::pow(n_neg1_to_1, dimension);
146 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device(
"source coordinates",
147 number_source_coords, 3);
148 Kokkos::View<double**>::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
151 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device (
"target coordinates", number_target_coords, 3);
152 Kokkos::View<double**>::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
156 int source_index = 0;
157 double this_coord[3] = {0,0,0};
158 for (
int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
159 this_coord[0] = i*h_spacing;
160 for (
int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
161 this_coord[1] = j*h_spacing;
162 for (
int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
163 this_coord[2] = k*h_spacing;
165 source_coords(source_index,0) = this_coord[0];
166 source_coords(source_index,1) = this_coord[1];
167 source_coords(source_index,2) = this_coord[2];
172 source_coords(source_index,0) = this_coord[0];
173 source_coords(source_index,1) = this_coord[1];
174 source_coords(source_index,2) = 0;
179 source_coords(source_index,0) = this_coord[0];
180 source_coords(source_index,1) = 0;
181 source_coords(source_index,2) = 0;
187 for(
int i=0; i<number_target_coords; i++){
190 double rand_dir[3] = {0,0,0};
192 for (
int j=0; j<dimension; ++j) {
194 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
198 for (
int j=0; j<dimension; ++j) {
199 target_coords(i,j) = rand_dir[j];
207 Kokkos::Profiling::popRegion();
208 Kokkos::Profiling::pushRegion(
"Creating Data");
214 Kokkos::deep_copy(source_coords_device, source_coords);
217 Kokkos::deep_copy(target_coords_device, target_coords);
220 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
221 source_coords_device.extent(0));
223 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device(
"samples of true gradient",
224 source_coords_device.extent(0), dimension);
226 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
227 (
"samples of true solution for divergence test", source_coords_device.extent(0), dimension);
229 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
230 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int i) {
233 double xval = source_coords_device(i,0);
234 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
235 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
238 sampling_data_device(i) =
trueSolution(xval, yval, zval, order, dimension);
241 double true_grad[3] = {0,0,0};
242 trueGradient(true_grad, xval, yval,zval, order, dimension);
244 for (
int j=0; j<dimension; ++j) {
245 gradient_sampling_data_device(i,j) = true_grad[j];
256 Kokkos::Profiling::popRegion();
257 Kokkos::Profiling::pushRegion(
"Neighbor Search");
266 double epsilon_multiplier = 1.5;
271 Kokkos::View<int*> neighbor_lists_device(
"neighbor lists",
273 Kokkos::View<int*>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
276 Kokkos::View<int*> number_of_neighbors_list_device(
"number of neighbor lists",
277 number_target_coords);
278 Kokkos::View<int*>::HostMirror number_of_neighbors_list = Kokkos::create_mirror_view(number_of_neighbors_list_device);
281 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device(
"h supports", number_target_coords);
282 Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
289 size_t storage_size = point_cloud_search.generateCRNeighborListsFromKNNSearch(
true , target_coords, neighbor_lists,
290 number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
293 Kokkos::resize(neighbor_lists_device, storage_size);
294 neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
297 point_cloud_search.generateCRNeighborListsFromKNNSearch(
false , target_coords, neighbor_lists,
298 number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
302 Kokkos::Profiling::popRegion();
313 Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
314 Kokkos::deep_copy(number_of_neighbors_list_device, number_of_neighbors_list);
315 Kokkos::deep_copy(epsilon_device, epsilon);
318 std::string solver_name;
319 if (solver_type == 0) {
321 }
else if (solver_type == 1) {
323 }
else if (solver_type == 2) {
328 std::string problem_name;
329 if (problem_type == 0) {
330 problem_name =
"STANDARD";
331 }
else if (problem_type == 1) {
332 problem_name =
"MANIFOLD";
336 std::string constraint_name;
337 if (constraint_type == 0) {
338 constraint_name =
"NO_CONSTRAINT";
339 }
else if (constraint_type == 1) {
340 constraint_name =
"NEUMANN_GRAD_SCALAR";
346 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
366 my_GMLS.
setProblemData(neighbor_lists_device, number_of_neighbors_list_device, source_coords_device, target_coords_device, epsilon_device);
369 std::vector<TargetOperation> lro(5);
377 my_GMLS.addTargets(lro);
383 my_GMLS.setWeightingPower(2);
386 my_GMLS.generateAlphas(number_of_batches, keep_coefficients );
391 double instantiation_time = timer.seconds();
392 std::cout <<
"Took " << instantiation_time <<
"s to complete alphas generation." << std::endl;
394 Kokkos::Profiling::pushRegion(
"Apply Alphas to Data");
426 decltype(output_curl) scalar_coefficients;
427 if (number_of_batches==1)
428 scalar_coefficients =
430 (sampling_data_device);
435 Kokkos::Profiling::popRegion();
437 Kokkos::Profiling::pushRegion(
"Comparison");
443 for (
int i=0; i<number_target_coords; i++) {
446 double GMLS_value = output_value(i);
449 double GMLS_Laplacian = output_laplacian(i);
455 double GMLS_GradX = (number_of_batches==1) ? scalar_coefficients(i,1)*1./epsilon(i) : output_gradient(i,0);
458 double GMLS_GradY = (dimension>1) ? output_gradient(i,1) : 0;
461 double GMLS_GradZ = (dimension>2) ? output_gradient(i,2) : 0;
464 double GMLS_Divergence = output_divergence(i);
467 double GMLS_CurlX = (dimension>1) ? output_curl(i,0) : 0;
468 double GMLS_CurlY = (dimension>1) ? output_curl(i,1) : 0;
469 double GMLS_CurlZ = (dimension>2) ? output_curl(i,2) : 0;
473 double xval = target_coords(i,0);
474 double yval = (dimension>1) ? target_coords(i,1) : 0;
475 double zval = (dimension>2) ? target_coords(i,2) : 0;
478 double actual_value =
trueSolution(xval, yval, zval, order, dimension);
479 double actual_Laplacian =
trueLaplacian(xval, yval, zval, order, dimension);
481 double actual_Gradient[3] = {0,0,0};
482 trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
484 double actual_Divergence;
485 actual_Divergence =
trueLaplacian(xval, yval, zval, order, dimension);
487 double actual_Curl[3] = {0,0,0};
498 if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
500 std::cout << i <<
" Failed Actual by: " << std::abs(actual_value - GMLS_value) << std::endl;
504 if(std::abs(actual_Laplacian - GMLS_Laplacian) > laplacian_failure_tolerance) {
506 std::cout << i <<
" Failed Laplacian by: " << std::abs(actual_Laplacian - GMLS_Laplacian) << std::endl;
510 if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
512 std::cout << i <<
" Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << std::endl;
514 if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
516 std::cout << i <<
" Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << std::endl;
520 if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
522 std::cout << i <<
" Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << std::endl;
528 if(std::abs(actual_Divergence - GMLS_Divergence) > failure_tolerance) {
530 std::cout << i <<
" Failed Divergence by: " << std::abs(actual_Divergence - GMLS_Divergence) << std::endl;
537 tmp_diff += std::abs(actual_Curl[0] - GMLS_CurlX) + std::abs(actual_Curl[1] - GMLS_CurlY);
539 tmp_diff += std::abs(actual_Curl[2] - GMLS_CurlZ);
540 if(std::abs(tmp_diff) > failure_tolerance) {
542 std::cout << i <<
" Failed Curl by: " << std::abs(tmp_diff) << std::endl;
551 Kokkos::Profiling::popRegion();
560 #ifdef COMPADRE_USE_MPI
566 fprintf(stdout,
"Passed test \n");
569 fprintf(stdout,
"Failed test \n");
Point evaluation of the curl of a vector (results in a vector)
Class handling Kokkos command line arguments and returning parameters.
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
PointCloudSearch< view_type > CreatePointCloudSearch(view_type src_view, const local_index_type dimensions=-1, const local_index_type max_leaf=-1)
CreatePointCloudSearch allows for the construction of an object of type PointCloudSearch with templat...
Kokkos::View< output_data_type, output_array_layout, output_memory_space > applyFullPolynomialCoefficientsBasisToDataAllComponents(view_type_input_data sampling_data, bool scalar_as_vector_if_needed=true) const
Generation of polynomial reconstruction coefficients by applying to data in GMLS. ...
int main(int argc, char *args[])
[Parse Command Line Arguments]
Point evaluation of a scalar.
static KOKKOS_INLINE_FUNCTION int getNP(const int m, const int dimension=3, const ReconstructionSpace r_space=ReconstructionSpace::ScalarTaylorPolynomial)
Returns size of the basis for a given polynomial order and dimension General to dimension 1...
KOKKOS_INLINE_FUNCTION void trueGradient(double *ans, double x, double y, double z, int order, int dimension)
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
Point evaluation of the laplacian of a scalar (could be on a manifold or not)
Point evaluation of the divergence of a vector (results in a scalar)
Point evaluation of the gradient of a scalar.
KOKKOS_INLINE_FUNCTION double trueLaplacian(double x, double y, double z, int order, int dimension)
Generalized Moving Least Squares (GMLS)
KOKKOS_INLINE_FUNCTION double curlTestSolution(double x, double y, double z, int component, int dimension)
KOKKOS_INLINE_FUNCTION double divergenceTestSamples(double x, double y, double z, int component, int dimension)
Kokkos::View< output_data_type, output_array_layout, output_memory_space > applyAlphasToDataAllComponentsAllTargetSites(view_type_input_data sampling_data, TargetOperation lro, const SamplingFunctional sro_in=PointSample, bool scalar_as_vector_if_needed=true, const int evaluation_site_local_index=0) const
Transformation of data under GMLS.
void setProblemData(view_type_1 neighbor_lists, view_type_2 source_coordinates, view_type_3 target_coordinates, view_type_4 epsilons)
Sets basic problem data (neighbor lists, source coordinates, and target coordinates) ...
struct SubviewND< T, T2, enable_if_t<(T::rank< 2)> >{T _data_in;T2 _data_original_view;bool _scalar_as_vector_if_needed;SubviewND(T data_in, T2 data_original_view, bool scalar_as_vector_if_needed){_data_in=data_in;_data_original_view=data_original_view;_scalar_as_vector_if_needed=scalar_as_vector_if_needed;}auto get1DView(const int column_num) -> decltype(Kokkos::subview(_data_in, Kokkos::ALL))
Creates 1D subviews of data from a 1D view, generally constructed with CreateNDSliceOnDeviceView.
KOKKOS_INLINE_FUNCTION double trueSolution(double x, double y, double z, int order, int dimension)
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.