Medusa  1.1
Coordinate Free Mehless Method implementation
test/interpolants/PUApproximant_test.cpp
#include <Eigen/Core>
#include "gtest/gtest.h"
namespace mm {
TEST(Interpolants, PUApproximant1D) {
// Domain.
BoxShape<Vec1d> b(0, 1);
DomainDiscretization<Vec1d> domain = b.discretizeBoundaryWithStep(0.01);
int support_size = 2;
domain.findSupport(FindClosest(support_size));
// Values.
Eigen::VectorXd values(2);
for (int i = 0; i < values.size(); i++) {
values[i] = static_cast<double>(i);
}
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec1d, ScaleToClosest> approx({3}, Monomials<Vec1d>(2));
// Test points.
Range<Vec1d> test_points(1);
test_points[0] = Vec1d(0.5);
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto value = PUApproximant<Vec1d>::evaluate(domain, values, test_points, effective_radius,
approx); // value = 0.5
EXPECT_NEAR(0.5, value[0], 1e-15);
}
TEST(Interpolants, PUApproximant2D) {
// Domain.
BoxShape<Vec2d> b(0, 1);
DomainDiscretization<Vec2d> domain = b.discretizeWithStep(0.05);
int support_size = 15;
domain.findSupport(FindClosest(support_size));
// Values.
Eigen::VectorXd values(domain.size());
values.setOnes();
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec2d, ScaleToClosest> approx({3}, Monomials<Vec2d>(2));
// Test points.
BoxShape<Vec2d> b_test(0, 1);
DomainDiscretization<Vec2d> domain_test = b_test.discretizeWithStep(0.1);
auto test_points = domain_test.positions();
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto interpolated_values =
PUApproximant<Vec2d>::evaluate(domain, values, test_points, effective_radius, approx);
for (auto value : interpolated_values) {
EXPECT_NEAR(1, value, 1e-15);
}
}
TEST(Interpolants, PUApproximant2DNotConstant) {
// Domain.
BoxShape<Vec2d> b(0, 1);
DomainDiscretization<Vec2d> domain = b.discretizeBoundaryWithStep(1);
int support_size = 4;
domain.findSupport(FindClosest(support_size));
// Values.
int N = domain.size();
Eigen::VectorXd values(N);
for (int i = 0; i < N; i++) {
values[i] = i;
}
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec2d, ScaleToClosest> approx({3}, Monomials<Vec2d>(2));
// Test points.
Range<Vec2d> test_points;
test_points.push_back(Vec2d(0.5, 0.5));
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto num_values =
PUApproximant<Vec2d>::evaluate(domain, values, test_points, effective_radius, approx);
EXPECT_NEAR(1, num_values[0], 1e-15);
}
TEST(Interpolants, PUApproximant3D) {
// Domain.
BoxShape<Vec3d> b(0, 1);
DomainDiscretization<Vec3d> domain = b.discretizeWithStep(0.1);
int support_size = 30;
domain.findSupport(FindClosest(support_size));
// Values.
Eigen::VectorXd values(domain.size());
values.setOnes();
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec3d, ScaleToClosest> approx({3}, Monomials<Vec3d>(2));
// Test points.
BoxShape<Vec3d> b_test(0, 1);
DomainDiscretization<Vec3d> domain_test = b_test.discretizeWithStep(0.2);
auto test_points = domain_test.positions();
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto interpolated_values =
PUApproximant<Vec3d>::evaluate(domain, values, test_points, effective_radius, approx);
for (auto value : interpolated_values) {
EXPECT_NEAR(1, value, 1e-15);
}
}
TEST(Interpolants, PUApproximantKDTree1D) {
// Domain.
BoxShape<Vec1d> b(0, 1);
DomainDiscretization<Vec1d> domain = b.discretizeBoundaryWithStep(0.01);
int support_size = 2;
domain.findSupport(FindClosest(support_size));
// Values.
Eigen::VectorXd values(2);
for (int i = 0; i < values.size(); i++) {
values[i] = static_cast<double>(i);
}
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec1d, ScaleToClosest> approx({3}, Monomials<Vec1d>(2));
// Test points.
Range<Vec1d> test_points(1);
test_points[0] = Vec1d(0.5);
KDTree<Vec1d> tree(test_points);
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto value = PUApproximant<Vec1d>::evaluate(domain, values, tree, effective_radius,
approx); // value = 0.5
EXPECT_NEAR(0.5, value[0], 1e-15);
}
TEST(Interpolants, PUApproximantKDTree2D) {
// Domain.
BoxShape<Vec2d> b(0, 1);
DomainDiscretization<Vec2d> domain = b.discretizeWithStep(0.05);
int support_size = 15;
domain.findSupport(FindClosest(support_size));
// Values.
Eigen::VectorXd values(domain.size());
values.setOnes();
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec2d, ScaleToClosest> approx({3}, Monomials<Vec2d>(2));
// Test points.
BoxShape<Vec2d> b_test(0, 1);
DomainDiscretization<Vec2d> domain_test = b_test.discretizeWithStep(0.1);
auto test_points = domain_test.positions();
KDTree<Vec2d> tree(test_points);
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto interpolated_values =
PUApproximant<Vec2d>::evaluate(domain, values, tree, effective_radius, approx);
for (auto value : interpolated_values) {
EXPECT_NEAR(1, value, 1e-15);
}
}
TEST(Interpolants, PUApproximantKDTree3D) {
// Domain.
BoxShape<Vec3d> b(0, 1);
DomainDiscretization<Vec3d> domain = b.discretizeWithStep(0.1);
int support_size = 30;
domain.findSupport(FindClosest(support_size));
// Values.
Eigen::VectorXd values(domain.size());
values.setOnes();
// Approximation engine.
RBFFD<Polyharmonic<double>, Vec3d, ScaleToClosest> approx({3}, Monomials<Vec3d>(2));
// Test points.
BoxShape<Vec3d> b_test(0, 1);
DomainDiscretization<Vec3d> domain_test = b_test.discretizeWithStep(0.2);
auto test_points = domain_test.positions();
KDTree<Vec3d> tree(test_points);
// Partition of unity approximation.
double effective_radius = 1.0; // Effective radius.
// Interpolated value.
auto interpolated_values =
PUApproximant<Vec3d>::evaluate(domain, values, tree, effective_radius, approx);
for (auto value : interpolated_values) {
EXPECT_NEAR(1, value, 1e-15);
}
}
} // namespace mm
mm
Root namespace for the whole library.
Definition: Gaussian.hpp:14
mm::PUApproximant::evaluate
static eigen_vt< scalar_t > evaluate(const DomainDiscretization< vec_t > &domain, const eigen_vt< scalar_t > &values, const Range< vec_t > &query_points, d_scalar_t radius_factor, const engine_t &engine)
Evaluates the Partion-of-unity approximant at points 'query_points' constructed on given pairs.
RBFFD.hpp
PUApproximant_fwd.hpp
numutils.hpp
FindClosest.hpp
Polyharmonic.hpp
mm::Vec1d
Vec< double, 1 > Vec1d
Convenience typedef for 1d vector of doubles.
Definition: Vec_fwd.hpp:33
mm::Vec3d
Vec< double, 3 > Vec3d
Convenience typedef for 3d vector of doubles.
Definition: Vec_fwd.hpp:35
Vec.hpp
BoxShape.hpp
mm::Vec2d
Vec< double, 2 > Vec2d
Convenience typedef for 2d vector of doubles.
Definition: Vec_fwd.hpp:34