Medusa  1.1
Coordinate Free Mehless Method implementation
WLS_fwd.hpp
Go to the documentation of this file.
1 #ifndef MEDUSA_BITS_APPROXIMATIONS_WLS_FWD_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_WLS_FWD_HPP_
3 
11 #include <medusa/Config.hpp>
13 #include <Eigen/Core>
14 #include <type_traits>
15 #include <iosfwd>
16 #include "WLSApproximant_fwd.hpp"
17 #include "Operators_fwd.hpp"
18 
19 namespace mm {
20 
21 template <typename scalar_t, int QRPreconditioner>
22 class JacobiSVDWrapper;
23 
24 template <typename vec_t>
25 class NoWeight;
26 
27 class NoScale;
28 
45 template <class basis_t, class weight_t = NoWeight<typename basis_t::vector_t>,
46  class scale_t = NoScale, class solver_t =
47  JacobiSVDWrapper<typename basis_t::scalar_t,
48  Eigen::ColPivHouseholderQRPreconditioner>>
49 class WLS {
50  public:
51  static_assert(std::is_same<typename basis_t::vector_t, typename weight_t::vector_t>::value,
52  "Basis and weight functions must have the same vector type.");
53 
55  typedef typename basis_t::scalar_t scalar_t;
57  typedef typename basis_t::vector_t vector_t;
59  enum { dim = vector_t::dim };
61  typedef typename Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> ei_matrix_t;
62 
63  private:
64  basis_t basis_;
65  weight_t weight_;
66 
70  solver_t solver_;
71 
72  public:
74  explicit WLS(const basis_t& basis) : basis_(basis), scale_(NaN) {}
80  WLS(const basis_t& basis, const weight_t& weight) :
82 
83  private:
85  ei_matrix_t getMatrix(const std::vector<vector_t>& local_coordinates) const;
86 
87  public:
98  void compute(const vector_t& point, const std::vector<vector_t>& support);
99 
107  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> getShape() const;
108 
119  template <class operator_t>
120  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> getShape(const operator_t& op) const;
121 
131  template <typename Derived>
132  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> getShapeFromRhs(
133  const Eigen::MatrixBase<Derived>& Lb) const;
134 
136  int requiredSupportSize() const { return basis_.size(); }
138  const basis_t& basis() const { return basis_; }
140  const weight_t& weight() const { return weight_; }
141 
146  vector_t center() const { return point_; }
147 
149  scalar_t scale() const { return scale_; }
150 
156 
162 
164  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> getWeights() const;
165 
170  const solver_t& solver() const { return solver_; }
171 
177  solver_t& solver() { return solver_; }
178 
189  template <typename Derived>
191  const std::vector<vector_t>& support, const Eigen::MatrixBase<Derived>& values) const;
192 
194  template <class B, class W, class S, class L>
195  friend std::ostream& operator<<(std::ostream& os, const WLS<B, W, S, L>& wls);
196 };
197 
198 } // namespace mm
199 
200 #endif // MEDUSA_BITS_APPROXIMATIONS_WLS_FWD_HPP_
mm::WLS::operator<<
friend std::ostream & operator<<(std::ostream &os, const WLS< B, W, S, L > &wls)
Output basic info about given approximation engine.
mm
Root namespace for the whole library.
Definition: Gaussian.hpp:14
scalar_t
Scalar scalar_t
Type of the elements, alias of Scalar.
Definition: MatrixBaseAddons.hpp:16
mm::WLS::WLS
WLS(const basis_t &basis)
Construct WLS with given basis and default constructed weight.
Definition: WLS_fwd.hpp:74
dim
@ dim
Number of elements of this matrix.
Definition: MatrixBaseAddons.hpp:14
mm::WLS::scale_
scalar_t scale_
Saved scale.
Definition: WLS_fwd.hpp:69
mm::NaN
static const double NaN
Not-a-number floating point value.
Definition: Config.hpp:46
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::WLS::center
vector_t center() const
Return the center point of the approximation.
Definition: WLS_fwd.hpp:146
mm::WLS::solver
const solver_t & solver() const
Read access to the linear solver used for matrix decomposition.
Definition: WLS_fwd.hpp:170
WLSApproximant_fwd.hpp
mm::WLS::solver_
solver_t solver_
Saved decomposition.
Definition: WLS_fwd.hpp:70
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
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::WLS
WLS(const basis_t &basis, const weight_t &weight)
Hold information about the type of approximation we want to have for the domain.
Definition: WLS_fwd.hpp:80
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...
Config.hpp
mm::WLS::vector_t
basis_t::vector_t vector_t
Vector data type.
Definition: WLS_fwd.hpp:57
Range_fwd.hpp
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
Operators_fwd.hpp
mm::WLS::local_coordinates_
Range< vector_t > local_coordinates_
Support points after translation and scaling.
Definition: WLS_fwd.hpp:67
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::weight
const weight_t & weight() const
Returns weight function used in the computation.
Definition: WLS_fwd.hpp:140
mm::WLS::scale
scalar_t scale() const
Return the scaling factor. Only initialized after the first call to compute.
Definition: WLS_fwd.hpp:149
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
mm::WLS::requiredSupportSize
int requiredSupportSize() const
Returns minimal support size required for this approximation to be performed.
Definition: WLS_fwd.hpp:136
mm::WLS::solver
solver_t & solver()
Read-write access to the linear solver used for matrix decomposition.
Definition: WLS_fwd.hpp:177
mm::Range< vector_t >
mm::WLS::localCoordinates
Range< vector_t > localCoordinates() const
Returns the local (translated and scaled) coordinates of the support nodes.
Definition: WLS_fwd.hpp:155
mm::WLS::basis
const basis_t & basis() const
Returns basis used in the computation.
Definition: WLS_fwd.hpp:138
mm::WLS::basis_
basis_t basis_
Basis functions.
Definition: WLS_fwd.hpp:64