1 #ifndef MEDUSA_BITS_APPROXIMATIONS_RBFFD_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_RBFFD_HPP_
18 template <
class RBFType,
class vec_t,
class scale_t,
class solver_t>
19 Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>
21 const std::vector<vector_t>& local_support)
const {
22 int n1 = local_support.size();
23 assert_msg(n1 > 0,
"Cannot construct a RBFFD approximation with no stencil nodes. "
24 "Did you forget to cal findSupport?.");
27 Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> M(N, N);
29 for (
int i = 0; i < n1; ++i) {
30 for (
int j = 0; j < i; ++j) {
31 M(i, j) = M(j, i) = rbf_((local_support[i]-local_support[j]).squaredNorm());
37 for (
int i = 0; i < n2; ++i) {
38 for (
int j = 0; j < n1; ++j) {
39 M(n1+i, j) = M(j, n1+i) = mon_.eval(i, local_support[j], local_support);
42 M.bottomRightCorner(n2, n2).setZero();
46 template <
class RBFType,
class vec_t,
class scale_t,
class solver_t>
48 const vec_t& point,
const std::vector<vec_t>& support) {
49 scale_ = scale_t::scale(point, support);
50 int n1 = support.size();
53 for (
int i = 0; i < n1; ++i) {
54 support_[i] = (support[i] - point) / scale_;
56 Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic> M = getMatrix(support_);
61 template <
class RBFType,
class vec_t,
class scale_t,
class solver_t>
62 Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, 1>
64 int n = support_.size();
65 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> rhs(n + mon_.size());
66 for (
int i = 0; i < n; ++i) {
67 rhs(i) = rbf_(support_[i].squaredNorm());
69 for (
int i = 0; i < mon_.size(); ++i) {
70 rhs(n+i) = mon_.evalAt0(i);
72 return getShapeFromRhs(rhs);
75 template <
class RBFType,
class vec_t,
class scale_t,
class solver_t>
76 template <
class operator_t>
77 Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, 1>
79 int n = support_.size();
80 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> rhs(n + mon_.size());
81 RBFBasis<RBFType, vec_t> rbf_basis(n, rbf_);
82 for (
int i = 0; i < n; ++i) {
83 rhs(i) = rbf_basis.evalOpAt0(i, op, support_, scale_);
85 for (
int i = 0; i < mon_.size(); ++i) {
86 rhs(n+i) = mon_.evalOpAt0(i, op, support_, scale_);
88 return getShapeFromRhs(rhs);
91 template <
class RBFType,
class vec_t,
class scale_t,
class solver_t>
92 template <
typename Derived>
93 Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, 1>
95 const Eigen::MatrixBase<Derived>& Lb)
const {
96 return solver_.solve(Lb).head(support_.size());
99 template <
class RBFType,
class vec_t,
class scale_t,
class solver_t>
100 template <
typename Derived>
101 RBFInterpolant<RBFType, vec_t>
103 const std::vector<vector_t>& support,
const Eigen::MatrixBase<Derived>& values)
const {
104 int n = support.size();
105 assert_msg(values.size() == n,
"Number of given values %d must be equal to support size %d.",
107 double scale = scale_t::scale(point, support);
109 std::vector<vec_t> local_support(n);
110 for (
int i = 0; i < n; ++i) {
111 local_support[i] = (support[i] - point) / scale;
113 Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic> M =
114 getMatrix(local_support);
117 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> rhs(n + mon_.size());
118 rhs.head(n) = values;
119 rhs.tail(mon_.size()).setZero();
120 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> coeff = solver.solve(rhs);
121 return RBFInterpolant<rbf_t, vec_t>(rbf_, mon_, point, support, scale, coeff);
127 #endif // MEDUSA_BITS_APPROXIMATIONS_RBFFD_HPP_