17 #include <Compadre_Config.h>
25 #ifdef COMPADRE_USE_MPI
29 #include <Kokkos_Timer.hpp>
30 #include <Kokkos_Core.hpp>
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;
65 const double failure_tolerance = 1e-9;
72 Kokkos::Profiling::pushRegion(
"Setup Point Data");
76 double h_spacing = 0.05;
77 int n_neg1_to_1 = 2*(1/h_spacing) + 1;
80 const int number_source_coords = std::pow(n_neg1_to_1, dimension);
83 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device(
"source coordinates",
84 number_source_coords, 3);
85 Kokkos::View<double**>::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
88 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device (
"target coordinates", number_target_coords, 3);
89 Kokkos::View<double**>::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
92 Kokkos::View<double***, Kokkos::DefaultExecutionSpace> tangent_bundles_device (
"tangent bundles", number_target_coords, dimension, dimension);
93 Kokkos::View<double***>::HostMirror tangent_bundles = Kokkos::create_mirror_view(tangent_bundles_device);
97 double this_coord[3] = {0,0,0};
98 for (
int i=-n_neg1_to_1/2; i<n_neg1_to_1/2+1; ++i) {
99 this_coord[0] = i*h_spacing;
100 for (
int j=-n_neg1_to_1/2; j<n_neg1_to_1/2+1; ++j) {
101 this_coord[1] = j*h_spacing;
102 for (
int k=-n_neg1_to_1/2; k<n_neg1_to_1/2+1; ++k) {
103 this_coord[2] = k*h_spacing;
105 source_coords(source_index,0) = this_coord[0];
106 source_coords(source_index,1) = this_coord[1];
107 source_coords(source_index,2) = this_coord[2];
112 source_coords(source_index,0) = this_coord[0];
113 source_coords(source_index,1) = this_coord[1];
114 source_coords(source_index,2) = 0;
119 source_coords(source_index,0) = this_coord[0];
120 source_coords(source_index,1) = 0;
121 source_coords(source_index,2) = 0;
127 for(
int i=0; i<number_target_coords; i++){
130 double rand_dir[3] = {0,0,0};
132 for (
int j=0; j<dimension; ++j) {
134 rand_dir[j] = ((double)rand() / (double) RAND_MAX) - 0.5;
138 for (
int j=0; j<dimension; ++j) {
139 target_coords(i,j) = rand_dir[j];
144 if (dimension == 2) {
145 tangent_bundles(i, 0, 0) = 0.0;
146 tangent_bundles(i, 0, 1) = 0.0;
147 tangent_bundles(i, 1, 0) = 1.0/(sqrt(2.0));
148 tangent_bundles(i, 1, 1) = 1.0/(sqrt(2.0));
149 }
else if (dimension == 3) {
150 tangent_bundles(i, 0, 0) = 0.0;
151 tangent_bundles(i, 0, 1) = 0.0;
152 tangent_bundles(i, 0, 2) = 0.0;
153 tangent_bundles(i, 1, 0) = 0.0;
154 tangent_bundles(i, 1, 1) = 0.0;
155 tangent_bundles(i, 1, 2) = 0.0;
156 tangent_bundles(i, 2, 0) = 1.0/(sqrt(3.0));
157 tangent_bundles(i, 2, 1) = 1.0/(sqrt(3.0));
158 tangent_bundles(i, 2, 2) = 1.0/(sqrt(3.0));
164 Kokkos::Profiling::popRegion();
165 Kokkos::Profiling::pushRegion(
"Creating Data");
171 Kokkos::deep_copy(source_coords_device, source_coords);
174 Kokkos::deep_copy(target_coords_device, target_coords);
177 Kokkos::deep_copy(tangent_bundles_device, tangent_bundles);
180 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
181 source_coords_device.extent(0));
183 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
184 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int i) {
187 double xval = source_coords_device(i,0);
188 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
189 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
192 sampling_data_device(i) =
trueSolution(xval, yval, zval, order, dimension);
197 Kokkos::Profiling::popRegion();
198 Kokkos::Profiling::pushRegion(
"Neighbor Search");
209 double epsilon_multiplier = 1.8;
210 int estimated_upper_bound_number_neighbors =
211 point_cloud_search.getEstimatedNumberNeighborsUpperBound(min_neighbors, dimension, epsilon_multiplier);
213 Kokkos::View<int**, Kokkos::DefaultExecutionSpace> neighbor_lists_device(
"neighbor lists",
214 number_target_coords, estimated_upper_bound_number_neighbors);
215 Kokkos::View<int**>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
218 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device(
"h supports", number_target_coords);
219 Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
224 point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , target_coords, neighbor_lists,
225 epsilon, min_neighbors, epsilon_multiplier);
229 Kokkos::Profiling::popRegion();
240 Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
241 Kokkos::deep_copy(epsilon_device, epsilon);
247 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
264 my_GMLS.
setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
285 double instantiation_time = timer.seconds();
286 std::cout <<
"Took " << instantiation_time <<
"s to complete normal vectors generation." << std::endl;
288 Kokkos::Profiling::pushRegion(
"Apply Alphas to Data");
308 Kokkos::Profiling::popRegion();
310 Kokkos::Profiling::pushRegion(
"Comparison");
315 for (
int i=0; i<number_target_coords; i++) {
317 double xval = target_coords(i,0);
318 double yval = (dimension>1) ? target_coords(i,1) : 0;
319 double zval = (dimension>2) ? target_coords(i,2) : 0;
322 int num_neigh_i = neighbor_lists(i, 0);
323 double b_i = my_GMLS.
getSolutionSetHost()->getAlpha0TensorTo0Tensor(lro, i, num_neigh_i);
326 double GMLS_value = output_value(i);
329 double actual_Laplacian =
trueLaplacian(xval, yval, zval, order, dimension);
332 double actual_Gradient[3] = {0,0,0};
333 trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
334 double g = (dimension == 3) ? (1.0/sqrt(3.0))*(actual_Gradient[0] + actual_Gradient[1] + actual_Gradient[2])
335 : (1.0/sqrt(2.0))*(actual_Gradient[0] + actual_Gradient[1]);
336 double adjusted_value = GMLS_value + b_i*g;
339 if(GMLS_value!=GMLS_value || std::abs(actual_Laplacian - adjusted_value) > failure_tolerance) {
341 std::cout << i <<
" Failed Actual by: " << std::abs(actual_Laplacian - adjusted_value) << std::endl;
348 Kokkos::Profiling::popRegion();
355 #ifdef COMPADRE_USE_MPI
361 fprintf(stdout,
"Passed test \n");
364 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 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)
Generalized Moving Least Squares (GMLS)
void addTargets(TargetOperation lro)
Adds a target to the vector of target functional to be applied to the reconstruction.
void setTangentBundle(view_type tangent_directions)
(OPTIONAL) Sets orthonormal tangent directions for reconstruction on a manifold.
void setWeightingParameter(int wp, int index=0)
Parameter for weighting kernel for GMLS problem index = 0 sets p paramater for weighting kernel index...
void generateAlphas(const int number_of_batches=1, const bool keep_coefficients=false, const bool clear_cache=true)
Meant to calculate target operations and apply the evaluations to the previously constructed polynomi...
decltype(_h_ss) * getSolutionSetHost(bool alpha_validity_check=true)
Get solution set on host.
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)
void setWeightingType(const std::string &wt)
Type for weighting kernel for GMLS problem.
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....
constexpr SamplingFunctional PointSample
Available sampling functionals.
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...
@ ScalarTaylorPolynomial
Scalar polynomial basis centered at the target site and scaled by sum of basis powers e....
TargetOperation
Available target functionals.
@ LaplacianOfScalarPointEvaluation
Point evaluation of the laplacian of a scalar (could be on a manifold or not)
std::string constraint_name