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;
62 bool keep_coefficients = number_of_batches==1;
66 const double failure_tolerance = 1e-9;
69 const double laplacian_failure_tolerance = 1e-9;
76 Kokkos::Profiling::pushRegion(
"Setup Point Data");
80 double h_spacing = 0.05;
81 int n_neg1_to_1 = 2*(1/h_spacing) + 1;
84 const int number_source_coords = std::pow(n_neg1_to_1, dimension);
87 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device(
"source coordinates",
88 number_source_coords, 3);
89 Kokkos::View<double**>::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
92 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device (
"target coordinates", number_target_coords, 3);
93 Kokkos::View<double**>::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
98 double this_coord[3] = {0,0,0};
99 for (
int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
100 this_coord[0] = i*h_spacing;
101 for (
int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
102 this_coord[1] = j*h_spacing;
103 for (
int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
104 this_coord[2] = k*h_spacing;
106 source_coords(source_index,0) = this_coord[0];
107 source_coords(source_index,1) = this_coord[1];
108 source_coords(source_index,2) = this_coord[2];
113 source_coords(source_index,0) = this_coord[0];
114 source_coords(source_index,1) = this_coord[1];
115 source_coords(source_index,2) = 0;
120 source_coords(source_index,0) = this_coord[0];
121 source_coords(source_index,1) = 0;
122 source_coords(source_index,2) = 0;
128 for(
int i=0; i<number_target_coords; i++){
131 double rand_dir[3] = {0,0,0};
133 for (
int j=0; j<dimension; ++j) {
135 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
139 for (
int j=0; j<dimension; ++j) {
140 target_coords(i,j) = rand_dir[j];
148 Kokkos::Profiling::popRegion();
149 Kokkos::Profiling::pushRegion(
"Creating Data");
155 Kokkos::deep_copy(source_coords_device, source_coords);
158 Kokkos::deep_copy(target_coords_device, target_coords);
161 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
162 source_coords_device.extent(0));
164 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device(
"samples of true gradient",
165 source_coords_device.extent(0), dimension);
167 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
168 (
"samples of true solution for divergence test", source_coords_device.extent(0), dimension);
170 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
171 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int i) {
174 double xval = source_coords_device(i,0);
175 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
176 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
179 sampling_data_device(i) =
trueSolution(xval, yval, zval, order, dimension);
182 double true_grad[3] = {0,0,0};
183 trueGradient(true_grad, xval, yval,zval, order, dimension);
185 for (
int j=0; j<dimension; ++j) {
186 gradient_sampling_data_device(i,j) = true_grad[j];
197 Kokkos::Profiling::popRegion();
198 Kokkos::Profiling::pushRegion(
"Neighbor Search");
207 double epsilon_multiplier = 1.4;
212 Kokkos::View<int*> neighbor_lists_device(
"neighbor lists",
214 Kokkos::View<int*>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
217 Kokkos::View<int*> number_of_neighbors_list_device(
"number of neighbor lists",
218 number_target_coords);
219 Kokkos::View<int*>::HostMirror number_of_neighbors_list = Kokkos::create_mirror_view(number_of_neighbors_list_device);
222 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device(
"h supports", number_target_coords);
223 Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
230 size_t storage_size = point_cloud_search.generateCRNeighborListsFromKNNSearch(
true , target_coords, neighbor_lists,
231 number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
234 Kokkos::resize(neighbor_lists_device, storage_size);
235 neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
238 point_cloud_search.generateCRNeighborListsFromKNNSearch(
false , target_coords, neighbor_lists,
239 number_of_neighbors_list, epsilon, min_neighbors, epsilon_multiplier);
243 Kokkos::Profiling::popRegion();
254 Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
255 Kokkos::deep_copy(number_of_neighbors_list_device, number_of_neighbors_list);
256 Kokkos::deep_copy(epsilon_device, epsilon);
261 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
281 my_GMLS.
setProblemData(neighbor_lists_device, number_of_neighbors_list_device, source_coords_device, target_coords_device, epsilon_device);
284 std::vector<TargetOperation> lro(5);
292 my_GMLS.addTargets(lro);
298 my_GMLS.setWeightingParameter(2);
301 my_GMLS.generateAlphas(number_of_batches, keep_coefficients );
306 double instantiation_time = timer.seconds();
307 std::cout <<
"Took " << instantiation_time <<
"s to complete alphas generation." << std::endl;
309 Kokkos::Profiling::pushRegion(
"Apply Alphas to Data");
341 decltype(output_curl) scalar_coefficients;
342 if (number_of_batches==1)
343 scalar_coefficients =
345 (sampling_data_device);
350 Kokkos::Profiling::popRegion();
352 Kokkos::Profiling::pushRegion(
"Comparison");
358 for (
int i=0; i<number_target_coords; i++) {
361 double GMLS_value = output_value(i);
364 double GMLS_Laplacian = output_laplacian(i);
370 double GMLS_GradX = (number_of_batches==1) ? scalar_coefficients(i,1)*1./epsilon(i) : output_gradient(i,0);
373 double GMLS_GradY = (dimension>1) ? output_gradient(i,1) : 0;
376 double GMLS_GradZ = (dimension>2) ? output_gradient(i,2) : 0;
379 double GMLS_Divergence = output_divergence(i);
382 double GMLS_CurlX = (dimension>1) ? output_curl(i,0) : 0;
383 double GMLS_CurlY = (dimension>1) ? output_curl(i,1) : 0;
384 double GMLS_CurlZ = (dimension>2) ? output_curl(i,2) : 0;
388 double xval = target_coords(i,0);
389 double yval = (dimension>1) ? target_coords(i,1) : 0;
390 double zval = (dimension>2) ? target_coords(i,2) : 0;
393 double actual_value =
trueSolution(xval, yval, zval, order, dimension);
394 double actual_Laplacian =
trueLaplacian(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_Curl[3] = {0,0,0};
413 if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
415 std::cout << i <<
" Failed Actual by: " << std::abs(actual_value - GMLS_value) << std::endl;
419 if(std::abs(actual_Laplacian - GMLS_Laplacian) > laplacian_failure_tolerance) {
421 std::cout << i <<
" Failed Laplacian by: " << std::abs(actual_Laplacian - GMLS_Laplacian) << std::endl;
425 if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
427 std::cout << i <<
" Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << std::endl;
429 if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
431 std::cout << i <<
" Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << std::endl;
435 if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
437 std::cout << i <<
" Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << std::endl;
443 if(std::abs(actual_Divergence - GMLS_Divergence) > failure_tolerance) {
445 std::cout << i <<
" Failed Divergence by: " << std::abs(actual_Divergence - GMLS_Divergence) << std::endl;
452 tmp_diff += std::abs(actual_Curl[0] - GMLS_CurlX) + std::abs(actual_Curl[1] - GMLS_CurlY);
454 tmp_diff += std::abs(actual_Curl[2] - GMLS_CurlZ);
455 if(std::abs(tmp_diff) > failure_tolerance) {
457 std::cout << i <<
" Failed Curl by: " << std::abs(tmp_diff) << std::endl;
466 Kokkos::Profiling::popRegion();
475 #ifdef COMPADRE_USE_MPI
481 fprintf(stdout,
"Passed test \n");
484 fprintf(stdout,
"Failed test \n");
Point evaluation of a scalar.
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
Lightweight Evaluator Helper This class is a lightweight wrapper for extracting and applying all rele...
Point evaluation of the gradient of a scalar.
Point evaluation of the laplacian of a scalar (could be on a manifold or not)
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)
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)
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 (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) ...
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.
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.