9 #include <Compadre_Config.h>
16 #ifdef COMPADRE_USE_MPI
20 #include <Kokkos_Timer.hpp>
21 #include <Kokkos_Core.hpp>
23 using namespace Compadre;
28 int main (
int argc,
char* args[]) {
31 #ifdef COMPADRE_USE_MPI
32 MPI_Init(&argc, &args);
36 Kokkos::initialize(argc, args);
39 bool all_passed =
true;
48 int constraint_type = 0;
50 int arg7toi = atoi(args[6]);
52 constraint_type = arg7toi;
62 int arg6toi = atoi(args[5]);
64 problem_type = arg6toi;
75 int arg5toi = atoi(args[4]);
77 solver_type = arg5toi;
85 int arg4toi = atoi(args[3]);
93 int number_target_coords = 200;
95 int arg3toi = atoi(args[2]);
97 number_target_coords = arg3toi;
105 int arg2toi = atoi(args[1]);
113 const double failure_tolerance = 1e-9;
120 Kokkos::Profiling::pushRegion(
"Setup Point Data");
124 double h_spacing = 0.05;
125 int n_neg1_to_1 = 2*(1/h_spacing) + 1;
128 const int number_source_coords = std::pow(n_neg1_to_1, dimension);
137 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_data(
"source coordinates",
138 number_source_coords, 3);
141 number_source_coords, 3);
142 scratch_matrix_left_type::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
146 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_data(
"target coordinates",
147 number_target_coords, 3);
150 scratch_matrix_right_type::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
154 int source_index = 0;
155 double this_coord[3] = {0,0,0};
156 for (
int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
157 this_coord[0] = i*h_spacing;
158 for (
int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
159 this_coord[1] = j*h_spacing;
160 for (
int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
161 this_coord[2] = k*h_spacing;
163 source_coords(source_index,0) = this_coord[0];
164 source_coords(source_index,1) = this_coord[1];
165 source_coords(source_index,2) = this_coord[2];
170 source_coords(source_index,0) = this_coord[0];
171 source_coords(source_index,1) = this_coord[1];
172 source_coords(source_index,2) = 0;
177 source_coords(source_index,0) = this_coord[0];
178 source_coords(source_index,1) = 0;
179 source_coords(source_index,2) = 0;
185 for(
int i=0; i<number_target_coords; i++){
188 double rand_dir[3] = {0,0,0};
190 for (
int j=0; j<dimension; ++j) {
192 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
196 for (
int j=0; j<dimension; ++j) {
197 target_coords(i,j) = rand_dir[j];
205 Kokkos::Profiling::popRegion();
206 Kokkos::Profiling::pushRegion(
"Creating Data");
212 Kokkos::deep_copy(source_coords_device, source_coords);
215 Kokkos::deep_copy(target_coords_device, target_coords);
218 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
219 source_coords_device.extent(0));
221 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device(
"samples of true gradient",
222 source_coords_device.extent(0), dimension);
224 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
225 (
"samples of true solution for divergence test", source_coords_device.extent(0), dimension);
227 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
228 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int i) {
231 double xval = source_coords_device(i,0);
232 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
233 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
236 sampling_data_device(i) =
trueSolution(xval, yval, zval, order, dimension);
239 double true_grad[3] = {0,0,0};
240 trueGradient(true_grad, xval, yval,zval, order, dimension);
242 for (
int j=0; j<dimension; ++j) {
243 gradient_sampling_data_device(i,j) = true_grad[j];
254 Kokkos::Profiling::popRegion();
255 Kokkos::Profiling::pushRegion(
"Neighbor Search");
266 double epsilon_multiplier = 1.5;
267 int estimated_upper_bound_number_neighbors =
268 point_cloud_search.getEstimatedNumberNeighborsUpperBound(min_neighbors, dimension, epsilon_multiplier);
270 Kokkos::View<int**, Kokkos::DefaultExecutionSpace> neighbor_lists_device(
"neighbor lists",
271 number_target_coords, estimated_upper_bound_number_neighbors);
272 Kokkos::View<int**>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
275 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device(
"h supports", number_target_coords);
276 Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
281 point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , target_coords, neighbor_lists,
282 epsilon, min_neighbors, epsilon_multiplier);
287 Kokkos::Profiling::popRegion();
298 Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
299 Kokkos::deep_copy(epsilon_device, epsilon);
302 std::string solver_name;
303 if (solver_type == 0) {
305 }
else if (solver_type == 1) {
307 }
else if (solver_type == 2) {
312 std::string problem_name;
313 if (problem_type == 0) {
314 problem_name =
"STANDARD";
315 }
else if (problem_type == 1) {
316 problem_name =
"MANIFOLD";
320 std::string constraint_name;
321 if (constraint_type == 0) {
322 constraint_name =
"NO_CONSTRAINT";
323 }
else if (constraint_type == 1) {
324 constraint_name =
"NEUMANN_GRAD_SCALAR";
330 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
347 my_GMLS.
setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
350 std::vector<TargetOperation> lro(5);
358 my_GMLS.addTargets(lro);
364 my_GMLS.setWeightingPower(2);
367 my_GMLS.generateAlphas(1,
true );
372 double instantiation_time = timer.seconds();
373 std::cout <<
"Took " << instantiation_time <<
"s to complete alphas generation." << std::endl;
375 Kokkos::Profiling::pushRegion(
"Apply Alphas to Data");
408 (sampling_data_device);
411 (gradient_sampling_data_device);
416 Kokkos::Profiling::popRegion();
418 Kokkos::Profiling::pushRegion(
"Comparison");
424 for (
int i=0; i<number_target_coords; i++) {
427 double GMLS_value = output_value(i);
433 double GMLS_GradX = scalar_coefficients(i,1)*1./epsilon(i);
436 double GMLS_GradY = (dimension>1) ? output_gradient(i,1) : 0;
439 double GMLS_GradZ = (dimension>2) ? output_gradient(i,2) : 0;
442 double GMLS_Divergence = output_divergence(i);
445 double GMLS_CurlX = (dimension>1) ? output_curl(i,0) : 0;
446 double GMLS_CurlY = (dimension>1) ? output_curl(i,1) : 0;
447 double GMLS_CurlZ = (dimension>2) ? output_curl(i,2) : 0;
449 auto NP = min_neighbors;
451 double GMLS_GradXX = output_hessian(i,0*dimension+0);
452 double GMLS_GradXY = (dimension>1) ? output_hessian(i,0*dimension+1) : 0;
453 double GMLS_GradXZ = (dimension>2) ? output_hessian(i,0*dimension+2) : 0;
454 double GMLS_GradYX = (dimension>1) ? output_hessian(i,1*dimension+0) : 0;
456 double GMLS_GradYY = (dimension>1) ? vector_coefficients(i,1*NP+2)*1./epsilon(i) : 0;
457 double GMLS_GradYZ = (dimension>2) ? output_hessian(i,1*dimension+2) : 0;
458 double GMLS_GradZX = (dimension>2) ? output_hessian(i,2*dimension+0) : 0;
460 double GMLS_GradZY = (dimension>2) ? vector_coefficients(i,2*NP+2)*1./epsilon(i) : 0;
461 double GMLS_GradZZ = (dimension>2) ? output_hessian(i,2*dimension+2) : 0;
464 double xval = target_coords(i,0);
465 double yval = (dimension>1) ? target_coords(i,1) : 0;
466 double zval = (dimension>2) ? target_coords(i,2) : 0;
469 double actual_value =
trueSolution(xval, yval, zval, order, dimension);
471 double actual_Gradient[3] = {0,0,0};
472 trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
474 double actual_Divergence;
475 actual_Divergence =
trueLaplacian(xval, yval, zval, order, dimension);
477 double actual_Hessian[9] = {0,0,0,0,0,0,0,0,0};
478 trueHessian(actual_Hessian, xval, yval, zval, order, dimension);
480 double actual_Curl[3] = {0,0,0};
491 if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
493 std::cout << i <<
" Failed Actual by: " << std::abs(actual_value - GMLS_value) << std::endl;
497 if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
499 std::cout << i <<
" Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << std::endl;
502 if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
504 std::cout << i <<
" Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << std::endl;
508 if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
510 std::cout << i <<
" Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << std::endl;
515 if(std::abs(actual_Divergence - GMLS_Divergence) > failure_tolerance) {
517 std::cout << i <<
" Failed Divergence by: " << std::abs(actual_Divergence - GMLS_Divergence) << std::endl;
521 if(std::abs(actual_Hessian[0] - GMLS_GradXX) > failure_tolerance) {
523 std::cout << i <<
" Failed GradXX by: " << std::abs(actual_Hessian[0] - GMLS_GradXX) << std::endl;
526 if(std::abs(actual_Hessian[1] - GMLS_GradXY) > failure_tolerance) {
528 std::cout << i <<
" Failed GradXY by: " << std::abs(actual_Hessian[1] - GMLS_GradXY) << std::endl;
530 if(std::abs(actual_Hessian[1*dimension+1] - GMLS_GradYY) > failure_tolerance) {
532 std::cout << i <<
" Failed GradYY by: " << std::abs(actual_Hessian[1*dimension+1] - GMLS_GradYY) << std::endl;
534 if(std::abs(actual_Hessian[1*dimension+0] - GMLS_GradYX) > failure_tolerance) {
536 std::cout << i <<
" Failed GradYX by: " << std::abs(actual_Hessian[1*dimension+0] - GMLS_GradYX) << std::endl;
540 if(std::abs(actual_Hessian[2] - GMLS_GradXZ) > failure_tolerance) {
542 std::cout << i <<
" Failed GradXZ by: " << std::abs(actual_Hessian[2] - GMLS_GradXZ) << std::endl;
544 if(std::abs(actual_Hessian[1*dimension+2] - GMLS_GradYZ) > failure_tolerance) {
546 std::cout << i <<
" Failed GradYZ by: " << std::abs(actual_Hessian[1*dimension+2] - GMLS_GradYZ) << std::endl;
548 if(std::abs(actual_Hessian[2*dimension+0] - GMLS_GradZX) > failure_tolerance) {
550 std::cout << i <<
" Failed GradZX by: " << std::abs(actual_Hessian[2*dimension+0] - GMLS_GradZX) << std::endl;
552 if(std::abs(actual_Hessian[2*dimension+1] - GMLS_GradZY) > failure_tolerance) {
554 std::cout << i <<
" Failed GradZY by: " << std::abs(actual_Hessian[2*dimension+1] - GMLS_GradZY) << std::endl;
556 if(std::abs(actual_Hessian[2*dimension+2] - GMLS_GradZZ) > failure_tolerance) {
558 std::cout << i <<
" Failed GradZZ by: " << std::abs(actual_Hessian[2*dimension+2] - GMLS_GradZZ) << std::endl;
566 tmp_diff += std::abs(actual_Curl[0] - GMLS_CurlX) + std::abs(actual_Curl[1] - GMLS_CurlY);
568 tmp_diff += std::abs(actual_Curl[2] - GMLS_CurlZ);
569 if(std::abs(tmp_diff) > failure_tolerance) {
571 std::cout << i <<
" Failed Curl by: " << std::abs(tmp_diff) << std::endl;
580 Kokkos::Profiling::popRegion();
589 #ifdef COMPADRE_USE_MPI
595 fprintf(stdout,
"Passed test \n");
598 fprintf(stdout,
"Failed test \n");
Point evaluation of the curl of a vector (results in a vector)
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
Point evaluation of a vector (reconstructs entire vector at once, requiring a ReconstructionSpace hav...
KOKKOS_INLINE_FUNCTION void trueHessian(double *ans, double x, double y, double z, int order, int dimension)
PointCloudSearch< view_type > CreatePointCloudSearch(view_type src_view, const local_index_type dimensions=-1, const local_index_type max_leaf=-1)
CreatePointCloudSearch allows for the construction of an object of type PointCloudSearch with templat...
Kokkos::View< output_data_type, output_array_layout, output_memory_space > applyFullPolynomialCoefficientsBasisToDataAllComponents(view_type_input_data sampling_data, bool scalar_as_vector_if_needed=true) const
Generation of polynomial reconstruction coefficients by applying to data in GMLS. ...
int main(int argc, char *args[])
[Parse Command Line Arguments]
Point evaluation of a scalar.
static KOKKOS_INLINE_FUNCTION int getNP(const int m, const int dimension=3, const ReconstructionSpace r_space=ReconstructionSpace::ScalarTaylorPolynomial)
Returns size of the basis for a given polynomial order and dimension General to dimension 1...
KOKKOS_INLINE_FUNCTION void trueGradient(double *ans, double x, double y, double z, int order, int dimension)
Kokkos::View< double **, layout_right, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_right_type
Kokkos::View< double **, layout_left, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_left_type
Point evaluation of the divergence of a vector (results in a scalar)
KOKKOS_INLINE_FUNCTION double trueLaplacian(double x, double y, double z, int order, int dimension)
Generalized Moving Least Squares (GMLS)
KOKKOS_INLINE_FUNCTION double curlTestSolution(double x, double y, double z, int component, int dimension)
KOKKOS_INLINE_FUNCTION double divergenceTestSamples(double x, double y, double z, int component, int dimension)
Kokkos::View< output_data_type, output_array_layout, output_memory_space > applyAlphasToDataAllComponentsAllTargetSites(view_type_input_data sampling_data, TargetOperation lro, const SamplingFunctional sro_in=PointSample, bool scalar_as_vector_if_needed=true, const int evaluation_site_local_index=0) const
Transformation of data under GMLS.
void setProblemData(view_type_1 neighbor_lists, view_type_2 source_coordinates, view_type_3 target_coordinates, view_type_4 epsilons)
Sets basic problem data (neighbor lists, source coordinates, and target coordinates) ...
Point evaluation of the gradient of a vector (results in a matrix, NOT CURRENTLY IMPLEMENTED) ...
KOKKOS_INLINE_FUNCTION double trueSolution(double x, double y, double z, int order, int dimension)
Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifol...
constexpr SamplingFunctional VectorPointSample
Point evaluations of the entire vector source function.