Medusa  1.1
Coordinate Free Mehless Method implementation
mm::WLS< basis_t, weight_t, scale_t, solver_t > Class Template Reference

#include <WLS_fwd.hpp>

Detailed Description

template<class basis_t, class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
class mm::WLS< basis_t, weight_t, scale_t, solver_t >

A class for generating approximations using Weighted Least Squares over local neighborhoods.

Template Parameters
basis_tType of basis used in approximation. Must satisfy the Basis function concept.
weight_tType of weight function used in approximation. Must satisfy the Weight function concept.
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: Jacobi SVD.

This class satisfies the Approximation engine concept.

Usage example:

WLS<Monomials<Vec2d>, NoWeight<Vec2d>, ScaleToFarthest, Eigen::PartialPivLU<Eigen::MatrixXd>>
appr(Monomials<Vec2d>::tensorBasis(2)); // full type specification with scale and solver
std::cout << appr << std::endl;
WLS<Monomials<Vec2d>> appr2(2); // using default parameters
// Local neighbourhood
double h = 0.1;
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.basis() << std::endl;
std::cout << appr.weight() << 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;
std::cout << appr.getWeights() << 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();
Examples
test/end2end/poisson_explicit.cpp.

Definition at line 49 of file WLS_fwd.hpp.

+ Collaboration diagram for mm::WLS< basis_t, weight_t, scale_t, solver_t >:

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_tlocalCoordinates () 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_tlocal_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...
 

Member Enumeration Documentation

◆ anonymous enum

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
anonymous enum

Store dimension of the domain.

Enumerator
dim 

Dimensionality of the domain.

Definition at line 59 of file WLS_fwd.hpp.

Constructor & Destructor Documentation

◆ WLS() [1/2]

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
mm::WLS< basis_t, weight_t, scale_t, solver_t >::WLS ( const basis_t &  basis)
inlineexplicit

Construct WLS with given basis and default constructed weight.

Definition at line 74 of file WLS_fwd.hpp.

◆ WLS() [2/2]

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
mm::WLS< basis_t, weight_t, scale_t, solver_t >::WLS ( const basis_t &  basis,
const weight_t &  weight 
)
inline

Hold information about the type of approximation we want to have for the domain.

Parameters
basisBasis to use in the algorithm.
weightWeight function to use.

Definition at line 80 of file WLS_fwd.hpp.

Member Function Documentation

◆ basis()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
const basis_t& mm::WLS< basis_t, weight_t, scale_t, solver_t >::basis ( ) const
inline

Returns basis used in the computation.

Definition at line 138 of file WLS_fwd.hpp.

◆ center()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
vector_t mm::WLS< basis_t, weight_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 146 of file WLS_fwd.hpp.

◆ compute()

template<class basis_t , class weight_t , class scale_t , class solver_t >
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.

Parameters
pointPoint around which to prepare for shape computation.
supportNeighbouring points where function values are expected.

Definition at line 19 of file WLS.hpp.

◆ getApproximant()

template<class basis_t , class weight_t , class scale_t , class solver_t >
template<typename Derived >
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 \]

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

Definition at line 123 of file WLS.hpp.

◆ getMatrix() [1/2]

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
ei_matrix_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::getMatrix ( ) const
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.

◆ getMatrix() [2/2]

template<class basis_t , class weight_t , class scale_t , class solver_t >
WLS< basis_t, weight_t, scale_t, solver_t >::ei_matrix_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::getMatrix ( const std::vector< vector_t > &  local_coordinates) const
private

Returns the collocation matrix.

Definition at line 37 of file WLS.hpp.

◆ getShape() [1/2]

template<class basis_t , class weight_t , class scale_t , class solver_t >
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.

Returns
Vector of nodal weights \(\chi\), such that \(u(p) \approx \sum_{i=1}^n \chi_i u(p_i)\) where \(p\) represents point and \(p_i\) are its support nodes.

Definition at line 54 of file WLS.hpp.

◆ getShape() [2/2]

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
template<class operator_t >
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.

Returns
Vector of nodal weights \(\chi\), such that \((\mathcal{L}u)(p) \approx \sum_{i=1}^n \chi_i u(p_i)\) where \(p\) represents 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.
Parameters
opOperator to be approximated.

◆ getShapeFromRhs()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
template<typename Derived >
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.

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

◆ getWeights()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
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.

◆ localCoordinates()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
Range<vector_t> mm::WLS< basis_t, weight_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 155 of file WLS_fwd.hpp.

◆ requiredSupportSize()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
int mm::WLS< basis_t, weight_t, scale_t, solver_t >::requiredSupportSize ( ) const
inline

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

Definition at line 136 of file WLS_fwd.hpp.

◆ scale()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
scalar_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::scale ( ) const
inline

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

Definition at line 149 of file WLS_fwd.hpp.

◆ solver() [1/2]

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
solver_t& mm::WLS< basis_t, weight_t, scale_t, solver_t >::solver ( )
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.

◆ solver() [2/2]

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
const solver_t& mm::WLS< basis_t, weight_t, scale_t, solver_t >::solver ( ) const
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.

◆ weight()

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
const weight_t& mm::WLS< basis_t, weight_t, scale_t, solver_t >::weight ( ) const
inline

Returns weight function used in the computation.

Definition at line 140 of file WLS_fwd.hpp.

Friends And Related Function Documentation

◆ operator<<

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
template<class B , class W , class S , class L >
std::ostream& operator<< ( std::ostream &  os,
const WLS< B, W, S, L > &  wls 
)
friend

Output basic info about given approximation engine.

Member Data Documentation

◆ basis_

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
basis_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::basis_
private

Basis functions.

Definition at line 64 of file WLS_fwd.hpp.

◆ local_coordinates_

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
Range<vector_t> mm::WLS< basis_t, weight_t, scale_t, solver_t >::local_coordinates_
private

Support points after translation and scaling.

Definition at line 67 of file WLS_fwd.hpp.

◆ point_

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
vector_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::point_
private

Saved center point.

Definition at line 68 of file WLS_fwd.hpp.

◆ scale_

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
scalar_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::scale_
private

Saved scale.

Definition at line 69 of file WLS_fwd.hpp.

◆ solver_

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
solver_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::solver_
private

Saved decomposition.

Definition at line 70 of file WLS_fwd.hpp.

◆ weight_

template<class basis_t , class weight_t = NoWeight<typename basis_t::vector_t>, class scale_t = NoScale, class solver_t = JacobiSVDWrapper<typename basis_t::scalar_t, Eigen::ColPivHouseholderQRPreconditioner>>
weight_t mm::WLS< basis_t, weight_t, scale_t, solver_t >::weight_
private

Weight function.

Definition at line 65 of file WLS_fwd.hpp.


The documentation for this class was generated from the following files:
mm::WLS::solver
const solver_t & solver() const
Read access to the linear solver used for matrix decomposition.
Definition: WLS_fwd.hpp:170
mm::Monomials::tensorBasis
static Monomials< vec_t > tensorBasis(int order)
Construct a tensor basis of monomials ranging from 0 up to order (inclusive) in each dimension.
Definition: Monomials.hpp:41