35 #define DEBUG_PRINT printf
36 #define DEBUG_MAT_PRINT MAT_PRINT
38 #define DEBUG_PRINT(...)
39 #define DEBUG_MAT_PRINT(...)
57 for (i = 0; i < n; i++) {
58 for (j = 0; j < (i + 1); j++) {
60 for (k = 0; k < j; k++) {
61 s += o[i][k] * o[j][k];
64 o[i][j] = sqrtf(in[i][i] - s);
67 o[i][j] = 1.0 / o[j][j] * (in[i][j] -
s);
94 float _z[
m][n], _z1[
m][n], _z2[
m][
m];
98 for (i = 0; i <
m; i++)
for (k = 0; k <
m; k++) { q[i][k] = &_q[i][k][0]; }
100 for (k = 0; k < n && k < m - 1; k++) {
101 float e[
m], x[
m], a, b;
105 if (in[k][k] > 0) { a = -a; }
106 for (i = 0; i <
m; i++) {
107 e[i] = (i == k) ? 1 : 0;
108 e[i] = x[i] + a * e[i];
116 for (i = 1; i < n && i < m - 1; i++) {
128 static inline float pythag(
float a,
float b)
134 if (absa == 0) {
return 0.0; }
135 return (absa * sqrtf(1.0 + (absb / absa) * (absb / absa)));
136 }
else if (absb == 0.0) {
139 return (absb * sqrtf(1.0 + (absa / absb) * (absa / absb)));
167 int flag, i, its, j, jj, k, l, NM;
168 float C, F,
H, S,
X,
Y,
Z, tmp;
174 for (i = 0; i < n; ++i) {
181 for (k = i; k <
m; ++k) {
182 Scale = Scale + fabsf(a[k][i]);
185 for (k = i; k <
m; ++k) {
186 a[k][i] = a[k][i] / Scale;
187 S = S + a[k][i] * a[k][i];
197 for (j = l; j < n; ++j) {
199 for (k = i; k <
m; ++k) {
200 S = S + a[k][i] * a[k][j];
207 for (k = i; k <
m; ++k) {
208 a[k][j] = a[k][j] + F * a[k][i];
212 for (k = i; k <
m; ++k) {
213 a[k][i] = Scale * a[k][i];
222 if ((i < m) && (i != (n - 1))) {
223 for (k = l; k < n; ++k) {
224 Scale = Scale + fabsf(a[i][k]);
227 for (k = l; k < n; ++k) {
228 a[i][k] = a[i][k] / Scale;
229 S = S + a[i][k] * a[i][k];
241 for (k = l; k < n; ++k) {
242 rv1[k] = a[i][k] /
H;
245 for (j = l; j <
m; ++j) {
247 for (k = l; k < n; ++k) {
248 S = S + a[j][k] * a[i][k];
250 for (k = l; k < n; ++k) {
251 a[j][k] = a[j][k] + S * rv1[k];
255 for (k = l; k < n; ++k) {
256 a[i][k] = Scale * a[i][k];
260 tmp = fabsf(w[i]) + fabsf(rv1[i]);
267 for (i = n - 1; i >= 0; --i) {
270 for (j = l; j < n; ++j) {
272 v[j][i] = (a[i][j] / a[i][l]) / G;
277 for (j = l; j < n; ++j) {
279 for (k = l; k < n; ++k) {
280 S = S + a[i][k] * v[k][j];
282 for (k = l; k < n; ++k) {
283 v[k][j] = v[k][j] + S * v[k][i];
287 for (j = l; j < n; ++j) {
298 for (i = n - 1; i >= 0; --i) {
302 for (j = l; j < n; ++j) {
309 for (j = l; j < n; ++j) {
311 for (k = l; k <
m; ++k) {
312 S = S + a[k][i] * a[k][j];
314 F = (S / a[i][i]) * G;
315 for (k = i; k <
m; ++k) {
316 a[k][j] = a[k][j] + F * a[k][i];
320 for (j = i; j <
m; ++j) {
321 a[j][i] = a[j][i] * G;
324 for (j = i; j <
m; ++j) {
328 a[i][i] = a[i][i] + 1.0;
333 for (k = (n - 1); k >= 0; --k) {
335 for (its = 1; its <= 30; ++its) {
339 for (l = k; l >= 0; --l) {
341 if ((fabsf(rv1[l]) + ANorm) == ANorm) {
344 }
else if ((fabsf(w[NM]) + ANorm) == ANorm) {
353 for (i = l; i <= k; ++i) {
355 if ((fabsf(F) + ANorm) != ANorm) {
367 for (j = 0; j <
m; ++j) {
370 a[j][NM] = (Y * C) + (Z * S);
371 a[j][i] = -(Y * S) + (Z * C);
382 for (j = 0; j < n; ++j) {
400 F = ((Y -
Z) * (Y + Z) + (G -
H) * (G + H)) / (2.0 * H * Y);
410 if ((F + tmp) != 0) {
411 F = ((X -
Z) * (X + Z) + H * ((Y / (F + tmp)) -
H)) /
X;
419 for (j = l; j <= NM; ++j) {
435 F = (X * C) + (G * S);
436 G = -(X * S) + (G * C);
439 for (jj = 0; jj < n; ++jj) {
442 v[jj][j] = (X * C) + (Z * S);
443 v[jj][i] = -(X * S) + (Z * C);
455 F = (C * G) + (S * Y);
456 X = -(S * G) + (C * Y);
457 for (jj = 0; jj <
m; ++jj) {
460 a[jj][j] = (Y * C) + (Z * S);
461 a[jj][i] = -(Y * S) + (Z * C);
496 for (k = 0; k < l; k++) {
497 for (j = 0; j < n; j++) {
500 for (i = 0; i <
m; i++) { s += u[i][j] * b[i][k]; }
505 for (j = 0; j < n; j++) {
507 for (jj = 0; jj < n; jj++) { s += v[j][jj] * tmp[jj]; }
540 n_samples = (n_samples < D_1) ? D_1 : n_samples;
542 n_samples = (n_samples < count) ? n_samples : count;
546 float _AA[count][D_1];
548 float _targets_all[count][1];
551 for (sam = 0; sam < count; sam++) {
552 for (d = 0; d <
D; d++) {
553 AA[sam][d] = samples[sam][d];
560 targets_all[sam][0] = targets[sam];
565 float _parameters[D_1][1];
584 MAT_MUL(count, D_1, 1, bb,
AA, parameters);
586 MAT_SUB(count, 1, C, bb, targets_all);
588 for (sam = 0; sam < count; sam++) {
589 *fit_error += fabsf(C[sam][0]);
595 for (d = 0; d < D_1; d++) {
596 params[d] = parameters[d][0];
615 float *params,
float *fit_error)
636 n_samples = (n_samples < D_1) ? D_1 : n_samples;
638 n_samples = (n_samples < count) ? n_samples : count;
645 float _AA[count][D_1];
647 float _AAT[D_1][count];
649 float _AATAA[D_1][D_1];
651 float _PRIOR[D_1][D_1];
653 for(d = 0; d <
D; d++) {
654 PRIOR[d][d] = priors[d];
657 PRIOR[
D][
D] = priors[
D];
663 float _targets_all[count][1];
666 for (sam = 0; sam < count; sam++) {
667 for (d = 0; d <
D; d++) {
668 AA[sam][d] = samples[sam][d];
675 targets_all[sam][0] = targets[sam];
699 float _INV_AATAA[D_1][D_1];
727 float _PINV[D_1][count];
730 MAT_MUL(D_1, D_1, count, PINV, INV_AATAA, AAT);
735 float _parameters[D_1][1];
737 MAT_MUL(D_1, count, 1, parameters, PINV, targets_all);
754 MAT_MUL(count, D_1, 1, bb,
AA, parameters);
757 MAT_SUB(count, 1, C, bb, targets_all);
760 for (sam = 0; sam < count; sam++) {
761 *fit_error += fabsf(C[sam][0]);
767 for (d = 0; d < D_1; d++) {
768 params[d] = parameters[d][0];
void pprz_cholesky_float(float **out, float **in, int n)
Cholesky decomposition.
static float float_vect_norm(const float *a, const int n)
||a||
Simple matrix helper macros.
#define MAKE_MATRIX_PTR(_ptr, _mat, _rows)
Make a pointer to a matrix of _rows lines.
void float_mat_invert(float **o, float **mat, int n)
Calculate inverse of any n x n matrix (passed as C array) o = mat^-1 Algorithm verified with Matlab...
static void float_mat_transpose_square(float **a, int n)
transpose square matrix
static float pythag(float a, float b)
Some SVD decomposition utility macros and functions.
static void float_mat_copy(float **a, float **b, int m, int n)
a = b
static void float_mat_col(float *o, float **a, int m, int c)
o = c-th column of matrix a[m x n]
static void float_mat_vmul(float **o, float *v, int n)
o = I - v v^T
void pprz_qr_float(float **Q, float **R, float **in, int m, int n)
QR decomposition.
static void float_vect_sdiv(float *o, const float *a, const float s, const int n)
o = a / s
static void float_mat_minor(float **o, float **a, int m, int n, int d)
matrix minor
Paparazzi floating point algebra.
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.
static void float_mat_zero(float **a, int m, int n)
a = 0
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.
Matrix decompositions in floating point.
#define DEBUG_MAT_PRINT(...)
int pprz_svd_float(float **a, float *w, float **v, int m, int n)
SVD decomposition.
static struct FloatVect3 H
#define MAT_SUB(_i, _j, C, A, B)
static void float_mat_transpose(float **o, float **a, int n, int m)
transpose non-square matrix
static void float_mat_sum(float **o, float **a, float **b, int m, int n)
o = a + b
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.
static void float_mat_mul(float **o, float **a, float **b, int m, int n, int l)
o = a * b
#define MAT_MUL(_i, _k, _j, C, A, B)