Paparazzi UAS
v7.0_unstable
Paparazzi is a free software Unmanned Aircraft System.
|
Matrix decompositions in floating point. More...
#include "math/pprz_simple_matrix.h"
#include "math/pprz_matrix_decomp_float.h"
#include "math/pprz_algebra_float.h"
#include <math.h>
#include <string.h>
Go to the source code of this file.
Macros | |
#define | DEBUG_PRINT(...) |
#define | DEBUG_MAT_PRINT(...) |
Functions | |
void | pprz_cholesky_float (float **out, float **in, int n) |
Cholesky decomposition. More... | |
void | pprz_qr_float (float **Q, float **R, float **in, int m, int n) |
QR decomposition. More... | |
static float | pythag (float a, float b) |
Some SVD decomposition utility macros and functions. More... | |
static int | imin (int num1, int num2) |
int | pprz_svd_float (float **a, float *w, float **v, int m, int n) |
SVD decomposition. More... | |
void | pprz_svd_solve_float (float **x, float **u, float *w, float **v, float **b, int m, int n, int l) |
SVD based linear solver. More... | |
void | fit_linear_model (float *targets, int D, float(*samples)[D], uint16_t count, bool use_bias, float *params, float *fit_error) |
Fit a linear model from samples to target values. More... | |
void | fit_linear_model_prior (float *targets, int D, float(*samples)[D], uint16_t count, bool use_bias, float *priors, float *params, float *fit_error) |
Fit a linear model from samples to target values with a prior. More... | |
Matrix decompositions in floating point.
Definition in file pprz_matrix_decomp_float.c.
#define DEBUG_MAT_PRINT | ( | ... | ) |
Definition at line 39 of file pprz_matrix_decomp_float.c.
#define DEBUG_PRINT | ( | ... | ) |
Definition at line 38 of file pprz_matrix_decomp_float.c.
void fit_linear_model | ( | float * | targets, |
int | D, | ||
float(*) | samples[D], | ||
uint16_t | count, | ||
bool | use_bias, | ||
float * | params, | ||
float * | fit_error | ||
) |
Fit a linear model from samples to target values.
Effectively a wrapper for the pprz_svd_float and pprz_svd_solve_float functions.
[in] | targets | The target values |
[in] | samples | The samples / feature vectors |
[in] | D | The dimensionality of the samples |
[in] | count | The number of samples |
[in] | use_bias | Whether to use the bias. Please note that params should always be of size D+1, but in case of no bias, the bias value is set to 0. |
[out] | parameters* | Parameters of the linear fit |
[out] | fit_error* | Total error of the fit |
Definition at line 529 of file pprz_matrix_decomp_float.c.
References D, MAKE_MATRIX_PTR, MAT_MUL, MAT_SUB, n_samples, parameters, pprz_svd_float(), and pprz_svd_solve_float().
Referenced by RANSAC_linear_model().
void fit_linear_model_prior | ( | float * | targets, |
int | D, | ||
float(*) | samples[D], | ||
uint16_t | count, | ||
bool | use_bias, | ||
float * | priors, | ||
float * | params, | ||
float * | fit_error | ||
) |
Fit a linear model from samples to target values with a prior.
Effectively a wrapper for the pprz_svd_float and pprz_svd_solve_float functions.
[in] | targets | The target values |
[in] | samples | The samples / feature vectors |
[in] | D | The dimensionality of the samples |
[in] | count | The number of samples |
[in] | use_bias | Whether to use the bias. Please note that params should always be of size D+1, but in case of no bias, the bias value is set to 0. |
[in] | priors | Prior per dimension. If use_bias, also for the dimension D+1. |
[out] | parameters* | Parameters of the linear fit |
[out] | fit_error* | Total error of the fit |
Definition at line 617 of file pprz_matrix_decomp_float.c.
References D, DEBUG_MAT_PRINT, DEBUG_PRINT, float_mat_invert(), float_mat_mul(), float_mat_sum(), float_mat_transpose(), MAKE_MATRIX_PTR, MAT_MUL, MAT_SUB, n_samples, and parameters.
|
static |
Definition at line 144 of file pprz_matrix_decomp_float.c.
Referenced by pprz_svd_float().
void pprz_cholesky_float | ( | float ** | out, |
float ** | in, | ||
int | n | ||
) |
Cholesky decomposition.
http://rosettacode.org/wiki/Cholesky_decomposition#C
out | pointer to the output array [n x n] |
in | pointer to the input array [n x n] |
n | dimension of the matrix |
Definition at line 50 of file pprz_matrix_decomp_float.c.
References float_mat_copy(), float_mat_zero(), MAKE_MATRIX_PTR, and s.
void pprz_qr_float | ( | float ** | Q, |
float ** | R, | ||
float ** | in, | ||
int | m, | ||
int | n | ||
) |
QR decomposition.
using Householder method
http://rosettacode.org/wiki/QR_decomposition#C
Q | square orthogonal matrix Q [m x m] |
R | upper triangular matrix R [m x n] |
in | pointer to the input array [m x n] |
m | number of rows of the input matrix |
n | number of column of the input matrix |
Definition at line 89 of file pprz_matrix_decomp_float.c.
References b, float_mat_col(), float_mat_copy(), float_mat_minor(), float_mat_mul(), float_mat_transpose_square(), float_mat_vmul(), float_vect_norm(), float_vect_sdiv(), simple_quad_sim::m, MAKE_MATRIX_PTR, z1, and z2.
int pprz_svd_float | ( | float ** | a, |
float * | w, | ||
float ** | v, | ||
int | m, | ||
int | n | ||
) |
SVD decomposition.
------------------------------------------------------------------— * Reference: "Numerical Recipes By W.H. Press, B. P. Flannery, * S.A. Teukolsky and W.T. Vetterling, Cambridge * University Press, 1986" [BIBLI 08]. * ------------------------------------------------------------------— *
Given a matrix a(m,n), this routine computes its singular value decomposition, A = U · W · Vt. The matrix U replaces a on output. The diagonal matrix of singular values W is output as a vector w(n). The matrix V (not the transpose Vt) is output as v(n,n).
a | input matrix [m x n] and output matrix U [m x n] |
w | output diagonal vector of matrix W [n] |
v | output square matrix V [n x n] |
m | number of rows of input the matrix |
n | number of columns of the input matrix |
Definition at line 167 of file pprz_matrix_decomp_float.c.
References H, imin(), simple_quad_sim::m, pythag(), mesonh.mesonh_atmosphere::X, mesonh.mesonh_atmosphere::Y, and mesonh.mesonh_atmosphere::Z.
Referenced by fit_linear_flow_field(), fit_linear_model(), fit_linear_model_OF(), and pprz_polyfit_float().
void pprz_svd_solve_float | ( | float ** | x, |
float ** | u, | ||
float * | w, | ||
float ** | v, | ||
float ** | b, | ||
int | m, | ||
int | n, | ||
int | l | ||
) |
SVD based linear solver.
Solves A · X = B for a vector X, where A is specified by the arrays u, w, v as returned by pprz_svd_float. m and n are the dimensions of a. b(m) is the input right-hand side. x(n) is the output solution vector. No input quantities are destroyed, so the routine may be called sequentially with different b's.
x | solution of the system ([n x l] matrix) |
u | U matrix from SVD decomposition |
w | diagonal of the W matrix from the SVD decomposition |
v | V matrrix from SVD decomposition |
b | right-hand side input matrix from system to solve (column vector [m x l]) |
m | number of rows of the matrix A |
n | number of columns of the matrix A |
l | number of columns of the matrix B |
Definition at line 494 of file pprz_matrix_decomp_float.c.
References b, simple_quad_sim::m, and s.
Referenced by fit_linear_flow_field(), fit_linear_model(), fit_linear_model_OF(), and pprz_polyfit_float().
|
inlinestatic |
Some SVD decomposition utility macros and functions.
Definition at line 128 of file pprz_matrix_decomp_float.c.
References b.
Referenced by pprz_svd_float().