Go to the source code of this file.
Functions | |
int | mbs_eig_0 (double **A, int n, double *eval_a, double *eval_b, double **evec_a, double **evec_b) |
Compute the eigen values and vector of a matrix with index starting a 0. More... | |
int | mbs_schur_0 (double **A, int n, double **T, double **Z, double *WR, double *WI, int reorder) |
A : matrix for the Eigenvalue problem in row major n : size of the matrix WR : Result of the Eigenvalue problem: real part of the n eigen values (ordered as they appear on the diagonal of the output Schur form T [n x 1] WI : Result of the Eigenvalue problem: real part of the n eigen values (ordered as they appear on the diagonal of the output Schur form T [n x 1] T : (Re-ordered) schur decomposition matrix Z : (Re-ordered) orthogonal matrix Z of Schur vectors [n X n] reorder : flag for re-ordering 1=re-ordered 0=non-ordered (default 1) More... | |
int | mbs_rank_0 (double **A, int m, int n) |
Compute the rank of matrix_in. More... | |
int | mbs_choldc_0 (double **A, int n) |
Compute the Cholesky decomposition of A, i.e. More... | |
int | mbs_cholsl_0 (double **A, int n, double **B, int nb) |
Solves a symmetric positive definite system of linear equations AX=B (B being the right hand side with possibly more than 1 column) using the Cholesky factorization computed by mbs_choldc_0 in COL MAJOR ! More... | |
int | mbs_invLU_0 (double **A, double **Am1, int n) |
Inverse the square matrix A using a LU decomposition. More... | |
int | mbs_svdDcmp_0 (double **A, int x, int y, double **U, double **S, double **VT) |
Realize the singular value decomposition. More... | |
int | mbs_over_under_determined (double **A, double *b, int x, int y) |
Solve a linear system of equations, A*x=b, using dgelss lapack function. More... | |
This c header file declares various functions about linear algebra.
Creation date: 2006
(c) Universite catholique de Louvain
int mbs_choldc_0 | ( | double ** | A, |
int | n | ||
) |
Compute the Cholesky decomposition of A, i.e.
the upper triangular matrix L
/!\/!\/!\ the upper triangular matrix is expressed in COL major /!\/!\/!\
WARNING: A is modified !
[in,out] | A | : IN : real symmetric positive definite matrix. OUT: in COLUMN major the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A |
int mbs_cholsl_0 | ( | double ** | A, |
int | n, | ||
double ** | B, | ||
int | nb | ||
) |
Solves a symmetric positive definite system of linear equations AX=B (B being the right hand side with possibly more than 1 column) using the Cholesky factorization computed by mbs_choldc_0 in COL MAJOR !
WARNING: B is modified and the solution is stored in the matrix !
/!\/!\/!\ the B and the X(solution) matrices are expressed in COL major /!\/!\/!\
[in,out] | A | COLUMN major the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A |
[in] | n | size of the A matrix |
[in,out] | B | On entry - In COLUMN major - the right hand side matrix B. On exit - In COLUMN major - the solution matrix X. |
[in] | nb | number of vectors for the RHS matrix B. |
int mbs_eig_0 | ( | double ** | A, |
int | n, | ||
double * | eval_a, | ||
double * | eval_b, | ||
double ** | evec_a, | ||
double ** | evec_b | ||
) |
Compute the eigen values and vector of a matrix with index starting a 0.
[in] | A | matrix for the Eigenvalue problem. |
[in] | n | size of the matrix |
[out] | eval_a | Result of the Eigenvalue problem: real part of the 2*nx eigen values [n X 1] |
[out] | eval_b | Result of the Eigenvalue problem: imag part of the 2*nx eigen values [n X 1] |
[out] | evec_a | Result of the Eigenvalue problem: real part of the eigen vectors (each column is a eigenvector) [n X n] |
[out] | evec_b | Result of the Eigenvalue problem: imag part of the eigen vectors (each column is a eigenvector) [n X n] |
int mbs_invLU_0 | ( | double ** | A, |
double ** | Am1, | ||
int | n | ||
) |
Inverse the square matrix A using a LU decomposition.
[in] | A | : the matrix to inverse |
[out] | Am1 | : the inverse of the matrix A |
[in] | n | : size of the matrix A |
int mbs_over_under_determined | ( | double ** | A, |
double * | b, | ||
int | x, | ||
int | y | ||
) |
Solve a linear system of equations, A*x=b, using dgelss lapack function.
This function mainly transfers the MBsysC arrays with index starting at 1 (see get_dmat_1() and get_dvec_1()) to Lapack-compatible arrays (index starting at 0, col. major memory).
The following explanations come from the Lapack documentation. This function solves overdetermined or underdetermined systems.
For overdetermined system it minimizes 2-norm(| b - A*x |).
For underdetermined system it minimizes 2-norm(| x |).
[in,out] | A | main matrix of size [x,y] with first line and row unused using row major memory. |
[in,out] | b | As input it contains the vector b of size [y] with index starting at 1. As output it contains the solution vector x. |
For overdetermined system, the residual sum-of-squares for the solution is the sum of the elements from b[x+1] to b[y+1].
[in] | x | the number of row of A, which is the number of equations. |
[in] | y | the number of col. of A, which is the number of unknowns. |
int mbs_rank_0 | ( | double ** | A, |
int | m, | ||
int | n | ||
) |
Compute the rank of matrix_in.
[in] | A | matrix for the rank computation |
[in] | m | matrix first dimension |
[in] | n | matrix second dimension |
using Lapacke libraries ...
Threshold for considering a singular value to be 0
int mbs_schur_0 | ( | double ** | A, |
int | n, | ||
double ** | T, | ||
double ** | Z, | ||
double * | WR, | ||
double * | WI, | ||
int | reorder | ||
) |
A : matrix for the Eigenvalue problem in row major n : size of the matrix WR : Result of the Eigenvalue problem: real part of the n eigen values (ordered as they appear on the diagonal of the output Schur form T [n x 1] WI : Result of the Eigenvalue problem: real part of the n eigen values (ordered as they appear on the diagonal of the output Schur form T [n x 1] T : (Re-ordered) schur decomposition matrix Z : (Re-ordered) orthogonal matrix Z of Schur vectors [n X n] reorder : flag for re-ordering 1=re-ordered 0=non-ordered (default 1)
int mbs_svdDcmp_0 | ( | double ** | A, |
int | x, | ||
int | y, | ||
double ** | U, | ||
double ** | S, | ||
double ** | VT | ||
) |
Realize the singular value decomposition.
A = U * S * transpose(V)
[in] | A | : the rectangular matrix |
[in] | x | : the number of row of A |
[in] | y | : the number of row of A |
[in] | U | : the x-by-x orthogonal matrix U; |
[in] | S | : the x-by-y matrix which is zero except for its min(x,y) diagonal elements |
[in] | VT | : the y-by-y orthogonal matrix V**T. |