Singular and eigen values¶
This section regroups all of the linear algebra functions concerning singular and eigen values which are currently implemented within castor. In order to access these functions, the user must include castor/linalg.hpp
. Two levels of interface are available. The high-level interface provides functions with simplified arguments while the low-level interface is much closer to the BLAS/LAPACK API.
The user will find here high-level linear algebra functions. Some examples of use may also be found at Linear algebra.
eig¶
-
auto castor::eig(matrix<std::complex<double>> const &A, matrix<std::complex<double>> const &B, std::string typ)¶
-
auto castor::eig(matrix<std::complex<float>> const &A, matrix<std::complex<float>> const &B, std::string typ)¶
-
auto castor::eig(matrix<float> const &A, std::string typ)¶
Eigenvalues and eigenvectors.
E = eig(A) produces a row vector E containing the eigenvalues of a square matrix A.
[E,U] = eig(A,”left”) produces a row vector E containing the eigenvalues of a square matrix A, associated to the left eigenvectors U as: U^t*A = diag(E)*U^t, where U^t is the conjugate transpose.
[E,V] = eig(A,”right”) produces a row vector E containing the eigenvalues of a square matrix A, associated to the right eigenvector V as: A*V = V*diag(E).
[E,V] = eig(A,”none”) is equivalent to E = eig(A).
E = eig(A,B) produces a row vector E containing the generalized eigenvalues of square matrices matrices A and B
[E,U] = eig(A,B,”left”) produces a row vector and a matrix U containing the generalized eigenvalues and the corresponding left eigenvectors such that U^t*A = diag(E)*U^t*B where U^t is the conjugate transpose.
[E,V] = eig(A,B,”right”) produces a row vector and a matrix V containing the generalized eigenvalues and the corresponding right eigenvectors such that A*V = B*V*diag(E)
[E,V] = eig(A,B,”none”) is equivalent to E = eig(A,B)
matrix<> A = 1-eye(4); matrix<std::complex<double>> E, V; std::tie(E,V) = eig(A,"right"); disp(E); disp(mtimes(A,V)-mtimes(V,diag(E)));
qrsvd¶
-
template<typename T>
auto castor::qrsvd(matrix<T> const &A, matrix<T> const &B, float tol)¶ Low-rank recompression.
[U,V] = qrsvd(A,B,TOL) provides a low-rank recompression using QR factorization and SVD.
matrix<> A = rand(5,3), B = rand(3,6), U, V; std::tie(U,V) = qrsvd(A,B,1e-3); disp(mtimes(A,B)-mtimes(U,V));
rank¶
Warning
doxygenfunction: Unable to resolve function “rank” with arguments () in doxygen xml output for project “castor” from directory: ../xml. Potential matches:
- template<typename T> std::size_t rank(matrix<T> const &A, float tol = 0)
svd¶
-
auto castor::svd(matrix<float> const &A, std::string typ)¶
Singular value decomposition.
S = svd(A) returns a row vector S containing the singular values of A, with nonnegative elements sorted in decreasing order.
[S,U,Vt] = svd(A,”vect”) returns singulars values and unitary matrices U and Vt (not V) so that : A = U*diag(S)*V^t
[U,S,Vt] = svd(A,”none”) is equivalent to S = svd(A).
matrix<> A = rand(3,3); matrix<> U, S, Vt; std::tie(S,U,Vt) = svd(A,"vect"); disp(mtimes(mtimes(U,diag(S)),Vt)-A);