Medusa  1.1
Coordinate Free Mehless Method implementation
WLS.hpp
Go to the documentation of this file.
1 #ifndef MEDUSA_BITS_APPROXIMATIONS_WLS_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_WLS_HPP_
3 
9 #include "WLS_fwd.hpp"
13 #include "JacobiSVDWrapper.hpp"
14 #include "WLSApproximant.hpp"
15 
16 namespace mm {
17 
18 template <class basis_t, class weight_t, class scale_t, class solver_t>
20  const std::vector<vector_t>& support) {
21  int n = support.size();
22  assert_msg(basis_.size() <= n,
23  "WLS cannot have more basis functions than points in support domain, but %d "
24  "support points and %d basis functions were given.", n, basis_.size());
25  assert_msg(basis_.size() > 0, "Cannot construct an approximation with an empty basis.");
26  scale_ = scale_t::scale(point, support);
27  local_coordinates_.resize(n);
28  for (int i = 0; i < n; ++i) {
29  local_coordinates_[i] = (support[i] - point) / scale_;
30  }
31  point_ = point;
32  solver_.compute(getMatrix(local_coordinates_));
33 }
34 
35 template <class basis_t, class weight_t, class scale_t, class solver_t>
38  const std::vector<vector_t>& local_coordinates) const {
39  // Evaluate basis functions in given points
40  int n = local_coordinates.size();
41  int m = basis_.size();
42  ei_matrix_t WB(m, n);
43  for (int j = 0; j < n; j++) {
44  scalar_t w = weight_(local_coordinates[j]);
45  for (int i = 0; i < m; i++) {
46  WB(i, j) = w * basis_.eval(i, local_coordinates[j], local_coordinates);
47  }
48  }
49  return WB;
50 }
51 
52 template <class basis_t, class weight_t, class scale_t, class solver_t>
53 Eigen::Matrix<typename basis_t::scalar_t, Eigen::Dynamic, 1>
55  int m = basis_.size();
56  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> b(m);
57  for (int i = 0; i < m; i++) {
58  b(i) = basis_.evalAt0(i, local_coordinates_);
59  }
60  return getShapeFromRhs(b);
61 }
62 
64 template <class basis_t, class weight_t, class scale_t, class solver_t>
65 template <class operator_t>
66 Eigen::Matrix<typename basis_t::scalar_t, Eigen::Dynamic, 1>
67 WLS<basis_t, weight_t, scale_t, solver_t>::getShape(const operator_t& op) const {
68  int m = basis_.size();
69  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> Lb(m);
70  for (int i = 0; i < m; i++) {
71  Lb(i) = basis_.evalOpAt0(i, op, local_coordinates_, scale_);
72  }
73  return getShapeFromRhs(Lb);
74 }
75 
76 template <class basis_t, class weight_t, class scale_t, class solver_t>
77 template <class Derived>
78 Eigen::Matrix<typename basis_t::scalar_t, Eigen::Dynamic, 1>
80  const Eigen::MatrixBase<Derived>& Lb) const {
81  int n = local_coordinates_.size();
82  // Solve (WB)^T x = Lb
83  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> result = solver_.solve(Lb);
84  // Multiply with weight
85  for (int j = 0; j < n; j++) {
86  result(j) *= weight_(local_coordinates_[j]);
87  }
88  return result;
89 }
90 
91 template <class basis_t, class weight_t, class scale_t, class solver_t>
92 Eigen::Matrix<typename basis_t::scalar_t, Eigen::Dynamic, 1>
94  int n = local_coordinates_.size();
95  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> w(n);
96  for (int j = 0; j < n; j++) {
97  w(j) = weight_(local_coordinates_[j]);
98  }
99  return w;
100 }
102 
104 template <class basis_t, class weight_t, class scale_t, class solver_t>
105 std::ostream& operator<<(std::ostream& os, const WLS<basis_t, weight_t, scale_t, solver_t>& wls) {
106  os << "WLS:\n"
107  << " dimension: " << wls.dim << '\n'
108  << " basis size: " << wls.basis_.size() << '\n'
109  << " basis: " << wls.basis_ << '\n'
110  << " weight: " << wls.weight_ << '\n'
111  << " scale: " << scale_t() << '\n';
112  if (wls.point_[0] != wls.point_[0]) {
113  os << " last point: not used yet";
114  } else {
115  os << " last point: " << wls.point_;
116  }
117  return os;
118 }
119 
120 template <class basis_t, class weight_t, class scale_t, class solver_t>
121 template <typename Derived>
122 WLSApproximant<basis_t>
124  const vector_t& point, const std::vector<vector_t>& support,
125  const Eigen::MatrixBase<Derived>& values) const {
126  int n = support.size();
127  assert_msg(values.size() == n, "Number of given values %d must be equal to support size %d.",
128  values.size(), n);
129  assert_msg(basis_.size() <= n,
130  "WLS cannot have more basis functions than points in support domain, but %d "
131  "support points and %d basis functions were given.", n, basis_.size());
132  double scale = scale_t::scale(point, support);
133  std::vector<vector_t> local_coordinates(n);
134  for (int i = 0; i < n; ++i) {
135  local_coordinates[i] = (support[i] - point) / scale;
136  }
137  solver_t solver;
138  ei_matrix_t matrix = getMatrix(local_coordinates).transpose();
139  solver.compute(matrix);
140 
141  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> w(n);
142  for (int j = 0; j < n; j++) {
143  w(j) = weight_(local_coordinates[j]);
144  }
145  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> coefficients = solver.solve(values.cwiseProduct(w));
146  scalar_t residual = (matrix*coefficients - values).squaredNorm();
147  return WLSApproximant<basis_t>(basis_, point, support, scale, coefficients, residual);
148 }
149 
150 } // namespace mm
151 
152 #endif // MEDUSA_BITS_APPROXIMATIONS_WLS_HPP_
mm
Root namespace for the whole library.
Definition: Gaussian.hpp:14
WLSApproximant.hpp
JacobiSVDWrapper.hpp
mm::WLS::getWeights
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getWeights() const
Returns the current support weights. Only initialized after the first call to compute.
mm::operator<<
std::ostream & operator<<(std::ostream &os, const Gaussian< S > &b)
Output basic information about given Gaussian RBF.
Definition: Gaussian.hpp:37
numutils.hpp
mm::WLS::ei_matrix_t
Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > ei_matrix_t
Eigen matrix data type.
Definition: WLS_fwd.hpp:61
assert_msg
#define assert_msg(cond,...)
Assert with better error reporting.
Definition: assert.hpp:75
mm::WLS::getShape
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getShape() const
Returns weights for approximation of function value at point point, knowing values in its support,...
Definition: WLS.hpp:54
mm::WLS::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...
mm::WLS::vector_t
basis_t::vector_t vector_t
Vector data type.
Definition: WLS_fwd.hpp:57
mm::WLS::compute
void compute(const vector_t &point, const std::vector< vector_t > &support)
Computes, decomposes and stores matrix for shape calculation at given point with given neighbours.
Definition: WLS.hpp:19
mm::WLS::getApproximant
WLSApproximant< basis_t > getApproximant(const vector_t &point, const std::vector< vector_t > &support, const Eigen::MatrixBase< Derived > &values) const
Construct a WLS approximant anchored at point with given values in given support nodes.
Definition: WLS.hpp:123
assert.hpp
mm::WLS::weight_
weight_t weight_
Weight function.
Definition: WLS_fwd.hpp:65
mm::WLS::dim
@ dim
Dimensionality of the domain.
Definition: WLS_fwd.hpp:59
mm::WLS::point_
vector_t point_
Saved center point.
Definition: WLS_fwd.hpp:68
mm::WLS::scalar_t
basis_t::scalar_t scalar_t
Numeric scalar data type, e.g. double.
Definition: WLS_fwd.hpp:52
mm::WLSApproximant
Class representing the function that is a WLS approximant using some basis function over some points.
Definition: WLSApproximant_fwd.hpp:27
mm::WLS::getMatrix
ei_matrix_t getMatrix() const
Returns the matrix used in the approximation.
Definition: WLS_fwd.hpp:161
mm::WLS
A class for generating approximations using Weighted Least Squares over local neighborhoods.
Definition: WLS_fwd.hpp:49
WLS_fwd.hpp
Range.hpp
mm::WLS::basis_
basis_t basis_
Basis functions.
Definition: WLS_fwd.hpp:64