Robotran C Documentation
Functions
mbs_linalg.h File Reference

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...
 

Detailed Description

This c header file declares various functions about linear algebra.

Creation date: 2006

Author
Robotran team

(c) Universite catholique de Louvain

Function Documentation

◆ mbs_choldc_0()

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 !

Parameters
[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

◆ mbs_cholsl_0()

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 /!\/!\/!\

Parameters
[in,out]ACOLUMN major the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A
[in]nsize of the A matrix
[in,out]BOn entry - In COLUMN major - the right hand side matrix B. On exit - In COLUMN major - the solution matrix X.
[in]nbnumber of vectors for the RHS matrix B.

◆ mbs_eig_0()

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.

Parameters
[in]Amatrix for the Eigenvalue problem.
[in]nsize of the matrix
[out]eval_aResult of the Eigenvalue problem: real part of the 2*nx eigen values [n X 1]
[out]eval_bResult of the Eigenvalue problem: imag part of the 2*nx eigen values [n X 1]
[out]evec_aResult of the Eigenvalue problem: real part of the eigen vectors (each column is a eigenvector) [n X n]
[out]evec_bResult of the Eigenvalue problem: imag part of the eigen vectors (each column is a eigenvector) [n X n]
Returns
Error status, <0 in case of failure.

◆ mbs_invLU_0()

int mbs_invLU_0 ( double **  A,
double **  Am1,
int  n 
)

Inverse the square matrix A using a LU decomposition.

Parameters
[in]A: the matrix to inverse
[out]Am1: the inverse of the matrix A
[in]n: size of the matrix A
Returns
0 for successful exit.

◆ mbs_over_under_determined()

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 |).

Parameters
[in,out]Amain matrix of size [x,y] with first line and row unused using row major memory.
[in,out]bAs 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].

Parameters
[in]xthe number of row of A, which is the number of equations.
[in]ythe number of col. of A, which is the number of unknowns.
Returns
error status, < 0 in case of failure.

◆ mbs_rank_0()

int mbs_rank_0 ( double **  A,
int  m,
int  n 
)

Compute the rank of matrix_in.

Parameters
[in]Amatrix for the rank computation
[in]mmatrix first dimension
[in]nmatrix second dimension

using Lapacke libraries ...

Threshold for considering a singular value to be 0

◆ mbs_schur_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)

Returns
0 for successful exit.

◆ mbs_svdDcmp_0()

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)

Parameters
[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.