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;
56 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> additional_target_coords_device (
"additional target coordinates", 2*number_target_coords , 3);
93 Kokkos::View<double**>::HostMirror additional_target_coords = Kokkos::create_mirror_view(additional_target_coords_device);
96 Kokkos::View<int**, Kokkos::DefaultExecutionSpace> additional_target_indices_device (
"additional target indices", number_target_coords, 4 );
97 Kokkos::View<int**>::HostMirror additional_target_indices = Kokkos::create_mirror_view(additional_target_indices_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];
152 int extra_evaluation_coordinates_count = 0;
153 for(
int i=0; i<number_target_coords; i++){
156 additional_target_indices(i,0) = (i%3)+1;
159 for (
int k=0; k<(i%3+1); ++k) {
160 for (
int j=0; j<dimension; ++j) {
161 additional_target_coords(extra_evaluation_coordinates_count,j) = target_coords(i,j) + (j==0)*1e-3 + (j==1)*1e-2 + (j==1)*(-1e-1);
163 additional_target_indices(i,k+1) = extra_evaluation_coordinates_count;
164 extra_evaluation_coordinates_count++;
171 Kokkos::Profiling::popRegion();
172 Kokkos::Profiling::pushRegion(
"Creating Data");
178 Kokkos::deep_copy(source_coords_device, source_coords);
181 Kokkos::deep_copy(target_coords_device, target_coords);
184 Kokkos::deep_copy(additional_target_coords_device, additional_target_coords);
187 Kokkos::deep_copy(additional_target_indices_device, additional_target_indices);
190 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
191 source_coords_device.extent(0));
193 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> gradient_sampling_data_device(
"samples of true gradient",
194 source_coords_device.extent(0), dimension);
196 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> divergence_sampling_data_device
197 (
"samples of true solution for divergence test", source_coords_device.extent(0), dimension);
199 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
200 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int i) {
203 double xval = source_coords_device(i,0);
204 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
205 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
208 sampling_data_device(i) =
trueSolution(xval, yval, zval, order, dimension);
211 double true_grad[3] = {0,0,0};
212 trueGradient(true_grad, xval, yval,zval, order, dimension);
214 for (
int j=0; j<dimension; ++j) {
215 gradient_sampling_data_device(i,j) = true_grad[j];
226 Kokkos::Profiling::popRegion();
227 Kokkos::Profiling::pushRegion(
"Neighbor Search");
238 double epsilon_multiplier = 1.5;
239 int estimated_upper_bound_number_neighbors =
240 point_cloud_search.getEstimatedNumberNeighborsUpperBound(min_neighbors, dimension, epsilon_multiplier);
242 Kokkos::View<int**, Kokkos::DefaultExecutionSpace> neighbor_lists_device(
"neighbor lists",
243 number_target_coords, estimated_upper_bound_number_neighbors);
244 Kokkos::View<int**>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
247 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device(
"h supports", number_target_coords);
248 Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
253 point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , target_coords, neighbor_lists,
254 epsilon, min_neighbors, epsilon_multiplier);
259 Kokkos::Profiling::popRegion();
270 Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
271 Kokkos::deep_copy(epsilon_device, epsilon);
276 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
293 my_GMLS.
setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
299 std::vector<TargetOperation> lro(2);
318 double instantiation_time = timer.seconds();
319 std::cout <<
"Took " << instantiation_time <<
"s to complete alphas generation." << std::endl;
321 Kokkos::Profiling::pushRegion(
"Apply Alphas to Data");
364 Kokkos::Profiling::popRegion();
366 Kokkos::Profiling::pushRegion(
"Comparison");
380 extra_evaluation_coordinates_count = 0;
381 for (
int i=0; i<number_target_coords; i++) {
383 for (
int k=0; k<(i%3)+1; ++k) {
386 GMLS_value = output_value1(i);
388 GMLS_GradX = output_gradient1(i,0);
390 GMLS_GradY = (dimension>1) ? output_gradient1(i,1) : 0;
392 GMLS_GradZ = (dimension>2) ? output_gradient1(i,2) : 0;
395 GMLS_value = output_value2(i);
397 GMLS_GradX = output_gradient2(i,0);
399 GMLS_GradY = (dimension>1) ? output_gradient2(i,1) : 0;
401 GMLS_GradZ = (dimension>2) ? output_gradient2(i,2) : 0;
404 GMLS_value = output_value3(i);
406 GMLS_GradX = output_gradient3(i,0);
408 GMLS_GradY = (dimension>1) ? output_gradient3(i,1) : 0;
410 GMLS_GradZ = (dimension>2) ? output_gradient3(i,2) : 0;
414 double xval = additional_target_coords(extra_evaluation_coordinates_count,0);
415 double yval = additional_target_coords(extra_evaluation_coordinates_count,1);
416 double zval = additional_target_coords(extra_evaluation_coordinates_count,2);
419 double actual_value =
trueSolution(xval, yval, zval, order, dimension);
421 double actual_Gradient[3] = {0,0,0};
422 trueGradient(actual_Gradient, xval, yval, zval, order, dimension);
425 if(GMLS_value!=GMLS_value || std::abs(actual_value - GMLS_value) > failure_tolerance) {
427 std::cout << i <<
" Failed Actual by: " << std::abs(actual_value - GMLS_value) <<
" for evaluation site: " << k << std::endl;
431 if(std::abs(actual_Gradient[0] - GMLS_GradX) > failure_tolerance) {
433 std::cout << i <<
" Failed GradX by: " << std::abs(actual_Gradient[0] - GMLS_GradX) <<
" for evaluation site: " << k << std::endl;
435 if(std::abs(actual_Gradient[1] - GMLS_GradY) > failure_tolerance) {
437 std::cout << i <<
" Failed GradY by: " << std::abs(actual_Gradient[1] - GMLS_GradY) <<
" for evaluation site: " << k << std::endl;
441 if(std::abs(actual_Gradient[2] - GMLS_GradZ) > failure_tolerance) {
443 std::cout << i <<
" Failed GradZ by: " << std::abs(actual_Gradient[2] - GMLS_GradZ) <<
" for evaluation site: " << k << std::endl;
447 extra_evaluation_coordinates_count++;
455 Kokkos::Profiling::popRegion();
464 #ifdef COMPADRE_USE_MPI
470 fprintf(stdout,
"Passed test \n");
473 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)
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 setAdditionalEvaluationSitesData(view_type_1 additional_evaluation_indices, view_type_2 additional_evaluation_coordinates)
(OPTIONAL) Sets additional evaluation sites for each target site
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...
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...
@ 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.
@ GradientOfScalarPointEvaluation
Point evaluation of the gradient of a scalar.
@ ScalarPointEvaluation
Point evaluation of a scalar.
std::string constraint_name