|
| FIR_filt () |
|
| ~FIR_filt () |
|
void | clear (void) |
|
void | set_coeffs (const arma::Mat< T2 > &_b) |
|
void | set_coeffs (const arma::Col< T2 > &_b_col) |
|
arma::Col< T2 > | get_coeffs () |
|
void | update_coeffs (const arma::Mat< T2 > &_b) |
|
T3 | operator() (const T1 &in) |
|
arma::Mat< T3 > | filter (const arma::Mat< T1 > &in) |
|
arma::Col< T3 > | filter (const arma::Col< T1 > &in) |
|
void | setup_lms (const arma::uword _N, const double _mu, const arma::uword _L=1) |
|
void | lms_adapt (const T3 _err) |
|
void | setup_nlms (const arma::uword _N, const double _mu, const T2 _c, const arma::uword _L=1) |
|
void | nlms_adapt (const T3 _err) |
|
void | setup_newt (const arma::uword _N, const double _mu, const T2 _c, const arma::uword _L=1) |
|
void | newt_adapt (const T3 _err) |
|
void | setup_rls (const arma::uword _N, const double _lmd, const double _P0) |
|
void | rls_adapt (const T3 _err) |
|
void | setup_kalman (const arma::uword _N, const double _P0, const double _Q0, const double _R0) |
|
void | kalman_adapt (const T3 _err) |
|
double | get_step_size (void) |
|
arma::Mat< T1 > | get_P (void) |
|
arma::Mat< T1 > | get_K (void) |
|
void | set_step_size (const double _mu) |
|
void | adapt_enable (void) |
|
void | adapt_disble (void) |
|
template<class T1, class T2, class T3>
class sp::FIR_filt< T1, T2, T3 >
Implements FIR/MA filter functions as
\[ y(n) = \sum_{k=0}^{M-1}{b_kx(n-k)}=b_0x(n)+b_1x(n-1)+...+b_{M-1}x(n-(M-1))\]
where M is the number of taps in the FIR filter. The filter order is M-1. Adaptive update of filter is possible with LMS or NLMS algorithms
- Examples:
- adaptive_filter.cpp, and fir_iir.cpp.
Definition at line 20 of file filter.h.
template<class T1, class T2, class T3>
void sp::FIR_filt< T1, T2, T3 >::kalman_adapt |
( |
const T3 |
_err | ) |
|
|
inline |
Kalman Filter update function.
The Kalman filter is updated as
\( \mathbf{b(n)} = \mathbf{b(n-1)}+\mathbf{Kx(n)}\)
where
\( \mathbf{K} =\frac{\mathbf{Px}}{R+\mathbf{x^TPx}} \)
and
\( \mathbf{P^+} =\mathbf{P^-+xP^-x^T }+Q \)
- Parameters
-
Definition at line 409 of file filter.h.
Referenced by main().
template<class T1, class T2, class T3>
LMS Filter update function.
The LMS filter is updated as
\( \mathbf{b(n)} = \mathbf{b(n-1)}+2\mu\mathbf{x(n)}err(n) \)
where
\( err(n) = d(n)-\mathbf{b(n-1)^Tx(n)} \)
- Parameters
-
Definition at line 172 of file filter.h.
template<class T1, class T2, class T3>
void sp::FIR_filt< T1, T2, T3 >::newt_adapt |
( |
const T3 |
_err | ) |
|
|
inline |
LMS-Newton Filter update function.
The LMS-Newton filter is updated as
\( \mathbf{b(n)} = \mathbf{b(n-1)}+2\mu\frac{\mathbf{x(n)}err(n)}{c+\mathbf{R_{xx}}} \)
where
\( err(n) = d(n)-\mathbf{b(n-1)^Tx(n)} \)
and \( \mathbf{R_{xx}} \) is the correlation matrix
- Parameters
-
Definition at line 290 of file filter.h.
template<class T1, class T2, class T3>
void sp::FIR_filt< T1, T2, T3 >::nlms_adapt |
( |
const T3 |
_err | ) |
|
|
inline |
NLMS Filter update function.
The NLMS filter is updated as
\( \mathbf{b(n)} = \mathbf{b(n-1)}+2\mu\frac{\mathbf{x(n)}err(n)}{c+\mathbf{x(n)^Tx(n)}} \)
where
\( err(n) = d(n)-\mathbf{b(n-1)^Tx(n)} \)
- Parameters
-
Definition at line 229 of file filter.h.
template<class T1, class T2, class T3>
RLS Filter update function.
The RLS filter is updated as
\( \mathbf{b(n)} = \mathbf{b(n-1)}+\mathbf{Kx(n)}\)
where
\( \mathbf{K} =\frac{\mathbf{Px}}{\lambda+\mathbf{x^TPx}} \)
and
\( \mathbf{P^+} =\frac{\mathbf{P^-+xP^-x^T }}{\lambda} \)
- Parameters
-
Definition at line 351 of file filter.h.