|
Medusa
1.1
Coordinate Free Mehless Method implementation
|
|
Go to the documentation of this file. 1 #ifndef MEDUSA_BITS_APPROXIMATIONS_WLS_FWD_HPP_
2 #define MEDUSA_BITS_APPROXIMATIONS_WLS_FWD_HPP_
14 #include <type_traits>
21 template <
typename scalar_t,
int QRPreconditioner>
22 class JacobiSVDWrapper;
24 template <
typename vec_t>
45 template <
class basis_t,
class weight_t = NoWeight<
typename basis_t::vector_t>,
46 class scale_t = NoScale,
class solver_t =
47 JacobiSVDWrapper<
typename basis_t::scalar_t,
48 Eigen::ColPivHouseholderQRPreconditioner>>
51 static_assert(std::is_same<typename basis_t::vector_t, typename weight_t::vector_t>::value,
52 "Basis and weight functions must have the same vector type.");
61 typedef typename Eigen::Matrix<scalar_t, Eigen::Dynamic, Eigen::Dynamic>
ei_matrix_t;
98 void compute(
const vector_t& point,
const std::vector<vector_t>& support);
107 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>
getShape()
const;
119 template <
class operator_t>
120 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>
getShape(
const operator_t& op)
const;
131 template <
typename Derived>
133 const Eigen::MatrixBase<Derived>& Lb)
const;
164 Eigen::Matrix<scalar_t, Eigen::Dynamic, 1>
getWeights()
const;
189 template <
typename Derived>
191 const std::vector<vector_t>& support,
const Eigen::MatrixBase<Derived>& values)
const;
194 template <
class B,
class W,
class S,
class L>
200 #endif // MEDUSA_BITS_APPROXIMATIONS_WLS_FWD_HPP_
friend std::ostream & operator<<(std::ostream &os, const WLS< B, W, S, L > &wls)
Output basic info about given approximation engine.
Root namespace for the whole library.
Scalar scalar_t
Type of the elements, alias of Scalar.
WLS(const basis_t &basis)
Construct WLS with given basis and default constructed weight.
@ dim
Number of elements of this matrix.
scalar_t scale_
Saved scale.
static const double NaN
Not-a-number floating point value.
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getWeights() const
Returns the current support weights. Only initialized after the first call to compute.
vector_t center() const
Return the center point of the approximation.
const solver_t & solver() const
Read access to the linear solver used for matrix decomposition.
solver_t solver_
Saved decomposition.
Eigen::Matrix< scalar_t, Eigen::Dynamic, Eigen::Dynamic > ei_matrix_t
Eigen matrix data type.
Eigen::Matrix< scalar_t, Eigen::Dynamic, 1 > getShape() const
Returns weights for approximation of function value at point point, knowing values in its support,...
WLS(const basis_t &basis, const weight_t &weight)
Hold information about the type of approximation we want to have for the domain.
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...
basis_t::vector_t vector_t
Vector data type.
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.
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.
Range< vector_t > local_coordinates_
Support points after translation and scaling.
weight_t weight_
Weight function.
@ dim
Dimensionality of the domain.
vector_t point_
Saved center point.
basis_t::scalar_t scalar_t
Numeric scalar data type, e.g. double.
Class representing the function that is a WLS approximant using some basis function over some points.
const weight_t & weight() const
Returns weight function used in the computation.
scalar_t scale() const
Return the scaling factor. Only initialized after the first call to compute.
ei_matrix_t getMatrix() const
Returns the matrix used in the approximation.
A class for generating approximations using Weighted Least Squares over local neighborhoods.
int requiredSupportSize() const
Returns minimal support size required for this approximation to be performed.
solver_t & solver()
Read-write access to the linear solver used for matrix decomposition.
Range< vector_t > localCoordinates() const
Returns the local (translated and scaled) coordinates of the support nodes.
const basis_t & basis() const
Returns basis used in the computation.
basis_t basis_
Basis functions.