17 #include <Compadre_Config.h>
25 #ifdef COMPADRE_USE_MPI
29 #include <Kokkos_Timer.hpp>
30 #include <Kokkos_Core.hpp>
50 double local_vec[3] = {0,0};
51 for (
int j=0; j<3; ++j) {
52 local_vec[0] += T(0,j) * u[j];
53 local_vec[1] += T(1,j) * u[j];
56 for (
int j=0; j<3; ++j) u[j] = 0;
57 for (
int j=0; j<3; ++j) {
58 u[j] += P(0, j) * local_vec[0];
59 u[j] += P(1, j) * local_vec[1];
67 int main (
int argc,
char* args[]) {
70 #ifdef COMPADRE_USE_MPI
71 MPI_Init(&argc, &args);
75 Kokkos::initialize(argc, args);
82 auto order = clp.
order;
96 Kokkos::Profiling::pushRegion(
"Setup Point Data");
101 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> source_coords_device(
"source coordinates",
102 1.25*N_pts_on_sphere, 3);
103 Kokkos::View<double**>::HostMirror source_coords = Kokkos::create_mirror_view(source_coords_device);
108 int number_source_coords;
114 double a = 4*
PI*r*r/N_pts_on_sphere;
115 double d = std::sqrt(a);
116 int M_theta = std::round(
PI/d);
117 double d_theta =
PI/M_theta;
118 double d_phi = a/d_theta;
119 for (
int i=0; i<M_theta; ++i) {
120 double theta =
PI*(i + 0.5)/M_theta;
121 int M_phi = std::round(2*
PI*std::sin(theta)/d_phi);
122 for (
int j=0; j<M_phi; ++j) {
123 double phi = 2*
PI*j/M_phi;
124 source_coords(N_count, 0) = theta;
125 source_coords(N_count, 1) = phi;
130 for (
int i=0; i<N_count; ++i) {
131 double theta = source_coords(i,0);
132 double phi = source_coords(i,1);
133 source_coords(i,0) = r*std::sin(theta)*std::cos(phi);
134 source_coords(i,1) = r*std::sin(theta)*std::sin(phi);
135 source_coords(i,2) = r*cos(theta);
138 number_source_coords = N_count;
142 Kokkos::resize(source_coords, number_source_coords, 3);
143 Kokkos::resize(source_coords_device, number_source_coords, 3);
146 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> target_coords_device (
"target coordinates",
147 number_target_coords, 3);
148 Kokkos::View<double**>::HostMirror target_coords = Kokkos::create_mirror_view(target_coords_device);
151 bool enough_pts_found =
false;
155 double a = 4.0*
PI*r*r/((double)(5.0*number_target_coords));
156 double d = std::sqrt(a);
157 int M_theta = std::round(
PI/d);
158 double d_theta =
PI/((double)M_theta);
159 double d_phi = a/d_theta;
160 for (
int i=0; i<M_theta; ++i) {
161 double theta =
PI*(i + 0.5)/M_theta;
162 int M_phi = std::round(2*
PI*std::sin(theta)/d_phi);
163 for (
int j=0; j<M_phi; ++j) {
164 double phi = 2*
PI*j/M_phi;
165 target_coords(N_count, 0) = theta;
166 target_coords(N_count, 1) = phi;
168 if (N_count == number_target_coords) {
169 enough_pts_found =
true;
173 if (enough_pts_found)
break;
176 for (
int i=0; i<N_count; ++i) {
177 double theta = target_coords(i,0);
178 double phi = target_coords(i,1);
179 target_coords(i,0) = r*std::sin(theta)*std::cos(phi);
180 target_coords(i,1) = r*std::sin(theta)*std::sin(phi);
181 target_coords(i,2) = r*cos(theta);
189 Kokkos::Profiling::popRegion();
191 Kokkos::Profiling::pushRegion(
"Creating Data");
197 Kokkos::deep_copy(source_coords_device, source_coords);
200 Kokkos::deep_copy(target_coords_device, target_coords);
207 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> sampling_data_device(
"samples of true solution",
208 source_coords_device.extent(0));
209 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> ones_data_device(
"samples of ones",
210 source_coords_device.extent(0));
211 Kokkos::deep_copy(ones_data_device, 1.0);
214 Kokkos::View<double**, Kokkos::DefaultExecutionSpace> sampling_vector_data_device(
"samples of vector true solution",
215 source_coords_device.extent(0), 3);
217 Kokkos::parallel_for(
"Sampling Manufactured Solutions", Kokkos::RangePolicy<Kokkos::DefaultExecutionSpace>
218 (0,source_coords.extent(0)), KOKKOS_LAMBDA(
const int i) {
221 double xval = source_coords_device(i,0);
222 double yval = (dimension>1) ? source_coords_device(i,1) : 0;
223 double zval = (dimension>2) ? source_coords_device(i,2) : 0;
228 for (
int j=0; j<3; ++j) {
229 double gradient[3] = {0,0,0};
231 sampling_vector_data_device(i,j) = gradient[j];
238 Kokkos::Profiling::popRegion();
239 Kokkos::Profiling::pushRegion(
"Neighbor Search");
250 double epsilon_multiplier = 1.7;
251 int estimated_upper_bound_number_neighbors =
252 point_cloud_search.getEstimatedNumberNeighborsUpperBound(min_neighbors, dimension, epsilon_multiplier);
254 Kokkos::View<int**, Kokkos::DefaultExecutionSpace> neighbor_lists_device(
"neighbor lists",
255 number_target_coords, estimated_upper_bound_number_neighbors);
256 Kokkos::View<int**>::HostMirror neighbor_lists = Kokkos::create_mirror_view(neighbor_lists_device);
259 Kokkos::View<double*, Kokkos::DefaultExecutionSpace> epsilon_device(
"h supports", number_target_coords);
260 Kokkos::View<double*>::HostMirror epsilon = Kokkos::create_mirror_view(epsilon_device);
265 point_cloud_search.generate2DNeighborListsFromKNNSearch(
false , target_coords, neighbor_lists,
266 epsilon, min_neighbors, epsilon_multiplier);
270 Kokkos::Profiling::popRegion();
281 Kokkos::deep_copy(neighbor_lists_device, neighbor_lists);
282 Kokkos::deep_copy(epsilon_device, epsilon);
286 GMLS my_GMLS_scalar(order, dimension,
287 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
304 my_GMLS_scalar.
setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
316 std::vector<TargetOperation> lro_scalar(3);
339 Kokkos::Profiling::pushRegion(
"Full Polynomial Basis GMLS Solution");
346 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
349 my_GMLS_vector.
setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
351 std::vector<TargetOperation> lro_vector(2);
360 Kokkos::Profiling::popRegion();
362 Kokkos::Profiling::pushRegion(
"Scalar Polynomial Basis Repeated to Form a Vector GMLS Solution");
387 solver_name.c_str(), problem_name.c_str(), constraint_name.c_str(),
390 my_GMLS_vector_of_scalar_clones.
setProblemData(neighbor_lists_device, source_coords_device, target_coords_device, epsilon_device);
392 std::vector<TargetOperation> lro_vector_of_scalar_clones(2);
395 my_GMLS_vector_of_scalar_clones.
addTargets(lro_vector_of_scalar_clones);
401 Kokkos::Profiling::popRegion();
406 double instantiation_time = timer.seconds();
407 std::cout <<
"Took " << instantiation_time <<
"s to complete alphas generation." << std::endl;
409 Kokkos::Profiling::pushRegion(
"Apply Alphas to Data");
424 Evaluator scalar_gmls_evaluator(&my_GMLS_scalar);
425 Evaluator vector_gmls_evaluator(&my_GMLS_vector);
426 Evaluator vector_gmls_evaluator_of_scalar_clones(&my_GMLS_vector_of_scalar_clones);
456 auto output_vector_of_scalar_clones =
496 Kokkos::Profiling::popRegion();
498 Kokkos::Profiling::pushRegion(
"Comparison");
502 double tangent_bundle_error = 0;
503 double tangent_bundle_norm = 0;
504 double values_error = 0;
505 double values_norm = 0;
508 double curl_ambient_error = 0;
509 double curl_ambient_norm = 0;
514 double vector_ambient_error = 0;
515 double vector_ambient_norm = 0;
518 double vector_of_scalar_clones_ambient_error = 0;
519 double vector_of_scalar_clones_ambient_norm = 0;
525 auto prestencil_weights = Kokkos::create_mirror_view(d_prestencil_weights);
526 Kokkos::deep_copy(prestencil_weights, d_prestencil_weights);
530 auto tangent_directions = Kokkos::create_mirror_view(d_tangent_directions);
531 Kokkos::deep_copy(tangent_directions, d_tangent_directions);
534 for (
int i=0; i<number_target_coords; i++) {
542 double GMLS_value = output_value(i);
545 double GMLS_gc = output_gaussian_curvature(i);
555 double xval = source_coords(neighbor_lists(i,1),0);
556 double yval = (dimension>1) ? source_coords(neighbor_lists(i,1),1) : 0;
557 double zval = (dimension>2) ? source_coords(neighbor_lists(i,1),2) : 0;
558 double coord[3] = {xval, yval, zval};
561 for (
int j=0; j<dimension-1; ++j) {
562 double tangent_inner_prod = 0;
563 for (
int k=0; k<std::min(dimension,3); ++k) {
564 tangent_inner_prod += coord[k] * prestencil_weights(0, i, 0 , j, k);
566 tangent_bundle_error += tangent_inner_prod * tangent_inner_prod;
568 double normal_inner_prod = 0;
569 for (
int k=0; k<dimension; ++k) {
570 normal_inner_prod += coord[k] * my_GMLS_scalar.
getTangentBundle(i, dimension-1, k);
573 double normal_error_to_sum = (normal_inner_prod > 0) ? normal_inner_prod - 1 : normal_inner_prod + 1;
574 tangent_bundle_error += normal_error_to_sum * normal_error_to_sum;
575 tangent_bundle_norm += 1;
579 double actual_gc = 1.0;
581 double actual_Gradient_ambient[3] = {0,0,0};
584 double actual_Curl_ambient[3] = {0,0,0};
587 values_error += (GMLS_value - actual_value)*(GMLS_value - actual_value);
588 values_norm += actual_value*actual_value;
590 gc_error += (GMLS_gc - actual_gc)*(GMLS_gc - actual_gc);
593 for (
int j=0; j<dimension; ++j) u[j] = output_curl(i,j);
595 for (
int j=0; j<dimension; ++j) output_curl(i,j) = u[j];
597 for (
int j=0; j<dimension; ++j) {
598 curl_ambient_error += (output_curl(i,j) - actual_Curl_ambient[j])*(output_curl(i,j) - actual_Curl_ambient[j]);
599 curl_ambient_norm += actual_Curl_ambient[j]*actual_Curl_ambient[j];
610 for (
int j=0; j<dimension; ++j) u[j] = output_vector(i,j);
612 for (
int j=0; j<dimension; ++j) output_vector(i,j) = u[j];
614 for (
int j=0; j<dimension; ++j) {
615 vector_ambient_error += (output_vector(i,j) - actual_Gradient_ambient[j])*(output_vector(i,j) - actual_Gradient_ambient[j]);
616 vector_ambient_norm += actual_Gradient_ambient[j]*actual_Gradient_ambient[j];
623 for (
int j=0; j<dimension; ++j) u[j] = output_vector_of_scalar_clones(i,j);
625 for (
int j=0; j<dimension; ++j) output_vector_of_scalar_clones(i,j) = u[j];
640 for (
int j=0; j<dimension; ++j) {
641 vector_of_scalar_clones_ambient_error += (output_vector_of_scalar_clones(i,j) - actual_Gradient_ambient[j])*(output_vector_of_scalar_clones(i,j) - actual_Gradient_ambient[j]);
642 vector_of_scalar_clones_ambient_norm += actual_Gradient_ambient[j]*actual_Gradient_ambient[j];
650 tangent_bundle_error /= number_target_coords;
651 tangent_bundle_error = std::sqrt(tangent_bundle_error);
652 tangent_bundle_norm /= number_target_coords;
653 tangent_bundle_norm = std::sqrt(tangent_bundle_norm);
655 values_error /= number_target_coords;
656 values_error = std::sqrt(values_error);
657 values_norm /= number_target_coords;
658 values_norm = std::sqrt(values_norm);
660 gc_error /= number_target_coords;
661 gc_error = std::sqrt(gc_error);
662 gc_norm /= number_target_coords;
663 gc_norm = std::sqrt(gc_norm);
665 curl_ambient_error /= number_target_coords;
666 curl_ambient_error = std::sqrt(curl_ambient_error);
667 curl_ambient_norm /= number_target_coords;
668 curl_ambient_norm = std::sqrt(curl_ambient_norm);
680 vector_ambient_error /= number_target_coords;
681 vector_ambient_error = std::sqrt(vector_ambient_error);
682 vector_ambient_norm /= number_target_coords;
683 vector_ambient_norm = std::sqrt(vector_ambient_norm);
690 vector_of_scalar_clones_ambient_error /= number_target_coords;
691 vector_of_scalar_clones_ambient_error = std::sqrt(vector_of_scalar_clones_ambient_error);
692 vector_of_scalar_clones_ambient_norm /= number_target_coords;
693 vector_of_scalar_clones_ambient_norm = std::sqrt(vector_of_scalar_clones_ambient_norm);
700 printf(
"Tangent Bundle Error: %g\n", tangent_bundle_error / tangent_bundle_norm);
701 printf(
"Point Value Error: %g\n", values_error / values_norm);
702 printf(
"Gaussian Curvature Error: %g\n", gc_error / gc_norm);
703 printf(
"Surface Curl (Ambient) Error: %g\n", curl_ambient_error / curl_ambient_norm);
706 printf(
"Surface Vector (VectorBasis) Error: %g\n", vector_ambient_error / vector_ambient_norm);
708 printf(
"Surface Vector (ScalarClones) Error: %g\n",
709 vector_of_scalar_clones_ambient_error / vector_of_scalar_clones_ambient_norm);
715 Kokkos::Profiling::popRegion();
724 #ifdef COMPADRE_USE_MPI
#define TO_GLOBAL(variable)
Kokkos::View< double **, layout_right, host_execution_space, Kokkos::MemoryTraits< Kokkos::Unmanaged > > host_scratch_matrix_right_type
KOKKOS_INLINE_FUNCTION double sphere_harmonic54(double x, double y, double z)
KOKKOS_INLINE_FUNCTION void gradient_sphereHarmonic54_ambient(double *gradient, double x, double y, double z)
KOKKOS_INLINE_FUNCTION void curl_sphere_harmonic54(double *curl, double x, double y, double z)
void AmbientLocalAmbient(XYZ &u, double *T_data, double *P_data)
[Ambient to Local Back To Ambient Helper Function]
int main(int argc, char *args[])
[Ambient to Local Back To Ambient Helper Function]
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 setReferenceOutwardNormalDirection(view_type outward_normal_directions, bool use_to_orient_surface=true)
(OPTIONAL) Sets outward normal direction.
void setCurvatureWeightingParameter(int wp, int index=0)
Parameter for weighting kernel for curvature index = 0 sets p paramater for weighting kernel index = ...
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....
decltype(_T) * getTangentDirections()
Get a view (device) of all tangent direction bundles.
void setCurvatureWeightingType(const std::string &wt)
Type for weighting kernel for curvature.
double getTangentBundle(const int target_index, const int direction, const int component) const
Get component of tangent or normal directions for manifold problems.
decltype(_prestencil_weights) getPrestencilWeights() const
Get a view (device) of all rank 2 preprocessing tensors This is a rank 5 tensor that is able to provi...
constexpr SamplingFunctional VaryingManifoldVectorPointSample
For integrating polynomial dotted with normal over an edge.
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...
@ VectorTaylorPolynomial
Vector polynomial basis having # of components _dimensions, or (_dimensions-1) in the case of manifol...
@ VectorOfScalarClonesTaylorPolynomial
Scalar basis reused as many times as there are components in the vector resulting in a much cheaper p...
@ CurlOfVectorPointEvaluation
Point evaluation of the curl of a vector (results in a vector)
@ GaussianCurvaturePointEvaluation
Point evaluation of Gaussian curvature.
@ VectorPointEvaluation
Point evaluation of a vector (reconstructs entire vector at once, requiring a ReconstructionSpace hav...
@ ScalarPointEvaluation
Point evaluation of a scalar.
std::string constraint_name