Medusa  1.1
Coordinate Free Mehless Method implementation
RBFFD_fwd.hpp
Go to the documentation of this file.
1 #ifndef MEDUSA_BITS_APPROXIMATIONS_RBFFD_FWD_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_RBFFD_FWD_HPP_
3 
11 #include <medusa/Config.hpp>
13 #include <Eigen/Core>
14 #include "RBFBasis_fwd.hpp"
15 #include "Monomials_fwd.hpp"
16 #include "ScaleFunction.hpp"
17 #include "RBFInterpolant_fwd.hpp"
18 
19 namespace mm {
20 
36 template <class RBFType, class vec_t, class scale_t = NoScale, class solver_t =
37  Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
38 class RBFFD {
39  static_assert(std::is_same<typename vec_t::scalar_t, typename RBFType::scalar_t>::value,
40  "Basis and underlying RBF must have the same scalar type.");
41  public:
42  typedef RBFType rbf_t;
43  typedef vec_t vector_t;
44  typedef typename vector_t::scalar_t scalar_t;
45  enum { dim = vec_t::dim };
47 
48  private:
51  solver_t solver_;
52 
56 
57  public:
59  explicit RBFFD(rbf_t rbf) : rbf_(rbf), mon_(), point_(NaN) {}
60 
62  RBFFD(rbf_t rbf, const Monomials<vec_t>& mon) : rbf_(rbf), mon_(mon), point_(NaN) {}
63 
65  const rbf_t& rbf() const { return rbf_; }
66 
68  const Monomials<vec_t>& monomials() const { return mon_; }
69 
71  int requiredSupportSize() const { return mon_.size(); }
72 
78  void compute(const vector_t& point, const std::vector<vector_t>& support);
79 
81  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> getShape() const;
82 
84  template <class operator_t>
85  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> getShape(const operator_t& op) const;
86 
97  template <typename Derived>
98  Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> getShapeFromRhs(
99  const Eigen::MatrixBase<Derived>& Lb) const;
100 
101  private:
103  Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> getMatrix(
104  const std::vector<vector_t>& local_coordinates) const;
105 
106  public:
111  vec_t center() const { return point_; }
112 
114  scalar_t scale() const { return scale_; }
115 
121 
126  Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> getMatrix() const {
127  return getMatrix(support_);
128  }
129 
134  const solver_t& solver() const { return solver_; }
135 
150  template <typename Derived>
152  const std::vector<vector_t>& support, const Eigen::MatrixBase<Derived>& values) const;
153 
155  template <class R, class V, class S, class L>
156  friend std::ostream& operator<<(std::ostream& os, const RBFFD<R, V, S, L>& e);
157 };
158 
160 template <class R, class V, class S, class L>
161 std::ostream& operator<<(std::ostream& os, const RBFFD<R, V, S, L>& e) {
162  os << "RBFFD approximation engine:\n"
163  << " dimension: " << e.dim << '\n';
164  if (e.point_[0] != e.point_[0]) {
165  os << " last point: not used yet\n";
166  } else {
167  os << " last point: " << e.point_ << '\n';
168  }
169  os << " RBF: " << e.rbf_ << '\n'
170  << " Augmentation: " << e.mon_ << '\n'
171  << " scale: " << S();
172  return os;
173 }
174 
175 } // namespace mm
176 
177 #endif // MEDUSA_BITS_APPROXIMATIONS_RBFFD_FWD_HPP_
mm::RBFFD::center
vec_t center() const
Return the center point of the approximation.
Definition: RBFFD_fwd.hpp:111
RBFBasis_fwd.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
scalar_t
Scalar scalar_t
Type of the elements, alias of Scalar.
Definition: MatrixBaseAddons.hpp:16
mm::RBFFD::scale
scalar_t scale() const
Return the scaling factor. Only initialized after the first call to compute.
Definition: RBFFD_fwd.hpp:114
dim
@ dim
Number of elements of this matrix.
Definition: MatrixBaseAddons.hpp:14
mm::RBFFD::solver_
solver_t solver_
Linear solver used in the approximation.
Definition: RBFFD_fwd.hpp:51
mm::NaN
static const double NaN
Not-a-number floating point value.
Definition: Config.hpp:46
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...
mm::RBFFD
Computes a RBF-FD approximation of given operator over local neighbourhood.
Definition: RBFFD_fwd.hpp:38
mm::operator<<
std::ostream & operator<<(std::ostream &os, const Gaussian< S > &b)
Output basic information about given Gaussian RBF.
Definition: Gaussian.hpp:37
mm::RBFFD::localCoordinates
Range< vec_t > localCoordinates() const
Returns the local (translated and scaled) coordinates of the support nodes.
Definition: RBFFD_fwd.hpp:120
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.
ScaleFunction.hpp
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.
mm::RBFFD::scale_
scalar_t scale_
Saved scale.
Definition: RBFFD_fwd.hpp:54
mm::RBFFD::monomials
const Monomials< vec_t > & monomials() const
Return the augmenting monomial basis.
Definition: RBFFD_fwd.hpp:68
Config.hpp
mm::RBFFD::RBFFD
RBFFD(rbf_t rbf, const Monomials< vec_t > &mon)
Construct a RBFFD engine from given RBF with additional monomial augmentation.
Definition: RBFFD_fwd.hpp:62
mm::RBFFD::rbf_
rbf_t rbf_
RBF used in the approximation.
Definition: RBFFD_fwd.hpp:49
mm::RBFFD::rbf_t
RBFType rbf_t
Radial basis function type.
Definition: RBFFD_fwd.hpp:40
mm::Monomials
A class representing Monomial basis.
Definition: Monomials_fwd.hpp:38
Range_fwd.hpp
mm::RBFFD::requiredSupportSize
int requiredSupportSize() const
Returns minimal support size required for this approximation to be performed.
Definition: RBFFD_fwd.hpp:71
RBFInterpolant_fwd.hpp
mm::RBFFD::RBFFD
RBFFD(rbf_t rbf)
Construct a RBFFD engine from given RBF using no monomial augmentation.
Definition: RBFFD_fwd.hpp:59
mm::RBFFD::point_
vector_t point_
Saved point.
Definition: RBFFD_fwd.hpp:53
mm::RBFFD::mon_
Monomials< vec_t > mon_
Monomials used in the approximation.
Definition: RBFFD_fwd.hpp:50
mm::RBFFD::solver
const solver_t & solver() const
Returns the stored decomposition of the augmented collocation matrix.
Definition: RBFFD_fwd.hpp:134
mm::RBFInterpolant
Class representing a RBF Interpolant over a set of nodes of the form.
Definition: RBFInterpolant_fwd.hpp:34
mm::RBFFD::operator<<
friend std::ostream & operator<<(std::ostream &os, const RBFFD< R, V, S, L > &e)
Output basic info about given approximation engine.
Definition: RBFFD_fwd.hpp:161
mm::RBFFD::vector_t
vec_t vector_t
Vector type.
Definition: RBFFD_fwd.hpp:43
Monomials_fwd.hpp
mm::RBFFD::scalar_t
vector_t::scalar_t scalar_t
Scalar type.
Definition: RBFFD_fwd.hpp:44
mm::RBFFD::rbf
const rbf_t & rbf() const
Return the RBF used in this approximation.
Definition: RBFFD_fwd.hpp:65
mm::RBFFD::support_
Range< vec_t > support_
Saved support scaled to local coordinate system.
Definition: RBFFD_fwd.hpp:55
mm::RBFFD::dim
@ dim
Dimensionality of the function domain.
Definition: RBFFD_fwd.hpp:46
mm::Range< vec_t >