#include <RBFFD_fwd.hpp>
Computes a RBF-FD approximation of given operator over local neighbourhood.
RBFType | Type of RBF used in approximation. Must satisfy the Radial Basis Function concept. |
vec_t | Vector type used in calculations. |
scale_t | Scale function to be used in approximation. Must satisfy the Scale function concept. |
solver_t | Which 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:
Definition at line 38 of file RBFFD_fwd.hpp.
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_t & | rbf () 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... | |
anonymous enum |
Store dimension of the domain.
Enumerator | |
---|---|
dim | Dimensionality of the function domain. |
Definition at line 46 of file RBFFD_fwd.hpp.
|
inlineexplicit |
Construct a RBFFD engine from given RBF using no monomial augmentation.
Definition at line 59 of file RBFFD_fwd.hpp.
|
inline |
Construct a RBFFD engine from given RBF with additional monomial augmentation.
Definition at line 62 of file RBFFD_fwd.hpp.
|
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.
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.
point | At which point to construct the approximation. |
support | Coordinates of local support points. |
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\).
point | Center point \(p\). |
support | Support values \(p_i\). |
values | Function values \(u_i\) at \(p_i\). |
|
inline |
Returns the augmented collocation matrix.
Only initialized after the first call to compute.
Definition at line 126 of file RBFFD_fwd.hpp.
|
private |
Returns the augmented collocation matrix.
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.
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.
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.
Lb | The 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. |
|
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.
|
inline |
Return the augmenting monomial basis.
Definition at line 68 of file RBFFD_fwd.hpp.
|
inline |
Return the RBF used in this approximation.
Definition at line 65 of file RBFFD_fwd.hpp.
|
inline |
Returns minimal support size required for this approximation to be performed.
Definition at line 71 of file RBFFD_fwd.hpp.
|
inline |
Return the scaling factor. Only initialized after the first call to compute.
Definition at line 114 of file RBFFD_fwd.hpp.
|
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.
|
friend |
Output basic info about given approximation engine.
Definition at line 161 of file RBFFD_fwd.hpp.
|
private |
Monomials used in the approximation.
Definition at line 50 of file RBFFD_fwd.hpp.
|
private |
Saved point.
Definition at line 53 of file RBFFD_fwd.hpp.
|
private |
RBF used in the approximation.
Definition at line 49 of file RBFFD_fwd.hpp.
|
private |
Saved scale.
Definition at line 54 of file RBFFD_fwd.hpp.
|
private |
Linear solver used in the approximation.
Definition at line 51 of file RBFFD_fwd.hpp.
|
private |
Saved support scaled to local coordinate system.
Definition at line 55 of file RBFFD_fwd.hpp.