17 #include <Compadre_Config.h>
25 #ifdef COMPADRE_USE_MPI
29 #include <Kokkos_Timer.hpp>
30 #include <Kokkos_Core.hpp>
32 using namespace Compadre;
37 int main (
int argc,
char* args[]) {
40 #ifdef COMPADRE_USE_MPI
41 MPI_Init(&argc, &args);
45 Kokkos::initialize(argc, args);
48 bool all_passed =
true;
55 auto order = clp.
order;
64 const double failure_tolerance = 1e-9;
71 Kokkos::Profiling::pushRegion(
"Setup Point Data");
75 double h_spacing = 0.05;
76 int n_neg1_to_1 = 2*(1/h_spacing) + 1;
79 const int number_source_coords = std::pow(n_neg1_to_1, dimension);
88 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_data(
"source coordinates",
89 number_source_coords, 3);
92 number_source_coords, 3);
93 scratch_matrix_left_type::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
97 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_data(
"target coordinates",
98 number_target_coords, 3);
101 scratch_matrix_right_type::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
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;
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];
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;
128 source_coords(source_index,0) = this_coord[0];
129 source_coords(source_index,1) = 0;
130 source_coords(source_index,2) = 0;
136 for(
int i=0; i<number_target_coords; i++){
139 double rand_dir[3] = {0,0,0};
141 for (
int j=0; j<dimension; ++j) {
143 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
147 for (
int j=0; j<dimension; ++j) {
148 target_coords(i,j) = rand_dir[j];
156 Kokkos::Profiling::popRegion();
157 Kokkos::Profiling::pushRegion(
"Creating Data");
163 Kokkos::deep_copy(source_coords_device, source_coords);
166 Kokkos::deep_copy(target_coords_device, target_coords);
169 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
170 source_coords_device.extent(0));
172 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device(
"samples of true gradient",
173 source_coords_device.extent(0), dimension);
175 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
176 (
"samples of true solution for divergence test", source_coords_device.extent(0), dimension);
178 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
179 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int 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;
187 sampling_data_device(i) =
trueSolution(xval, yval, zval, order, dimension);
190 double true_grad[3] = {0,0,0};
191 trueGradient(true_grad, xval, yval,zval, order, dimension);
193 for (
int j=0; j<dimension; ++j) {
194 gradient_sampling_data_device(i,j) = true_grad[j];
205 Kokkos::Profiling::popRegion();
206 Kokkos::Profiling::pushRegion(
"Neighbor Search");
217 double epsilon_multiplier = 1.5;
218 int estimated_upper_bound_number_neighbors =
219 point_cloud_search.getEstimatedNumberNeighborsUpperBound(min_neighbors, dimension, epsilon_multiplier);
221 Kokkos::View<int**, Kokkos::DefaultExecutionSpace> neighbor_lists_device(
"neighbor lists",
222 number_target_coords, estimated_upper_bound_number_neighbors);
223 Kokkos::View<int**>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
226 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device(
"h supports", number_target_coords);
227 Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
232 point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , target_coords, neighbor_lists,
233 epsilon, min_neighbors, epsilon_multiplier);
238 Kokkos::Profiling::popRegion();
249 Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
250 Kokkos::deep_copy(epsilon_device, epsilon);
255 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
272 my_GMLS.
setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
275 std::vector<TargetOperation> lro(5);
283 my_GMLS.addTargets(lro);
289 my_GMLS.setWeightingParameter(2);
292 my_GMLS.generateAlphas(1,
true );
297 double instantiation_time = timer.seconds();
298 std::cout <<
"Took " << instantiation_time <<
"s to complete alphas generation." << std::endl;
300 Kokkos::Profiling::pushRegion(
"Apply Alphas to Data");
333 (sampling_data_device);
336 (gradient_sampling_data_device);
341 Kokkos::Profiling::popRegion();
343 Kokkos::Profiling::pushRegion(
"Comparison");
349 for (
int i=0; i<number_target_coords; i++) {
352 double GMLS_value = output_value(i);
358 double GMLS_GradX = scalar_coefficients(i,1)*1./epsilon(i);
361 double GMLS_GradY = (dimension>1) ? output_gradient(i,1) : 0;
364 double GMLS_GradZ = (dimension>2) ? output_gradient(i,2) : 0;
367 double GMLS_Divergence = output_divergence(i);
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;
374 auto NP = min_neighbors;
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;
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;
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;
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;
394 double actual_value =
trueSolution(xval, yval, zval, order, dimension);
396 double actual_Gradient[3] = {0,0,0};
397 trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
399 double actual_Divergence;
400 actual_Divergence =
trueLaplacian(xval, yval, zval, order, dimension);
402 double actual_Hessian[9] = {0,0,0,0,0,0,0,0,0};
403 trueHessian(actual_Hessian, xval, yval, zval, order, dimension);
405 double actual_Curl[3] = {0,0,0};
416 if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
418 std::cout << i <<
" Failed Actual by: " << std::abs(actual_value - GMLS_value) << std::endl;
422 if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
424 std::cout << i <<
" Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << std::endl;
427 if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
429 std::cout << i <<
" Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << std::endl;
433 if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
435 std::cout << i <<
" Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << std::endl;
440 if(std::abs(actual_Divergence - GMLS_Divergence) > failure_tolerance) {
442 std::cout << i <<
" Failed Divergence by: " << std::abs(actual_Divergence - GMLS_Divergence) << std::endl;
446 if(std::abs(actual_Hessian[0] - GMLS_GradXX) > failure_tolerance) {
448 std::cout << i <<
" Failed GradXX by: " << std::abs(actual_Hessian[0] - GMLS_GradXX) << std::endl;
451 if(std::abs(actual_Hessian[1] - GMLS_GradXY) > failure_tolerance) {
453 std::cout << i <<
" Failed GradXY by: " << std::abs(actual_Hessian[1] - GMLS_GradXY) << std::endl;
455 if(std::abs(actual_Hessian[1*dimension+1] - GMLS_GradYY) > failure_tolerance) {
457 std::cout << i <<
" Failed GradYY by: " << std::abs(actual_Hessian[1*dimension+1] - GMLS_GradYY) << std::endl;
459 if(std::abs(actual_Hessian[1*dimension+0] - GMLS_GradYX) > failure_tolerance) {
461 std::cout << i <<
" Failed GradYX by: " << std::abs(actual_Hessian[1*dimension+0] - GMLS_GradYX) << std::endl;
465 if(std::abs(actual_Hessian[2] - GMLS_GradXZ) > failure_tolerance) {
467 std::cout << i <<
" Failed GradXZ by: " << std::abs(actual_Hessian[2] - GMLS_GradXZ) << std::endl;
469 if(std::abs(actual_Hessian[1*dimension+2] - GMLS_GradYZ) > failure_tolerance) {
471 std::cout << i <<
" Failed GradYZ by: " << std::abs(actual_Hessian[1*dimension+2] - GMLS_GradYZ) << std::endl;
473 if(std::abs(actual_Hessian[2*dimension+0] - GMLS_GradZX) > failure_tolerance) {
475 std::cout << i <<
" Failed GradZX by: " << std::abs(actual_Hessian[2*dimension+0] - GMLS_GradZX) << std::endl;
477 if(std::abs(actual_Hessian[2*dimension+1] - GMLS_GradZY) > failure_tolerance) {
479 std::cout << i <<
" Failed GradZY by: " << std::abs(actual_Hessian[2*dimension+1] - GMLS_GradZY) << std::endl;
481 if(std::abs(actual_Hessian[2*dimension+2] - GMLS_GradZZ) > failure_tolerance) {
483 std::cout << i <<
" Failed GradZZ by: " << std::abs(actual_Hessian[2*dimension+2] - GMLS_GradZZ) << std::endl;
491 tmp_diff += std::abs(actual_Curl[0] - GMLS_CurlX) + std::abs(actual_Curl[1] - GMLS_CurlY);
493 tmp_diff += std::abs(actual_Curl[2] - GMLS_CurlZ);
494 if(std::abs(tmp_diff) > failure_tolerance) {
496 std::cout << i <<
" Failed Curl by: " << std::abs(tmp_diff) << std::endl;
505 Kokkos::Profiling::popRegion();
514 #ifdef COMPADRE_USE_MPI
520 fprintf(stdout,
"Passed test \n");
523 fprintf(stdout,
"Failed test \n");
Point evaluation of a scalar.
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
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...
int main(int argc, char **argv)
Point evaluation of the curl of a vector (results in a vector)
Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifol...
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...
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
Point evaluation of the divergence of a vector (results in a scalar)
Kokkos::View< double **, layout_left, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_left_type
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 (allocates memory for output)
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) ...
std::string constraint_name
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.
Point evaluation of the gradient of a vector (results in a matrix, NOT CURRENTLY IMPLEMENTED) ...
Point evaluation of a vector (reconstructs entire vector at once, requiring a ReconstructionSpace hav...