24 #include <Compadre_Config.h>
32 #ifdef COMPADRE_USE_MPI
36 #include <Kokkos_Timer.hpp>
37 #include <Kokkos_Core.hpp>
44 int main (
int argc,
char* args[]) {
47 #ifdef COMPADRE_USE_MPI
48 MPI_Init(&argc, &args);
52 Kokkos::initialize(argc, args);
55 bool all_passed =
true;
62 auto order = clp.
order;
71 const double failure_tolerance = 1e-9;
74 const double laplacian_failure_tolerance = 1e-9;
84 double h_spacing = 0.05;
85 int n_neg1_to_1 = 2*(1/h_spacing) + 1;
88 const int number_source_coords = std::pow(n_neg1_to_1, dimension);
91 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device(
"source coordinates",
92 number_source_coords, 3);
93 Kokkos::View<double**>::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
96 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device (
"target coordinates", number_target_coords, 3);
97 Kokkos::View<double**>::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
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;
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];
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;
124 source_coords(source_index,0) = this_coord[0];
125 source_coords(source_index,1) = 0;
126 source_coords(source_index,2) = 0;
132 for(
int i=0; i<number_target_coords; i++){
135 double rand_dir[3] = {0,0,0};
137 for (
int j=0; j<dimension; ++j) {
139 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
143 for (
int j=0; j<dimension; ++j) {
144 target_coords(i,j) = rand_dir[j];
157 Kokkos::deep_copy(source_coords_device, source_coords);
160 Kokkos::deep_copy(target_coords_device, target_coords);
163 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
164 source_coords_device.extent(0));
166 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device(
"samples of true gradient",
167 source_coords_device.extent(0), dimension);
169 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
170 (
"samples of true solution for divergence test", source_coords_device.extent(0), dimension);
172 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
173 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int 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;
181 sampling_data_device(i) =
trueSolution(xval, yval, zval, order, dimension);
184 double true_grad[3] = {0,0,0};
185 trueGradient(true_grad, xval, yval,zval, order, dimension);
187 for (
int j=0; j<dimension; ++j) {
188 gradient_sampling_data_device(i,j) = true_grad[j];
205 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
209 std::vector<TargetOperation> lro(5);
217 my_GMLS.addTargets(lro);
223 my_GMLS.setWeightingParameter(2);
226 my_GMLS.setSourceSites(source_coords_device);
234 for (
int i=0; i<number_target_coords; i++) {
253 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> single_target_coords_device (
"single target coordinates", 1, 3);
254 Kokkos::View<double**>::HostMirror 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);
259 Kokkos::deep_copy(single_target_coords_device, single_target_coords);
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> single_neighbor_lists_device(
"neighbor lists",
271 1, estimated_upper_bound_number_neighbors);
272 Kokkos::View<int**>::HostMirror single_neighbor_lists = Kokkos::create_mirror_view(single_neighbor_lists_device);
275 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> single_epsilon_device(
"h supports", 1);
276 Kokkos::View<double*>::HostMirror single_epsilon = Kokkos::create_mirror_view(single_epsilon_device);
281 point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , single_target_coords,
282 single_neighbor_lists, single_epsilon, min_neighbors, epsilon_multiplier);
291 Kokkos::deep_copy(single_neighbor_lists_device, single_neighbor_lists);
292 Kokkos::deep_copy(single_epsilon_device, single_epsilon);
296 my_GMLS.setNeighborLists(single_neighbor_lists_device);
297 my_GMLS.setTargetSites(single_target_coords_device);
298 my_GMLS.setWindowSizes(single_epsilon_device);
302 my_GMLS.generateAlphas(1,
true );
307 double instantiation_time = timer.seconds();
308 std::cout <<
"Took " << instantiation_time <<
"s to complete alphas generation." << std::endl;
343 (sampling_data_device);
353 double GMLS_value = output_value(0);
356 double GMLS_Laplacian = output_laplacian(0);
362 double GMLS_GradX = scalar_coefficients(0,1)*1./single_epsilon(0);
366 double GMLS_GradY = (dimension>1) ? output_gradient(0,1) : 0;
369 double GMLS_GradZ = (dimension>2) ? output_gradient(0,2) : 0;
372 double GMLS_Divergence = output_divergence(0);
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;
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;
386 double actual_value =
trueSolution(xval, yval, zval, order, dimension);
387 double actual_Laplacian =
trueLaplacian(xval, yval, zval, order, dimension);
389 double actual_Gradient[3] = {0,0,0};
390 trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
392 double actual_Divergence;
393 actual_Divergence =
trueLaplacian(xval, yval, zval, order, dimension);
395 double actual_Curl[3] = {0,0,0};
406 if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
408 std::cout << i <<
" Failed Actual by: " << std::abs(actual_value - GMLS_value) << std::endl;
412 if(std::abs(actual_Laplacian - GMLS_Laplacian) > laplacian_failure_tolerance) {
414 std::cout << i <<
" Failed Laplacian by: " << std::abs(actual_Laplacian - GMLS_Laplacian) << std::endl;
418 if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
420 std::cout << i <<
" Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) << std::endl;
422 if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
424 std::cout << i <<
" Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) << std::endl;
428 if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
430 std::cout << i <<
" Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) << std::endl;
436 if(std::abs(actual_Divergence - GMLS_Divergence) > failure_tolerance) {
438 std::cout << i <<
" Failed Divergence by: " << std::abs(actual_Divergence - GMLS_Divergence) << std::endl;
445 tmp_diff += std::abs(actual_Curl[0] - GMLS_CurlX) + std::abs(actual_Curl[1] - GMLS_CurlY);
447 tmp_diff += std::abs(actual_Curl[2] - GMLS_CurlZ);
448 if(std::abs(tmp_diff) > failure_tolerance) {
450 std::cout << i <<
" Failed Curl by: " << std::abs(tmp_diff) << std::endl;
466 #ifdef COMPADRE_USE_MPI
472 fprintf(stdout,
"Passed test \n");
475 fprintf(stdout,
"Failed test \n");
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.
std::string constraint_name