Lapack functions
[Linear Algebra]

Collaboration diagram for Lapack functions:

Wrappers or bindings for Lapack functions are provided. More...

Classes

struct  mimas::gglse_t< T >
 Function-object for gglse. More...
struct  mimas::inv_t< T >
 Function-object for inv. More...
struct  mimas::gels_t< T >
 Function-object for gels. More...
struct  mimas::gesvd_t< T >
 Function-object for gesvd. More...
struct  mimas::gesdd_t< T >
 Function-wrapper for gesdd. More...
struct  mimas::syev_t< T >
 Function-wrapper for syev. More...
struct  mimas::pptrf_t< T >
 Function-object for pptrf. More...

Functions

template<typename T>
boost::numeric::ublas::vector<
T > 
mimas::gglse (const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &A, const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &B, const boost::numeric::ublas::vector< T > &c, const boost::numeric::ublas::vector< T > &d)
 Solve the linear equality-constrained least squares (LSE) problem.
template<typename T>
boost::numeric::ublas::matrix<
T, boost::numeric::ublas::column_major > 
mimas::inv (const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &A)
 Compute matrix inverse by using LU factorization.
template<typename T>
boost::numeric::ublas::vector<
T > 
mimas::gels (const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &A, const boost::numeric::ublas::vector< T > &b)
 solve overdetermined linear systems using QR factorization.
template<typename T>
boost::numeric::ublas::banded_matrix<
T > 
mimas::gesvd (const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &A, boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *U, boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *Vt)
 Compute the singular value decomposition (SVD) of a matrix A.
template<typename T>
boost::numeric::ublas::banded_matrix<
T > 
mimas::gesdd (const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &A, boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *U, boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *Vt)
 Compute the singular value decomposition (SVD) of a matrix A (fast).
template<typename T>
boost::numeric::ublas::diagonal_matrix<
T > 
mimas::syev (const boost::numeric::ublas::symmetric_matrix< T, boost::numeric::ublas::upper > &A, boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *E)
 Compute all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.
template<typename T>
boost::numeric::ublas::triangular_matrix<
T, boost::numeric::ublas::upper,
boost::numeric::ublas::column_major > 
mimas::pptrf (const boost::numeric::ublas::symmetric_matrix< T, boost::numeric::ublas::upper, boost::numeric::ublas::column_major > &A)
 Compute the Cholesky decomposition of a symmetric matrix A.

Detailed Description

Wrappers or bindings for Lapack functions are provided.

The matrix- and vector-classes of boost are used as datatypes.

Note, that lapack uses matrices in column-major orientation. You are forced to convert to

    boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major >
The default orientation of boost matrices is row-major.

Author:
Bala Amavasai (bala@amavasai.org)

Jan Wedekind (jan at wedesoft.de)

Todo:
Extend the solver-classes to complex values.
Date:
Fri May 31 12:07:52 UTC 2002

Function Documentation

template<typename T>
boost::numeric::ublas::vector< T > mimas::gels ( const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &  A,
const boost::numeric::ublas::vector< T > &  b 
)

solve overdetermined linear systems using QR factorization.

Solve the minimum least square problem. $A$ and $\vec{b}$ are given and $\vec{x}$ is the desired solution for

\[\vec{x}:=\displaystyle\mathop{argmin}_{\vec{x}}|\vec{b}-A\,\vec{x}|\]

This is much faster than computing $\vec{x}:=(A^\top\,A)^{-1}\,A^\top\,\vec{b}$ with the matrix-classes.

See manpage dgels(1) for more documentation. The parametrisation for underdetermined systems wasn't working as expected and was not included in the interface therefore.

See also:
gels_t

Definition at line 218 of file linalg.h.

template<typename T>
boost::numeric::ublas::banded_matrix< T > mimas::gesdd ( const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &  A,
boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *  U,
boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *  Vt 
)

Compute the singular value decomposition (SVD) of a matrix A (fast).

Same as mimas::gesvd but faster algorithm (divide-and-conquer approach).

Note that the singular value decomposition doesn't necessarily converge. In the later case an exception is thrown.

See manpage dgesdd(1) for more documentation.

Note, that, in contrast to the interface of mimas::gesvd, U and Vt must both be defined or both be NULL!

