Functions | |
bool | itpp::chol (const mat &X, mat &F) |
Cholesky factorisation of real symmetric and positive definite matrix. | |
mat | itpp::chol (const mat &X) |
Cholesky factorisation of real symmetric and positive definite matrix. | |
bool | itpp::chol (const cmat &X, cmat &F) |
Cholesky factorisation of complex hermitian and positive-definite matrix. | |
cmat | itpp::chol (const cmat &X) |
Cholesky factorisation of complex hermitian and positive-definite matrix. | |
bool | itpp::eig_sym (const mat &A, vec &d, mat &V) |
Calculates the eigenvalues and eigenvectors of a symmetric real matrix. | |
bool | itpp::eig_sym (const mat &A, vec &d) |
Calculates the eigenvalues of a symmetric real matrix. | |
vec | itpp::eig_sym (const mat &A) |
Calculates the eigenvalues of a symmetric real matrix. | |
bool | itpp::eig_sym (const cmat &A, vec &d, cmat &V) |
Calculates the eigenvalues and eigenvectors of a hermitian complex matrix. | |
bool | itpp::eig_sym (const cmat &A, vec &d) |
Calculates the eigenvalues of a hermitian complex matrix. | |
vec | itpp::eig_sym (const cmat &A) |
Calculates the eigenvalues of a hermitian complex matrix. | |
bool | itpp::eig (const mat &A, cvec &d, cmat &V) |
Calculates the eigenvalues and eigenvectors of a real non-symmetric matrix. | |
bool | itpp::eig (const mat &A, cvec &d) |
Calculates the eigenvalues of a real non-symmetric matrix. | |
cvec | itpp::eig (const mat &A) |
Calculates the eigenvalues of a real non-symmetric matrix. | |
bool | itpp::eig (const cmat &A, cvec &d, cmat &V) |
Calculates the eigenvalues and eigenvectors of a complex non-hermitian matrix. | |
bool | itpp::eig (const cmat &A, cvec &d) |
Calculates the eigenvalues of a complex non-hermitian matrix. | |
cvec | itpp::eig (const cmat &A) |
Calculates the eigenvalues of a complex non-hermitian matrix. | |
bool | itpp::lu (const mat &X, mat &L, mat &U, ivec &p) |
LU factorisation of real matrix. | |
bool | itpp::lu (const cmat &X, cmat &L, cmat &U, ivec &p) |
LU factorisation of real matrix. | |
void | itpp::interchange_permutations (vec &b, const ivec &p) |
Makes swapping of vector b according to the interchange permutation vector p. | |
bmat | itpp::permutation_matrix (const ivec &p) |
Make permutation matrix P from the interchange permutation vector p. | |
bool | itpp::qr (const mat &A, mat &Q, mat &R) |
QR factorisation of real matrix. | |
bool | itpp::qr (const mat &A, mat &Q, mat &R, bmat &P) |
QR factorisation of real matrix with pivoting. | |
bool | itpp::qr (const cmat &A, cmat &Q, cmat &R) |
QR factorisation of a complex matrix. | |
bool | itpp::qr (const cmat &A, cmat &Q, cmat &R, bmat &P) |
QR factorisation of a complex matrix with pivoting. | |
bool | itpp::schur (const mat &A, mat &U, mat &T) |
Schur decomposition of a real matrix. | |
mat | itpp::schur (const mat &A) |
Schur decomposition of a real matrix. | |
bool | itpp::schur (const cmat &A, cmat &U, cmat &T) |
Schur decomposition of a complex matrix. | |
cmat | itpp::schur (const cmat &A) |
Schur decomposition of a complex matrix. | |
bool | itpp::svd (const mat &A, vec &s) |
Get singular values s of a real matrix A using SVD. | |
bool | itpp::svd (const cmat &A, vec &s) |
Get singular values s of a complex matrix A using SVD. | |
vec | itpp::svd (const mat &A) |
Return singular values of a real matrix A using SVD. | |
vec | itpp::svd (const cmat &A) |
Return singular values of a complex matrix A using SVD. | |
bool | itpp::svd (const mat &A, mat &U, vec &s, mat &V) |
Perform Singular Value Decomposition (SVD) of a real matrix A . | |
bool | itpp::svd (const cmat &A, cmat &U, vec &s, cmat &V) |
Perform Singular Value Decomposition (SVD) of a complex matrix A . |
Cholesky factorisation of real symmetric and positive definite matrix.
The Cholesky factorisation of a real symmetric positive-definite matrix of size
is given by
where is an upper triangular
matrix.
Returns true if calculation succeeded. False otherwise.
Cholesky factorisation of real symmetric and positive definite matrix.
The Cholesky factorisation of a real symmetric positive-definite matrix of size
is given by
where is an upper triangular
matrix.
Cholesky factorisation of complex hermitian and positive-definite matrix.
The Cholesky factorisation of a hermitian positive-definite matrix of size
is given by
where is an upper triangular
matrix.
Returns true if calculation succeeded. False otherwise.
If X
is positive definite, true is returned and F=chol
(X) produces an upper triangular F
. If also X
is symmetric then F'*F
= X. If X
is not positive definite, false is returned.
Cholesky factorisation of complex hermitian and positive-definite matrix.
The Cholesky factorisation of a hermitian positive-definite matrix of size
is given by
where is an upper triangular
matrix.
Calculates the eigenvalues and eigenvectors of a symmetric real matrix.
The Eigenvalues and the eigenvectors
of the real and symmetric
matrix
satisfies
The eigenvectors are the columns of the matrix V. True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine DSYEV.
Calculates the eigenvalues of a symmetric real matrix.
The Eigenvalues and the eigenvectors
of the real and symmetric
matrix
satisfies
True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine DSYEV.
Calculates the eigenvalues of a symmetric real matrix.
The Eigenvalues and the eigenvectors
of the real and symmetric
matrix
satisfies
Uses the LAPACK routine DSYEV.
Calculates the eigenvalues and eigenvectors of a hermitian complex matrix.
The Eigenvalues and the eigenvectors
of the complex and hermitian
matrix
satisfies
The eigenvectors are the columns of the matrix V. True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine ZHEEV.
Calculates the eigenvalues of a hermitian complex matrix.
The Eigenvalues and the eigenvectors
of the complex and hermitian
matrix
satisfies
True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine ZHEEV.
Calculates the eigenvalues of a hermitian complex matrix.
The Eigenvalues and the eigenvectors
of the complex and hermitian
matrix
satisfies
Uses the LAPACK routine ZHEEV.
Calculates the eigenvalues and eigenvectors of a real non-symmetric matrix.
The Eigenvalues and the eigenvectors
of the real
matrix
satisfies
The eigenvectors are the columns of the matrix V. True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine DGEEV.
Calculates the eigenvalues of a real non-symmetric matrix.
The Eigenvalues and the eigenvectors
of the real
matrix
satisfies
True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine DGEEV.
Calculates the eigenvalues of a real non-symmetric matrix.
The Eigenvalues and the eigenvectors
of the real
matrix
satisfies
Uses the LAPACK routine DGEEV.
Calculates the eigenvalues and eigenvectors of a complex non-hermitian matrix.
The Eigenvalues and the eigenvectors
of the complex
matrix
satisfies
The eigenvectors are the columns of the matrix V. True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine ZGEEV.
Calculates the eigenvalues of a complex non-hermitian matrix.
The Eigenvalues and the eigenvectors
of the complex
matrix
satisfies
True is returned if the calculation was successful. Otherwise false.
Uses the LAPACK routine ZGEEV.
Calculates the eigenvalues of a complex non-hermitian matrix.
The Eigenvalues and the eigenvectors
of the complex
matrix
satisfies
Uses the LAPACK routine ZGEEV.
LU factorisation of real matrix.
The LU factorization of the real matrix of size
is given by
where and
are lower and upper triangular matrices and
is a permutation matrix.
The interchange permutation vector p is such that k and p(k) should be changed for all k. Given this vector a permutation matrix can be constructed using the function
bmat permutation_matrix(const ivec &p)
If X is an n by n matrix lu(X,L,U,p) computes the LU decomposition. L is a lower triangular, U an upper triangular matrix. p is the interchange permutation vector such that k and p(k) should be changed for all k.
Returns true is calculation succeeds. False otherwise.
LU factorisation of real matrix.
The LU factorization of the complex matrix of size
is given by
where and
are lower and upper triangular matrices and
is a permutation matrix.
The interchange permutation vector p is such that k and p(k) should be changed for all k. Given this vector a permutation matrix can be constructed using the function
bmat permutation_matrix(const ivec &p)
If X is an n by n matrix lu(X,L,U,p) computes the LU decomposition. L is a lower triangular, U an upper triangular matrix. p is the interchange permutation vector such that elements k and row p(k) should be interchanged.
Returns true is calculation succeeds. False otherwise.
QR factorisation of real matrix.
The QR factorization of the real matrix of size
is given by
where is an
orthogonal matrix and
is an
upper triangular matrix.
Returns true is calculation succeeds. False otherwise. Uses the LAPACK routine DGEQRF and DORGQR.
QR factorisation of real matrix with pivoting.
The QR factorization of the real matrix of size
is given by
where is an
orthogonal matrix,
is an
upper triangular matrix and
is an
permutation matrix.
Returns true is calculation succeeds. False otherwise. Uses the LAPACK routines DGEQP3 and DORGQR.
QR factorisation of a complex matrix.
The QR factorization of the complex matrix of size
is given by
where is an
unitary matrix and
is an
upper triangular matrix.
Returns true is calculation succeeds. False otherwise. Uses the LAPACK routines ZGEQRF and ZUNGQR.
QR factorisation of a complex matrix with pivoting.
The QR factorization of the complex matrix of size
is given by
where is an
unitary matrix,
is an
upper triangular matrix and
is an
permutation matrix.
Returns true is calculation succeeds. False otherwise. Uses the LAPACK routines ZGEQP3 and ZUNGQR.
Schur decomposition of a real matrix.
This function computes the Schur form of a square real matrix . The Schur decomposition satisfies the following equation:
where: is a unitary,
is upper quasi-triangular, and
is the transposed
matrix.
The upper quasi-triangular matrix may have blocks on its diagonal.
Uses the LAPACK routine DGEES.
Schur decomposition of a real matrix.
This function computes the Schur form of a square real matrix . The Schur decomposition satisfies the following equation:
where: is a unitary,
is upper quasi-triangular, and
is the transposed
matrix.
The upper quasi-triangular matrix may have blocks on its diagonal.
uses the LAPACK routine DGEES.
Schur decomposition of a complex matrix.
This function computes the Schur form of a square complex matrix . The Schur decomposition satisfies the following equation:
where: is a unitary,
is upper triangular, and
is the Hermitian transposition of the
matrix.
Uses the LAPACK routine ZGEES.
Schur decomposition of a complex matrix.
This function computes the Schur form of a square complex matrix . The Schur decomposition satisfies the following equation:
where: is a unitary,
is upper triangular, and
is the Hermitian transposition of the
matrix.
Uses the LAPACK routine ZGEES.
Get singular values s
of a real matrix A
using SVD.
This function calculates singular values from the SVD decomposition of a real matrix
. The SVD algorithm computes the decomposition of a real
matrix
so that
where are the singular values of
. Or put differently:
where
Referenced by orth().
Get singular values s
of a complex matrix A
using SVD.
This function calculates singular values from the SVD decomposition of a complex matrix
. The SVD algorithm computes the decomposition of a complex
matrix
so that
where are the singular values of
. Or put differently:
where
Return singular values of a real matrix A
using SVD.
This function returns singular values from the SVD decomposition of a real matrix . The SVD algorithm computes the decomposition of a real
matrix
so that
where are the singular values of
. Or put differently:
where
Return singular values of a complex matrix A
using SVD.
This function returns singular values from the SVD decomposition of a complex matrix . The SVD algorithm computes the decomposition of a complex
matrix
so that
where are the singular values of
. Or put differently:
where
Perform Singular Value Decomposition (SVD) of a real matrix A
.
This function returns two orthonormal matrices and
and a vector of singular values
. The SVD algorithm computes the decomposition of a real
matrix
so that
where the elements of ,
are the singular values of
. Or put differently:
where
Perform Singular Value Decomposition (SVD) of a complex matrix A
.
This function returns two orthonormal matrices and
and a vector of singular values
. The SVD algorithm computes the decomposition of a complex
matrix
so that
where the elements of ,
are the singular values of
. Or put differently:
where
Generated on Wed Jan 20 23:03:08 2010 for IT++ by Doxygen 1.6.2