|
Medusa
1.1
Coordinate Free Mehless Method implementation
|
|
Go to the documentation of this file. 1 #ifndef MEDUSA_BITS_APPROXIMATIONS_RBFFD_FWD_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_RBFFD_FWD_HPP_
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>>>
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.");
42 typedef RBFType
rbf_t;
78 void compute(
const vector_t& point,
const std::vector<vector_t>& support);
81 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>
getShape()
const;
84 template <
class operator_t>
85 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>
getShape(
const operator_t& op)
const;
97 template <
typename Derived>
99 const Eigen::MatrixBase<Derived>& Lb)
const;
103 Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic>
getMatrix(
104 const std::vector<vector_t>& local_coordinates)
const;
126 Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic>
getMatrix()
const {
150 template <
typename Derived>
152 const std::vector<vector_t>& support,
const Eigen::MatrixBase<Derived>& values)
const;
155 template <
class R,
class V,
class S,
class L>
160 template <
class R,
class V,
class S,
class L>
162 os <<
"RBFFD approximation engine:\n"
163 <<
" dimension: " << e.
dim <<
'\n';
165 os <<
" last point: not used yet\n";
167 os <<
" last point: " << e.
point_ <<
'\n';
169 os <<
" RBF: " << e.
rbf_ <<
'\n'
170 <<
" Augmentation: " << e.
mon_ <<
'\n'
171 <<
" scale: " << S();
177 #endif // MEDUSA_BITS_APPROXIMATIONS_RBFFD_FWD_HPP_
vec_t center() const
Return the center point of the approximation.
Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > getMatrix() const
Returns the augmented collocation matrix.
Root namespace for the whole library.
Scalar scalar_t
Type of the elements, alias of Scalar.
scalar_t scale() const
Return the scaling factor. Only initialized after the first call to compute.
@ dim
Number of elements of this matrix.
solver_t solver_
Linear solver used in the approximation.
static const double NaN
Not-a-number floating point value.
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...
Computes a RBF-FD approximation of given operator over local neighbourhood.
std::ostream & operator<<(std::ostream &os, const Gaussian< S > &b)
Output basic information about given Gaussian RBF.
Range< vec_t > localCoordinates() const
Returns the local (translated and scaled) coordinates of the support nodes.
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.
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getShape() const
Return shape function (stencil weights) for value approximation.
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.
scalar_t scale_
Saved scale.
const Monomials< vec_t > & monomials() const
Return the augmenting monomial basis.
RBFFD(rbf_t rbf, const Monomials< vec_t > &mon)
Construct a RBFFD engine from given RBF with additional monomial augmentation.
rbf_t rbf_
RBF used in the approximation.
RBFType rbf_t
Radial basis function type.
A class representing Monomial basis.
int requiredSupportSize() const
Returns minimal support size required for this approximation to be performed.
RBFFD(rbf_t rbf)
Construct a RBFFD engine from given RBF using no monomial augmentation.
vector_t point_
Saved point.
Monomials< vec_t > mon_
Monomials used in the approximation.
const solver_t & solver() const
Returns the stored decomposition of the augmented collocation matrix.
Class representing a RBF Interpolant over a set of nodes of the form.
friend std::ostream & operator<<(std::ostream &os, const RBFFD< R, V, S, L > &e)
Output basic info about given approximation engine.
vec_t vector_t
Vector type.
vector_t::scalar_t scalar_t
Scalar type.
const rbf_t & rbf() const
Return the RBF used in this approximation.
Range< vec_t > support_
Saved support scaled to local coordinate system.
@ dim
Dimensionality of the function domain.