Parameters:
A Matrix to be decomposed.
U Returns matrix with left hand singular vectors $U$ (may be NULL, if it is of no interest).
Vt Returns transposed matrix with right hand singular vectors `$V^\top$ (may be NULL, if it is of no interest).
Returns:
Diagonal matrix $\Sigma$.
See also:
gesvd

gesdd_t

Definition at line 308 of file linalg.h.

template<typename T>
boost::numeric::ublas::banded_matrix< T > mimas::gesvd ( const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &  A,
boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *  U,
boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *  Vt 
)

Compute the singular value decomposition (SVD) of a matrix A.

$A$ is given. $U$, $\Sigma$ and $V^\top$ have to be determined, such that: $A=U\,\Sigma\,V^\top$ and $\Sigma=\mathop{diag}(\sigma_1,\sigma_2,\ldots,\sigma_n)$ diagonal matrix with the singular values $\sigma_1,\sigma_2,\ldots,\sigma_n$ ($\sigma_i\ge 0$) in descending order. The singular values are sorted in descending order.

Note that the singular value decomposition doesn't necessarily converge. In the later case an exception is thrown.

See mimas::gesdd, which is a faster algorithm.

See manpage dgesvd(1) for more documentation.

Parameters:
A Matrix to be decomposed.
U Returns matrix with left hand singular vectors $U$ (may be NULL, if it is of no interest).
Vt Returns transposed matrix with right hand singular vectors $V^\top$ (may be NULL, if it is of no interest).
Returns:
Diagonal matrix $\Sigma$.
Todo:
Does not converge on AMD64-processor.
See also:
gesdd

gesvd_t

Definition at line 265 of file linalg.h.

template<typename T>
boost::numeric::ublas::vector< T > mimas::gglse ( const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &  A,
const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &  B,
const boost::numeric::ublas::vector< T > &  c,
const boost::numeric::ublas::vector< T > &  d 
)

Solve the linear equality-constrained least squares (LSE) problem.

$A$, $\vec{c}$, $B$ and $\vec{d}$ are given and $\vec{x}$ is sought:

minimize $\vec{x}:=\displaystyle\mathop{argmin}_{\vec{x}}|\vec{c}-A\,\vec{x}|$ subject to $B\,\vec{x}=\vec{d}$.

See manpage dgglse(1) for more documentation.

See also:
gglse_t

Definition at line 161 of file linalg.h.

template<typename T>
boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > mimas::inv ( const boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > &  A  ) 

Compute matrix inverse by using LU factorization.

The matrix $A\in\mathbf{K}^{n\times n}$ is given and the inverse $A^{-1}$ has to be computed.

See manpages dgetrf(1) and dgetri(1) for more documentation.

See also:
inv_t

Definition at line 188 of file linalg.h.

template<typename T>
boost::numeric::ublas::triangular_matrix< T, boost::numeric::ublas::upper, boost::numeric::ublas::column_major > mimas::pptrf ( const boost::numeric::ublas::symmetric_matrix< T, boost::numeric::ublas::upper, boost::numeric::ublas::column_major > &  A  ) 

Compute the Cholesky decomposition of a symmetric matrix A.

Given a symmetric positive-definite matrix$A\in\mathbf{K}^{n\times n}$, the Cholesky decomposition returns the upper triangular matrix $U\in\mathbf{K}^{n\times n}$ such that:

\[A=U^\top\,U\]

Parameters:
A Symmetric matrix (with upper and column-based storage) to compute Cholesky decomposition for.
Returns:
Upper triangular matrix (with column-based storage).
Author:
Julien Faucher (faucherj@gmail.com)
See also:
pptrf_t

Definition at line 396 of file linalg.h.

template<typename T>
boost::numeric::ublas::diagonal_matrix< T > mimas::syev ( const boost::numeric::ublas::symmetric_matrix< T, boost::numeric::ublas::upper > &  A,
boost::numeric::ublas::matrix< T, boost::numeric::ublas::column_major > *  E 
)

Compute all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

The eigenvalues and eigenvectors of a symmetric matrix $A\in\mathbf{K}^{n\times n}$ can be computed.

The eigenvectors are orthogonal by nature and normalised by definition. If the matrix $E$ contains the eigenvectors of $A$, then the eigentransform can be written as ($E$ is orthonormal $\Rightarrow$ $E^{-1}=E^\top$):

\[A=E\,\Lambda\,E^\top\]

Where $\Lambda=\mathop{diag}(\lambda_1,\lambda_2,\ldots,\lambda_n)$ is a diagonal matrix with the eigenvalues $\lambda_1,\lambda_2,\ldots,\lambda_n$ in ascending order.

See manpage dsyev(1) for more documentation.

Compute eigenvalues and, optionally, the eigenvectors of A.

Parameters:
A Symmetric matrix (with upper storage) to compute eigentransform for.
E Pointer to matrix, where the eigenvectors have to be written to (may be NULL, if it is of no interest).
Returns:
Diagonal matrix with eigenvalues in ascending order.
See also:
syev_t

Definition at line 357 of file linalg.h.


[GNU/Linux] [Qt] [Mesa] [STL] [Lapack] [Boost] [Magick++] [Xalan-C and Xerces-C] [doxygen] [graphviz] [FFTW] [popt] [xine] [Gnuplot] [gnu-arch] [gcc] [gstreamer] [autoconf/automake/make] [freshmeat.net] [opensource.org] [sourceforge.net] [MMVL]
mimas 2.1 - Copyright Mon Oct 30 11:31:22 2006, Bala Amavasai, Stuart Meikle, Arul Selvan, Fabio Caparrelli, Jan Wedekind, Manuel Boissenin, ...