Implementation of the Reduced Heteroscedastic Errors In Variables surface fitting algorithm [2].
More...
|
using | x_t = Eigen::Matrix< double, k, 1 > |
| Input vector type .
|
|
using | xs_t = Eigen::Matrix< double, k, -1 > |
| Container type for the input data array.
|
|
using | u_t = Eigen::Matrix< double, l, 1 > |
| Transformed input vector type .
|
|
using | P_t = Eigen::Matrix< double, k, k > |
| Covariance type of the input vector .
|
|
using | Ps_t = std::vector< P_t > |
| Container type for covariances P of the input data array.
|
|
using | z_t = Eigen::Matrix< double, lr, 1 > |
| Type of a reduced transformed input vector .
|
|
using | zs_t = Eigen::Matrix< double, lr, -1 > |
| Container type for an array of the reduced transformed input vectors z .
|
|
using | f_z_t = typename std::function< zs_t(const xs_t &)> |
| Function signature of the mapping function.
|
|
using | dzdx_t = Eigen::Matrix< double, lr, k > |
| Type of the jacobian matrix , evaluated at .
|
|
using | dzdxs_t = std::vector< dzdx_t > |
| Contained type for an array of the jacobian matrices.
|
|
using | f_dzdx_t = typename std::function< dzdx_t(const xs_t &)> |
| Function signature of the jacobian .
|
|
using | theta_t = Eigen::Matrix< double, l, 1 > |
| Parameter vector type .
|
|
using | eta_t = z_t |
| Reduced parameter vector type .
|
|
|
| RHEIV () |
| Convenience default constructor. More...
|
|
| RHEIV (const f_z_t &f_z, const f_dzdx_t &f_dzdx, const double min_dtheta=1e-15, const unsigned max_its=100) |
| The main constructor. More...
|
|
template<typename time_t > |
| RHEIV (const f_z_t &f_z, const f_dzdx_t &f_dzdx, const double min_dtheta=1e-15, const unsigned max_its=100, const time_t &timeout=std::chrono::duration_cast< time_t >(ms_t::zero()), const int debug_nth_it=-1) |
| The main constructor. More...
|
|
| RHEIV (const f_z_t &f_z, const dzdx_t &dzdx, const double min_dtheta=1e-15, const unsigned max_its=100) |
| A convenience constructor constructor. More...
|
|
template<typename time_t > |
| RHEIV (const f_z_t &f_z, const dzdx_t &dzdx, const double min_dtheta=1e-15, const unsigned max_its=100, const time_t &timeout=std::chrono::duration_cast< time_t >(ms_t::zero()), const int debug_nth_it=-1) |
| A convenience constructor constructor. More...
|
|
theta_t | fit (const xs_t &xs, const Ps_t &Ps) |
| Fit the defined model to the provided data. More...
|
|
template<class T_it1 , class T_it2 > |
theta_t | fit (const T_it1 &xs_begin, const T_it1 &xs_end, const T_it2 &Ps_begin, const T_it2 &Ps_end) |
| Fit the defined model to the provided data. More...
|
|
theta_t | fit_ALS (const xs_t &xs) |
| Fit the defined model to the provided data using Algebraic Least Squares (not RHEIV). More...
|
|
theta_t | get_last_estimate () const |
| Returns the last valid estimate of . More...
|
|
theta_t | get_ALS_estimate () const |
| Returns the Algebraic Least Squares estimate of . More...
|
|
template<int n_states, int n_params>
class mrs_lib::RHEIV< n_states, n_params >
Implementation of the Reduced Heteroscedastic Errors In Variables surface fitting algorithm [2].
This class estimates a vector of parameters of a model, which is defined by an equation , where is a data sample vector and is a vector, obtained by transforming in a problem-dependent manner.
Note that must fulfill two conditions:
- each element of is a quadratic form of a vector , and
- the last element of is equal to one.
Such model can be used to describe e.g. various surfaces, such as planes or conics, which makes this class a useful tool for computer vision. The model is described in the context of the RHEIV class by the length of the vector , the length of the vector , and the Jacobian matrix of , where is defined by the formula .
For more information, see [2].
- Examples
- rheiv/example.cpp.
template<int n_states, int n_params>
template<typename time_t >
mrs_lib::RHEIV< n_states, n_params >::RHEIV |
( |
const f_z_t & |
f_z, |
|
|
const dzdx_t & |
dzdx, |
|
|
const double |
min_dtheta = 1e-15 , |
|
|
const unsigned |
max_its = 100 , |
|
|
const time_t & |
timeout = std::chrono::duration_cast<time_t>(ms_t::zero()) , |
|
|
const int |
debug_nth_it = -1 |
|
) |
| |
|
inline |
template<int n_states, int n_params>
template<class T_it1 , class T_it2 >
theta_t mrs_lib::RHEIV< n_states, n_params >::fit |
( |
const T_it1 & |
xs_begin, |
|
|
const T_it1 & |
xs_end, |
|
|
const T_it2 & |
Ps_begin, |
|
|
const T_it2 & |
Ps_end |
|
) |
| |
|
inline |