#include <WLS_fwd.hpp>
A class for generating approximations using Weighted Least Squares over local neighborhoods.
basis_t | Type of basis used in approximation. Must satisfy the Basis function concept. |
weight_t | Type of weight function used in approximation. Must satisfy the Weight function concept. |
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: Jacobi SVD. |
This class satisfies the Approximation engine concept.
Usage example:
Definition at line 49 of file WLS_fwd.hpp.
Public Member Functions | |
WLS (const basis_t &basis) | |
Construct WLS with given basis and default constructed weight. More... | |
WLS (const basis_t &basis, const weight_t &weight) | |
Hold information about the type of approximation we want to have for the domain. More... | |
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. More... | |
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > | getShape () const |
Returns weights for approximation of function value at point point , knowing values in its support , given by a previous call to compute. More... | |
template<class operator_t > | |
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > | getShape (const operator_t&op) const |
Returns weights for approximation of differentiation operator, given by der at point point , knowing values in its support , given by a previous call to compute. 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. More... | |
int | requiredSupportSize () const |
Returns minimal support size required for this approximation to be performed. More... | |
const basis_t & | basis () const |
Returns basis used in the computation. More... | |
const weight_t & | weight () const |
Returns weight function used in the computation. More... | |
vector_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< vector_t > | localCoordinates () const |
Returns the local (translated and scaled) coordinates of the support nodes. More... | |
ei_matrix_t | getMatrix () const |
Returns the matrix used in the approximation. More... | |
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > | getWeights () const |
Returns the current support weights. Only initialized after the first call to compute. More... | |
const solver_t & | solver () const |
Read access to the linear solver used for matrix decomposition. More... | |
solver_t & | solver () |
Read-write access to the linear solver used for matrix decomposition. More... | |
template<typename Derived > | |
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. More... | |
Public Types | |
enum | { dim = vector_t::dim } |
Store dimension of the domain. More... | |
typedef basis_t::scalar_t | scalar_t |
Numeric scalar data type, e.g. double. More... | |
typedef basis_t::vector_t | vector_t |
Vector data type. More... | |
typedef Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > | ei_matrix_t |
Eigen matrix data type. More... | |
Friends | |
template<class B , class W , class S , class L > | |
std::ostream & | operator<< (std::ostream &os, const WLS< B, W, S, L > &wls) |
Output basic info about given approximation engine. More... | |
Private Member Functions | |
ei_matrix_t | getMatrix (const std::vector< vector_t > &local_coordinates) const |
Returns the collocation matrix. More... | |
Private Attributes | |
basis_t | basis_ |
Basis functions. More... | |
weight_t | weight_ |
Weight function. More... | |
Range< vector_t > | local_coordinates_ |
Support points after translation and scaling. More... | |
vector_t | point_ |
Saved center point. More... | |
scalar_t | scale_ |
Saved scale. More... | |
solver_t | solver_ |
Saved decomposition. More... | |
anonymous enum |
Store dimension of the domain.
Enumerator | |
---|---|
dim | Dimensionality of the domain. |
Definition at line 59 of file WLS_fwd.hpp.
|
inlineexplicit |
Construct WLS with given basis and default constructed weight.
Definition at line 74 of file WLS_fwd.hpp.
|
inline |
Hold information about the type of approximation we want to have for the domain.
basis | Basis to use in the algorithm. |
weight | Weight function to use. |
Definition at line 80 of file WLS_fwd.hpp.
|
inline |
Returns basis used in the computation.
Definition at line 138 of file WLS_fwd.hpp.
|
inline |
Return the center point of the approximation.
Only initialized after the first call to compute.
Definition at line 146 of file WLS_fwd.hpp.
void mm::WLS< basis_t, weight_t, scale_t, solver_t >::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.
Given basis and weight function are used and evaluated in the local coordinate system with origin in point
and scales by scale_. Decomposition given by solver_ is used. This function must be called before any calls to getShape are made for this point configuration.
point | Point around which to prepare for shape computation. |
support | Neighbouring points where function values are expected. |
WLSApproximant< basis_t > mm::WLS< basis_t, weight_t, scale_t, solver_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.
The approximation function \(\hat{u}(x) = \sum_{j=1}^m \alpha_j b_j\left(\frac{x - p}{s}\right)\) minimises
\[ J = \sum_{i=1}^n w\left( \frac{p_i-p}{s}\right)^2 (u_i - \hat{u}(p_i))^2 \]
point | Center point \(p\). |
support | Support values \(p_i\). |
values | Function values \(u_i\) at \(p_i\). |
|
inline |
Returns the matrix used in the approximation.
Only initialized after the first call to compute.
Definition at line 161 of file WLS_fwd.hpp.
|
private |
Eigen::Matrix< typename basis_t::scalar_t, Eigen::Dynamic, 1 > mm::WLS< basis_t, weight_t, scale_t, solver_t >::getShape |
Returns weights for approximation of function value at point point
, knowing values in its support
, given by a previous call to compute.
point
and \(p_i\) are its support nodes. Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> mm::WLS< basis_t, weight_t, scale_t, solver_t >::getShape | ( | const operator_t& | op | ) | const |
Returns weights for approximation of differentiation operator, given by der
at point point
, knowing values in its support
, given by a previous call to compute.
point
, \(p_i\) are its support nodes and \(\mathcal{L}\) is given by der
as \(\mathcal{L} = \dpar{^{\sum_{i=1}^d m_i}}{x_1^{m_1} \cdots \partial x_d^{m_d}}\), where \(m_i\) are elements of der
. op | Operator to be approximated. |
Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> mm::WLS< basis_t, weight_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.
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 at the point where compute was called. |
Eigen::Matrix<scalar_t, Eigen::Dynamic, 1> mm::WLS< basis_t, weight_t, scale_t, solver_t >::getWeights | ( | ) | const |
Returns the current support weights. Only initialized after the first call to compute.
|
inline |
Returns the local (translated and scaled) coordinates of the support nodes.
Only initialized after the first call to compute.
Definition at line 155 of file WLS_fwd.hpp.
|
inline |
Returns minimal support size required for this approximation to be performed.
Definition at line 136 of file WLS_fwd.hpp.
|
inline |
Return the scaling factor. Only initialized after the first call to compute.
Definition at line 149 of file WLS_fwd.hpp.
|
inline |
Read-write access to the linear solver used for matrix decomposition.
Can be used to set additional solver parameters. Only initialized after the first call to compute.
Definition at line 177 of file WLS_fwd.hpp.
|
inline |
Read access to the linear solver used for matrix decomposition.
Only initialized after the first call to compute.
Definition at line 170 of file WLS_fwd.hpp.
|
inline |
Returns weight function used in the computation.
Definition at line 140 of file WLS_fwd.hpp.
|
friend |
Output basic info about given approximation engine.
|
private |
Basis functions.
Definition at line 64 of file WLS_fwd.hpp.
|
private |
Support points after translation and scaling.
Definition at line 67 of file WLS_fwd.hpp.
|
private |
Saved center point.
Definition at line 68 of file WLS_fwd.hpp.
|
private |
Saved scale.
Definition at line 69 of file WLS_fwd.hpp.
|
private |
Saved decomposition.
Definition at line 70 of file WLS_fwd.hpp.
|
private |
Weight function.
Definition at line 65 of file WLS_fwd.hpp.