Medusa  1.1
Coordinate Free Mehless Method implementation
RBFFD.hpp
Go to the documentation of this file.
1 #ifndef MEDUSA_BITS_APPROXIMATIONS_RBFFD_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_RBFFD_HPP_
3 
9 #include "RBFFD_fwd.hpp"
10 #include <vector>
12 #include <Eigen/LU>
13 #include "RBFInterpolant.hpp"
14 
15 namespace mm {
16 
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?.");
25  int n2 = mon_.size();
26  int N = n1 + n2;
27  Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> M(N, N);
28  // RBF part
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());
32  }
33  M(i, i) = rbf_(0);
34  }
35 
36  // Monomial part
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);
40  }
41  }
42  M.bottomRightCorner(n2, n2).setZero();
43  return M;
44 }
45 
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();
51  // Local scaled support
52  support_.resize(n1);
53  for (int i = 0; i < n1; ++i) {
54  support_[i] = (support[i] - point) / scale_;
55  }
56  Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic> M = getMatrix(support_);
57  point_ = point;
58  solver_.compute(M);
59 }
60 
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());
68  }
69  for (int i = 0; i < mon_.size(); ++i) {
70  rhs(n+i) = mon_.evalAt0(i);
71  }
72  return getShapeFromRhs(rhs);
73 }
74 
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>
78 RBFFD<RBFType, vec_t, scale_t, solver_t>::getShape(const operator_t& op) const {
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_);
84  }
85  for (int i = 0; i < mon_.size(); ++i) {
86  rhs(n+i) = mon_.evalOpAt0(i, op, support_, scale_);
87  }
88  return getShapeFromRhs(rhs);
89 }
90 
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());
97 }
98 
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.",
106  values.size(), n);
107  double scale = scale_t::scale(point, support);
108  // Local scaled 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;
112  }
113  Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic> M =
114  getMatrix(local_support);
115  solver_t solver;
116  solver.compute(M);
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);
122 }
124 
125 } // namespace mm
126 
127 #endif // MEDUSA_BITS_APPROXIMATIONS_RBFFD_HPP_
mm::RBFFD::getMatrix
Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > getMatrix() const
Returns the augmented collocation matrix.
Definition: RBFFD_fwd.hpp:126
mm
Root namespace for the whole library.
Definition: Gaussian.hpp:14
RBFInterpolant.hpp
mm::RBFFD::getShapeFromRhs
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getShapeFromRhs(const Eigen::MatrixBase< Derived > &Lb) const
Returns weights for approximation with a given right hand side, which is expected to contains the val...
numutils.hpp
mm::RBFFD::compute
void compute(const vector_t &point, const std::vector< vector_t > &support)
Setup this engine to approximate operators at given point over given local neighbourhood.
assert_msg
#define assert_msg(cond,...)
Assert with better error reporting.
Definition: assert.hpp:75
mm::RBFFD::getShape
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getShape() const
Return shape function (stencil weights) for value approximation.
mm::RBFFD::getApproximant
RBFInterpolant< rbf_t, vec_t > getApproximant(const vector_t &point, const std::vector< vector_t > &support, const Eigen::MatrixBase< Derived > &values) const
Construct a RBF interpolant of form.
RBFFD_fwd.hpp