Medusa  1.1
Coordinate Free Mehless Method implementation
mm::RBFFD< RBFType, vec_t, scale_t, solver_t > Class Template Reference

#include <RBFFD_fwd.hpp>

Detailed Description

template<class RBFType, class vec_t, class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
class mm::RBFFD< RBFType, vec_t, scale_t, solver_t >

Computes a RBF-FD approximation of given operator over local neighbourhood.

Template Parameters
RBFTypeType of RBF used in approximation. Must satisfy the Radial Basis Function concept.
vec_tVector type used in calculations.
scale_tScale function to be used in approximation. Must satisfy the Scale function concept.
solver_tWhich linear solver to use in weight computations. Must satisfy the Linear solver concept. Default: LU with partial pivoting.

This class satisfies the Approximation engine concept.

Usage example:

double h = 0.1123;
double s = 1.5;
Gaussian<double> g(s);
// Gaussian RBF's augmented with a constant
RBFFD<Gaussian<double>, Vec2d> appr(g, Monomials<Vec2d>(0));
std::cout << appr << std::endl;
// Local neighbourhood
Range<Vec2d> support = {{0, 0}, {0, h}, {h, 0}, {0, -h}, {-h, 0},
{-h, h}, {-h, -h}, {h, -h}, {h, h}};
// compute approximations at point `{0.0, 0.0}`.
appr.compute({0.0, 0.0}, support);
// Get info about the computation.
std::cout << appr.rbf() << std::endl;
std::cout << appr.monomials() << std::endl;
std::cout << appr.center() << std::endl;
std::cout << appr.scale() << std::endl;
std::cout << appr.localCoordinates() << std::endl;
std::cout << appr.getMatrix() << std::endl;
Eigen::PartialPivLU<Eigen::MatrixXd> solver = appr.solver();
// Get shape (stencil weights) for approximation of Laplacian.
Eigen::VectorXd shape = appr.getShape(Lap<2>());
// Get shape (stencil weights) for approximation of value
shape = appr.getShape();

Definition at line 38 of file RBFFD_fwd.hpp.

+ Collaboration diagram for mm::RBFFD< RBFType, vec_t, scale_t, solver_t >:

Public Member Functions

 RBFFD (rbf_t rbf)
 Construct a RBFFD engine from given RBF using no monomial augmentation. More...
 
 RBFFD (rbf_t rbf, const Monomials< vec_t > &mon)
 Construct a RBFFD engine from given RBF with additional monomial augmentation. More...
 
const rbf_trbf () const
 Return the RBF used in this approximation. More...
 
const Monomials< vec_t > & monomials () const
 Return the augmenting monomial basis. More...
 
int requiredSupportSize () const
 Returns minimal support size required for this approximation to be performed. More...
 
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. More...
 
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getShape () const
 Return shape function (stencil weights) for value approximation. More...
 
template<class operator_t >
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getShape (const operator_t&op) const
 Return shape function (stencil weights) for given derivative operator. More...
 
template<typename Derived >
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 values of an operator applied to your basis functions (both RBFs and monomials). More...
 
vec_t center () const
 Return the center point of the approximation. More...
 
scalar_t scale () const
 Return the scaling factor. Only initialized after the first call to compute. More...
 
Range< vec_t > localCoordinates () const
 Returns the local (translated and scaled) coordinates of the support nodes. More...
 
Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > getMatrix () const
 Returns the augmented collocation matrix. More...
 
const solver_t & solver () const
 Returns the stored decomposition of the augmented collocation matrix. More...
 
template<typename Derived >
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. More...
 

Public Types

enum  { dim = vec_t::dim }
 Store dimension of the domain. More...
 
typedef RBFType rbf_t
 Radial basis function type. More...
 
typedef vec_t vector_t
 Vector type. More...
 
typedef vector_t::scalar_t scalar_t
 Scalar type. More...
 

Friends

template<class R , class V , class S , class L >
std::ostream & operator<< (std::ostream &os, const RBFFD< R, V, S, L > &e)
 Output basic info about given approximation engine. More...
 

Private Member Functions

Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > getMatrix (const std::vector< vector_t > &local_coordinates) const
 Returns the augmented collocation matrix. More...
 

Private Attributes

rbf_t rbf_
 RBF used in the approximation. More...
 
Monomials< vec_t > mon_
 Monomials used in the approximation. More...
 
solver_t solver_
 Linear solver used in the approximation. More...
 
vector_t point_
 Saved point. More...
 
scalar_t scale_
 Saved scale. More...
 
Range< vec_t > support_
 Saved support scaled to local coordinate system. More...
 

Member Enumeration Documentation

◆ anonymous enum

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
anonymous enum

Store dimension of the domain.

Enumerator
dim 

Dimensionality of the function domain.

Definition at line 46 of file RBFFD_fwd.hpp.

Constructor & Destructor Documentation

◆ RBFFD() [1/2]

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::RBFFD ( rbf_t  rbf)
inlineexplicit

Construct a RBFFD engine from given RBF using no monomial augmentation.

Definition at line 59 of file RBFFD_fwd.hpp.

◆ RBFFD() [2/2]

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::RBFFD ( rbf_t  rbf,
const Monomials< vec_t > &  mon 
)
inline

Construct a RBFFD engine from given RBF with additional monomial augmentation.

Definition at line 62 of file RBFFD_fwd.hpp.

Member Function Documentation

◆ center()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
vec_t mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::center ( ) const
inline

Return the center point of the approximation.

Only initialized after the first call to compute.

Definition at line 111 of file RBFFD_fwd.hpp.

◆ compute()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
void mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::compute ( const vector_t point,
const std::vector< vector_t > &  support 
)

Setup this engine to approximate operators at given point over given local neighbourhood.

Parameters
pointAt which point to construct the approximation.
supportCoordinates of local support points.

◆ getApproximant()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
template<typename Derived >
RBFInterpolant<rbf_t, vec_t> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::getApproximant ( const vector_t point,
const std::vector< vector_t > &  support,
const Eigen::MatrixBase< Derived > &  values 
) const

Construct a RBF interpolant of form.

\[ u(p) = \sum_{i=1}^n \alpha_i \phi\left(\left\|\frac{p - p_i}{s}\right\|\right) + \sum_{i=1}^s \beta_i q_i\left(\frac{p - p_c}{s} \right), \]

where \(\phi\) is the RBF and \(q_i\) are monomials. Additionally, \(s\) constraints are imposed for \(\beta_i\), \(\sum_{i=1}^n q_j(p_i) = 0\), for \(j = 1, \ldots, s\).

Parameters
pointCenter point \(p\).
supportSupport values \(p_i\).
valuesFunction values \(u_i\) at \(p_i\).
Returns
An object representing the interpolant.

◆ getMatrix() [1/2]

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::getMatrix ( ) const
inline

Returns the augmented collocation matrix.

Only initialized after the first call to compute.

Definition at line 126 of file RBFFD_fwd.hpp.

◆ getMatrix() [2/2]

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::getMatrix ( const std::vector< vector_t > &  local_coordinates) const
private

Returns the augmented collocation matrix.

◆ getShape() [1/2]

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::getShape ( ) const

Return shape function (stencil weights) for value approximation.

◆ getShape() [2/2]

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
template<class operator_t >
Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::getShape ( const operator_t&  op) const

Return shape function (stencil weights) for given derivative operator.

◆ getShapeFromRhs()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
template<typename Derived >
Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::getShapeFromRhs ( const Eigen::MatrixBase< Derived > &  Lb) const

Returns weights for approximation with a given right hand side, which is expected to contains the values of an operator applied to your basis functions (both RBFs and monomials).

This method is useful if values of the operator are computed differently to the procedure used by getShape.

Parameters
LbThe right hand side of the approximation system representing an values of an operator applied to basis functions (both RBFs and monomials) at the point where compute was called.
See also
getShape

◆ localCoordinates()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
Range<vec_t> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::localCoordinates ( ) const
inline

Returns the local (translated and scaled) coordinates of the support nodes.

Only initialized after the first call to compute.

Definition at line 120 of file RBFFD_fwd.hpp.

◆ monomials()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
const Monomials<vec_t>& mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::monomials ( ) const
inline

Return the augmenting monomial basis.

Definition at line 68 of file RBFFD_fwd.hpp.

◆ rbf()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
const rbf_t& mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::rbf ( ) const
inline

Return the RBF used in this approximation.

Definition at line 65 of file RBFFD_fwd.hpp.

◆ requiredSupportSize()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
int mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::requiredSupportSize ( ) const
inline

Returns minimal support size required for this approximation to be performed.

Definition at line 71 of file RBFFD_fwd.hpp.

◆ scale()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
scalar_t mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::scale ( ) const
inline

Return the scaling factor. Only initialized after the first call to compute.

Definition at line 114 of file RBFFD_fwd.hpp.

◆ solver()

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
const solver_t& mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::solver ( ) const
inline

Returns the stored decomposition of the augmented collocation matrix.

Only initialized after the first call to compute.

Definition at line 134 of file RBFFD_fwd.hpp.

Friends And Related Function Documentation

◆ operator<<

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
template<class R , class V , class S , class L >
std::ostream& operator<< ( std::ostream &  os,
const RBFFD< R, V, S, L > &  e 
)
friend

Output basic info about given approximation engine.

Definition at line 161 of file RBFFD_fwd.hpp.

Member Data Documentation

◆ mon_

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
Monomials<vec_t> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::mon_
private

Monomials used in the approximation.

Definition at line 50 of file RBFFD_fwd.hpp.

◆ point_

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
vector_t mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::point_
private

Saved point.

Definition at line 53 of file RBFFD_fwd.hpp.

◆ rbf_

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
rbf_t mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::rbf_
private

RBF used in the approximation.

Definition at line 49 of file RBFFD_fwd.hpp.

◆ scale_

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
scalar_t mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::scale_
private

Saved scale.

Definition at line 54 of file RBFFD_fwd.hpp.

◆ solver_

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
solver_t mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::solver_
private

Linear solver used in the approximation.

Definition at line 51 of file RBFFD_fwd.hpp.

◆ support_

template<class RBFType , class vec_t , class scale_t = NoScale, class solver_t = Eigen::PartialPivLU<Eigen::Matrix<typename vec_t::scalar_t, Eigen::Dynamic, Eigen::Dynamic>>>
Range<vec_t> mm::RBFFD< RBFType, vec_t, scale_t, solver_t >::support_
private

Saved support scaled to local coordinate system.

Definition at line 55 of file RBFFD_fwd.hpp.


The documentation for this class was generated from the following file:
mm::RBFFD::solver
const solver_t & solver() const
Returns the stored decomposition of the augmented collocation matrix.
Definition: RBFFD_fwd.hpp:134
mm::Vec2d
Vec< double, 2 > Vec2d
Convenience typedef for 2d vector of doubles.
Definition: Vec_fwd.hpp:34