![]() |
Octopus
|
Data Types | |
interface | lalg_cholesky |
interface | lalg_determinant |
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine. More... | |
interface | lalg_eigensolve |
interface | lalg_eigensolve_nonh |
interface | lalg_eigensolve_parallel |
interface | lalg_eigensolve_tridiagonal |
interface | lalg_geneigensolve |
interface | lalg_inverse |
interface | lalg_least_squares |
interface | lalg_linsyssolve |
interface | lalg_lowest_eigensolve |
interface | lalg_lowest_geneigensolve |
interface | lalg_matrix_function |
interface | lalg_matrix_norm2 |
interface | lalg_pseudo_inverse |
interface | lalg_singular_value_decomp |
interface | lalg_svd_inverse |
interface | lapack_geev |
Functions/Subroutines | |
real(real64) function | sfmin () |
Auxiliary function. More... | |
subroutine | lalg_dgeev (jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info) |
subroutine | lalg_zgeev (jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info) |
subroutine, public | zlalg_exp (nn, pp, aa, ex, hermitian) |
subroutine, public | zlalg_phi (nn, pp, aa, ex, hermitian) |
subroutine, public | lalg_zeigenderivatives (n, mat, zeigenvec, zeigenval, zmat) |
subroutine | lalg_zpseudoinverse (n, mat, imat) |
Computes the Moore-Penrose pseudoinverse of a complex matrix. More... | |
subroutine, public | lalg_check_zeigenderivatives (n, mat) |
complex(real64) function, public | lalg_zdni (eigenvec, alpha, beta) |
complex(real64) function, public | lalg_zduialpha (eigenvec, mmatrix, alpha, gamma, delta) |
complex(real64) function, public | lalg_zd2ni (eigenvec, mmatrix, alpha, beta, gamma, delta) |
pure real(real64) function | pseudoinverse_default_tolerance (m, n, sg_values) |
Computes the default Moore-Penrose pseudoinverse tolerance for zeroing. More... | |
real(real64) function, dimension(1:n, 1:n), public | lalg_remove_rotation (n, A) |
Remove rotation from affine transformation A by computing the polar decomposition and discarding the rotational part. The polar decomposition of A is given by A = U P with P = sqrt(A^T A), where U is a rotation matrix and P is a scaling matrix. This function returns P. More... | |
subroutine | zcholesky (n, a, bof, err_code) |
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a, dim(a) = n x n. On return a = u^T u with u upper triangular matrix. More... | |
subroutine | zgeneigensolve (n, a, b, e, preserve_mat, bof, err_code) |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form \( Ax=\lambda Bx \). B is also positive definite. More... | |
subroutine | zeigensolve_nonh (n, a, e, err_code, side, sort_eigenvectors) |
Computes all the eigenvalues and the right (left) eigenvectors of a real or complex (non-Hermitian) eigenproblem, of the form A*x=(lambda)*x. More... | |
subroutine | zlowest_geneigensolve (k, n, a, b, e, v, preserve_mat, bof, err_code) |
Computes the k lowest eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form A*x=(lambda)*B*x. B is also positive definite. More... | |
subroutine | zeigensolve (n, a, e, bof, err_code) |
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A. More... | |
subroutine | zeigensolve_tridiagonal (n, a, e, bof, err_code) |
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian case, the matrix is assumed to be real symmetric (even if stored in complex) More... | |
subroutine | zlowest_eigensolve (k, n, a, e, v, preserve_mat) |
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem, of the form A*x=(lambda)*x. Here A is assumed to be symmetric. More... | |
complex(real64) function | zdeterminant (n, a, preserve_mat) |
Invert a real symmetric or complex Hermitian square matrix a. More... | |
subroutine | zdirect_inverse (n, a, det) |
Invert a real symmetric or complex Hermitian square matrix a. More... | |
subroutine | zmatrix_norm2 (m, n, a, norm_l2, preserve_mat) |
Norm of a 2D matrix. More... | |
subroutine | zsym_inverse (uplo, n, a) |
Invert a real/complex symmetric square matrix a. More... | |
subroutine | zlinsyssolve (n, nrhs, a, b, x) |
compute the solution to a complex system of linear equations A*X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices. More... | |
subroutine | zsingular_value_decomp (m, n, a, u, vt, sg_values, preserve_mat) |
Computes the singular value decomposition of a complex MxN matrix a. More... | |
subroutine | zsvd_inverse (m, n, a, threshold) |
Computes inverse of a complex MxN matrix, a, using the SVD decomposition. More... | |
subroutine | zlalg_pseudo_inverse (a, threshold) |
Invert a matrix with the Moore-Penrose pseudo-inverse. More... | |
subroutine | zupper_triangular_inverse (n, a) |
Calculate the inverse of a real/complex upper triangular matrix (in unpacked storage). (lower triangular would be a trivial variant of this) More... | |
subroutine | zleast_squares_vec (nn, aa, bb, xx, preserve_mat) |
subroutine | zeigensolve_parallel (n, a, e, bof, err_code) |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenproblem in parallel using ScaLAPACK or ELPA on all processors n: dimension of matrix a: input matrix, on exit: contains eigenvectors e: eigenvalues. More... | |
subroutine | zinverse (n, a, method, det, threshold, uplo) |
An interface to different method to invert a matrix. More... | |
subroutine, public | zlalg_matrix_function (n, factor, a, fun_a, fun, hermitian, tridiagonal) |
This routine calculates a function of a matrix by using an eigenvalue decomposition. More... | |
subroutine | dcholesky (n, a, bof, err_code) |
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a, dim(a) = n x n. On return a = u^T u with u upper triangular matrix. More... | |
subroutine | dgeneigensolve (n, a, b, e, preserve_mat, bof, err_code) |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form \( Ax=\lambda Bx \). B is also positive definite. More... | |
subroutine | deigensolve_nonh (n, a, e, err_code, side, sort_eigenvectors) |
Computes all the eigenvalues and the right (left) eigenvectors of a real or complex (non-Hermitian) eigenproblem, of the form A*x=(lambda)*x. More... | |
subroutine | dlowest_geneigensolve (k, n, a, b, e, v, preserve_mat, bof, err_code) |
Computes the k lowest eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form A*x=(lambda)*B*x. B is also positive definite. More... | |
subroutine | deigensolve (n, a, e, bof, err_code) |
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A. More... | |
subroutine | deigensolve_tridiagonal (n, a, e, bof, err_code) |
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian case, the matrix is assumed to be real symmetric (even if stored in complex) More... | |
subroutine | dlowest_eigensolve (k, n, a, e, v, preserve_mat) |
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem, of the form A*x=(lambda)*x. Here A is assumed to be symmetric. More... | |
real(real64) function | ddeterminant (n, a, preserve_mat) |
Invert a real symmetric or complex Hermitian square matrix a. More... | |
subroutine | ddirect_inverse (n, a, det) |
Invert a real symmetric or complex Hermitian square matrix a. More... | |
subroutine | dmatrix_norm2 (m, n, a, norm_l2, preserve_mat) |
Norm of a 2D matrix. More... | |
subroutine | dsym_inverse (uplo, n, a) |
Invert a real/complex symmetric square matrix a. More... | |
subroutine | dlinsyssolve (n, nrhs, a, b, x) |
compute the solution to a real system of linear equations A*X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices. More... | |
subroutine | dsingular_value_decomp (m, n, a, u, vt, sg_values, preserve_mat) |
Computes the singular value decomposition of a real M x N matrix a. More... | |
subroutine | dsvd_inverse (m, n, a, threshold) |
Computes the inverse of a real M x N matrix, a, using the SVD decomposition. More... | |
subroutine | dlalg_pseudo_inverse (a, threshold) |
Invert a matrix with the Moore-Penrose pseudo-inverse. More... | |
subroutine | dupper_triangular_inverse (n, a) |
Calculate the inverse of a real/complex upper triangular matrix (in unpacked storage). (lower triangular would be a trivial variant of this) More... | |
subroutine | dleast_squares_vec (nn, aa, bb, xx, preserve_mat) |
subroutine | deigensolve_parallel (n, a, e, bof, err_code) |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenproblem in parallel using ScaLAPACK or ELPA on all processors n: dimension of matrix a: input matrix, on exit: contains eigenvectors e: eigenvalues. More... | |
subroutine | dinverse (n, a, method, det, threshold, uplo) |
An interface to different method to invert a matrix. More... | |
subroutine | dlalg_matrix_function (n, factor, a, fun_a, fun, hermitian, tridiagonal) |
This routine calculates a function of a matrix by using an eigenvalue decomposition. More... | |
|
private |
Auxiliary function.
Definition at line 246 of file lalg_adv.F90.
|
private |
[in,out] | a | a(lda,n) |
[out] | w | w(n) |
[out] | vr | vl(ldvl,n), vl(ldvr,n) |
[out] | rwork | rwork(max(1,2n)) |
[out] | work | work(lwork) |
Definition at line 258 of file lalg_adv.F90.
|
private |
[in,out] | a | a(lda,n) |
[out] | w | w(n) |
[out] | vr | vl(ldvl,n), vl(ldvr,n) |
[out] | rwork | rwork(max(1,2n)) |
[out] | work | work(lwork) |
Definition at line 285 of file lalg_adv.F90.
subroutine, public lalg_adv_oct_m::zlalg_exp | ( | integer, intent(in) | nn, |
complex(real64), intent(in) | pp, | ||
complex(real64), dimension(:, :), intent(in) | aa, | ||
complex(real64), dimension(:, :), intent(inout) | ex, | ||
logical, intent(in) | hermitian | ||
) |
This routine calculates the exponential of a matrix by using an eigenvalue decomposition.
For the hermitian case:
A = V D V^T => exp(A) = V exp(D) V^T
and in general
A = V D V^-1 => exp(A) = V exp(D) V^-1
This is slow but it is simple to implement, and for the moment it does not affect performance.
Definition at line 319 of file lalg_adv.F90.
subroutine, public lalg_adv_oct_m::zlalg_phi | ( | integer, intent(in) | nn, |
complex(real64), intent(in) | pp, | ||
complex(real64), dimension(:, :), intent(in) | aa, | ||
complex(real64), dimension(:, :), intent(inout) | ex, | ||
logical, intent(in) | hermitian | ||
) |
This routine calculates phi(pp*A), where A is a matrix, pp is any complex number, and phi is the function:
phi(x) = (e^x - 1)/x
For the Hermitian case, for any function f:
A = V D V^T => f(A) = V f(D) V^T
and in general
A = V D V^-1 => f(A) = V f(D) V^-1
Definition at line 397 of file lalg_adv.F90.
subroutine, public lalg_adv_oct_m::lalg_zeigenderivatives | ( | integer, intent(in) | n, |
complex(real64), dimension(:, :), intent(in), contiguous | mat, | ||
complex(real64), dimension(:, :), intent(out), contiguous | zeigenvec, | ||
complex(real64), dimension(:), intent(out), contiguous | zeigenval, | ||
complex(real64), dimension(:, :, :), intent(out), contiguous | zmat | ||
) |
Computes the necessary ingredients to obtain, later, the first and second derivatives of the eigenvalues of a Hermitean complex matrix zmat, and the first derivatives of the eigenvectors.
This follows the scheme of J. R. Magnus, Econometric Theory 1, 179 (1985), restricted to Hermitean matrices, although probably this can be
Definition at line 472 of file lalg_adv.F90.
|
private |
Computes the Moore-Penrose pseudoinverse of a complex matrix.
Definition at line 524 of file lalg_adv.F90.
subroutine, public lalg_adv_oct_m::lalg_check_zeigenderivatives | ( | integer, intent(in) | n, |
complex(real64), dimension(:, :), intent(in) | mat | ||
) |
The purpose of this routine is to check that "lalg_zeigenderivatives" is working properly, and therefore, it is not really called anywhere in the code. It is here only for debugging purposes (perhaps it will
Definition at line 583 of file lalg_adv.F90.
complex(real64) function, public lalg_adv_oct_m::lalg_zdni | ( | complex(real64), dimension(2) | eigenvec, |
integer, intent(in) | alpha, | ||
integer, intent(in) | beta | ||
) |
Definition at line 701 of file lalg_adv.F90.
complex(real64) function, public lalg_adv_oct_m::lalg_zduialpha | ( | complex(real64), dimension(2) | eigenvec, |
complex(real64), dimension(2, 2) | mmatrix, | ||
integer, intent(in) | alpha, | ||
integer, intent(in) | gamma, | ||
integer, intent(in) | delta | ||
) |
Definition at line 707 of file lalg_adv.F90.
complex(real64) function, public lalg_adv_oct_m::lalg_zd2ni | ( | complex(real64), dimension(2) | eigenvec, |
complex(real64), dimension(2, 2) | mmatrix, | ||
integer, intent(in) | alpha, | ||
integer, intent(in) | beta, | ||
integer, intent(in) | gamma, | ||
integer, intent(in) | delta | ||
) |
Definition at line 713 of file lalg_adv.F90.
|
private |
Computes the default Moore-Penrose pseudoinverse tolerance for zeroing.
We use here the value suggested here https:
Definition at line 725 of file lalg_adv.F90.
real(real64) function, dimension(1:n, 1:n), public lalg_adv_oct_m::lalg_remove_rotation | ( | integer, intent(in) | n, |
real(real64), dimension(1:n, 1:n), intent(in) | A | ||
) |
Remove rotation from affine transformation A by computing the polar decomposition and discarding the rotational part. The polar decomposition of A is given by A = U P with P = sqrt(A^T A), where U is a rotation matrix and P is a scaling matrix. This function returns P.
Definition at line 737 of file lalg_adv.F90.
|
private |
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a, dim(a) = n x n. On return a = u^T u with u upper triangular matrix.
[in,out] | a | (n,n) |
[in,out] | bof | Bomb on failure. |
Definition at line 822 of file lalg_adv.F90.
|
private |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form \( Ax=\lambda Bx \). B is also positive definite.
For optimal performances, this uses the divide and conquer algoritm
[in,out] | a | (n,n) |
[in,out] | b | (n,n) |
[out] | e | (n) |
[in] | preserve_mat | If true, the matrix a and b on exit are the same |
[in,out] | bof | Bomb on failure. |
Definition at line 877 of file lalg_adv.F90.
|
private |
Computes all the eigenvalues and the right (left) eigenvectors of a real or complex (non-Hermitian) eigenproblem, of the form A*x=(lambda)*x.
[in,out] | a | (n,n) |
[out] | e | (n) |
[in] | side | which eigenvectors ('L' or 'R') |
[in] | sort_eigenvectors | only applies to complex version, sorts by real part |
Definition at line 988 of file lalg_adv.F90.
|
private |
Computes the k lowest eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form A*x=(lambda)*B*x. B is also positive definite.
[in,out] | a | (n, n) |
[in,out] | b | (n, n) |
[out] | e | (n) |
[out] | v | (n, n) |
[in] | preserve_mat | If true, the matrix a and b on exit are the same |
[in,out] | bof | Bomb on failure. |
Definition at line 1107 of file lalg_adv.F90.
|
private |
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A.
[in,out] | a | (n,n) |
[out] | e | (n) |
[in,out] | bof | Bomb on failure. |
Definition at line 1226 of file lalg_adv.F90.
|
private |
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian case, the matrix is assumed to be real symmetric (even if stored in complex)
[in,out] | a | (n,n) |
[out] | e | (n) |
[in,out] | bof | Bomb on failure. |
Definition at line 1296 of file lalg_adv.F90.
|
private |
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem, of the form A*x=(lambda)*x. Here A is assumed to be symmetric.
[in] | k | Number of eigenvalues requested |
[in] | n | Dimensions of a |
[in,out] | a | (n, n) |
[out] | e | (n) The first k elements contain the selected eigenvalues in ascending order. |
[out] | v | (n, k) |
[in] | preserve_mat | If true, the matrix a and b on exit are the same |
Definition at line 1375 of file lalg_adv.F90.
|
private |
Invert a real symmetric or complex Hermitian square matrix a.
[in,out] | a | (n,n) |
Definition at line 1466 of file lalg_adv.F90.
|
private |
Invert a real symmetric or complex Hermitian square matrix a.
[in,out] | a | (n,n) |
Definition at line 1513 of file lalg_adv.F90.
|
private |
Norm of a 2D matrix.
The spectral norm of a matrix \(A\) is the largest singular value of \(A\). i.e., the largest eigenvalue of the matrix \(\sqrt{\dagger{A}A}\).
[in] | m,n | Dimensions of A |
[in] | a | 2D matrix |
[out] | norm2 | L2 norm of A |
Definition at line 1581 of file lalg_adv.F90.
|
private |
Invert a real/complex symmetric square matrix a.
[in,out] | a | (n,n) |
Definition at line 1615 of file lalg_adv.F90.
|
private |
compute the solution to a complex system of linear equations A*X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
[in,out] | a | (n, n) |
[in,out] | b | (n, nrhs) |
[out] | x | (n, nrhs) |
Definition at line 1666 of file lalg_adv.F90.
|
private |
Computes the singular value decomposition of a complex MxN matrix a.
[in,out] | a | (m,n) |
[out] | vt | (n,n) and (m,m) |
[out] | sg_values | (n) |
Definition at line 1744 of file lalg_adv.F90.
|
private |
Computes inverse of a complex MxN matrix, a, using the SVD decomposition.
[in,out] | a | Input (m,n) |
Definition at line 1827 of file lalg_adv.F90.
|
private |
Invert a matrix with the Moore-Penrose pseudo-inverse.
SVD is used to find the U, V and Sigma matrices:
\[ A = U \Sigma V^\dagger \]
Diagonal terms in Sigma <= threshold
are set to zero, and the inverse is constructed as:
\[] A^{-1} \approx V \Sigma U^\dagger \]
[in,out] | a | Input: (m, n) |
Definition at line 1893 of file lalg_adv.F90.
|
private |
Calculate the inverse of a real/complex upper triangular matrix (in unpacked storage). (lower triangular would be a trivial variant of this)
[in,out] | a | (n,n) |
Definition at line 1957 of file lalg_adv.F90.
|
private |
Definition at line 2005 of file lalg_adv.F90.
|
private |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenproblem in parallel using ScaLAPACK or ELPA on all processors n: dimension of matrix a: input matrix, on exit: contains eigenvectors e: eigenvalues.
[in,out] | a | (n,n) |
[out] | e | (n) |
[in,out] | bof | Bomb on failure. |
Definition at line 2080 of file lalg_adv.F90.
|
private |
An interface to different method to invert a matrix.
The possible methods are: svd, dir, sym, upp For the SVD, an optional argument threshold an be specified For the direct inverse, an optional output determinant can be obtained For the symmetric matrix case, the optional argument uplo must be specified
[in,out] | a | (n,n) |
[out] | det | Determinant of the matrix. Direct inversion only |
[in] | threshold | Threshold for the SVD pseudoinverse |
[in] | uplo | Is the symmetric matrix stored in the upper or lower part? |
Definition at line 2213 of file lalg_adv.F90.
subroutine, public lalg_adv_oct_m::zlalg_matrix_function | ( | integer, intent(in) | n, |
complex(real64), intent(in) | factor, | ||
complex(real64), dimension(:, :), intent(in) | a, | ||
complex(real64), dimension(:, :), intent(inout) | fun_a, | ||
fun, | |||
logical, intent(in) | hermitian, | ||
logical, intent(in), optional | tridiagonal | ||
) |
This routine calculates a function of a matrix by using an eigenvalue decomposition.
For the hermitian case:
\[ A = V D V^T \implies fun(A) = V fun(D) V^T \]
and in general
\[ A = V D V^-1 \implies fun(A) = V fun(D) V^-1 \]
where \(V\) are the eigenvectors, and \(D\) is a diagonal matrix containing the eigenvalues.
In addition, this function can compute \(fun(factor*A)\) for a complex factor.
This is slow but it is simple to implement, and for the moment it does not affect performance.
[in] | n | dimension of the matrix A |
[in] | factor | complex factor |
[in] | a | matrix A |
[in,out] | fun_a | fun(A) |
[in] | hermitian | is the matrix hermitian? |
[in] | tridiagonal | is the matrix tridiagonal? |
Definition at line 2261 of file lalg_adv.F90.
|
private |
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a, dim(a) = n x n. On return a = u^T u with u upper triangular matrix.
[in,out] | a | (n,n) |
[in,out] | bof | Bomb on failure. |
Definition at line 2443 of file lalg_adv.F90.
|
private |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form \( Ax=\lambda Bx \). B is also positive definite.
For optimal performances, this uses the divide and conquer algoritm
[in,out] | a | (n,n) |
[in,out] | b | (n,n) |
[out] | e | (n) |
[in] | preserve_mat | If true, the matrix a and b on exit are the same |
[in,out] | bof | Bomb on failure. |
Definition at line 2498 of file lalg_adv.F90.
|
private |
Computes all the eigenvalues and the right (left) eigenvectors of a real or complex (non-Hermitian) eigenproblem, of the form A*x=(lambda)*x.
[in,out] | a | (n,n) |
[out] | e | (n) |
[in] | side | which eigenvectors ('L' or 'R') |
[in] | sort_eigenvectors | only applies to complex version, sorts by real part |
Definition at line 2604 of file lalg_adv.F90.
|
private |
Computes the k lowest eigenvalues and the eigenvectors of a real symmetric or complex Hermitian generalized definite eigenproblem, of the form A*x=(lambda)*B*x. B is also positive definite.
[in,out] | a | (n, n) |
[in,out] | b | (n, n) |
[out] | e | (n) |
[out] | v | (n, n) |
[in] | preserve_mat | If true, the matrix a and b on exit are the same |
[in,out] | bof | Bomb on failure. |
Definition at line 2723 of file lalg_adv.F90.
|
private |
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A.
[in,out] | a | (n,n) |
[out] | e | (n) |
[in,out] | bof | Bomb on failure. |
Definition at line 2843 of file lalg_adv.F90.
|
private |
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian case, the matrix is assumed to be real symmetric (even if stored in complex)
[in,out] | a | (n,n) |
[out] | e | (n) |
[in,out] | bof | Bomb on failure. |
Definition at line 2913 of file lalg_adv.F90.
|
private |
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem, of the form A*x=(lambda)*x. Here A is assumed to be symmetric.
[in] | k | Number of eigenvalues requested |
[in] | n | Dimensions of a |
[in,out] | a | (n, n) |
[out] | e | (n) The first k elements contain the selected eigenvalues in ascending order. |
[out] | v | (n, k) |
[in] | preserve_mat | If true, the matrix a and b on exit are the same |
Definition at line 2992 of file lalg_adv.F90.
|
private |
Invert a real symmetric or complex Hermitian square matrix a.
[in,out] | a | (n,n) |
Definition at line 3081 of file lalg_adv.F90.
|
private |
Invert a real symmetric or complex Hermitian square matrix a.
[in,out] | a | (n,n) |
Definition at line 3128 of file lalg_adv.F90.
|
private |
Norm of a 2D matrix.
The spectral norm of a matrix \(A\) is the largest singular value of \(A\). i.e., the largest eigenvalue of the matrix \(\sqrt{\dagger{A}A}\).
[in] | m,n | Dimensions of A |
[in] | a | 2D matrix |
[out] | norm2 | L2 norm of A |
Definition at line 3196 of file lalg_adv.F90.
|
private |
Invert a real/complex symmetric square matrix a.
[in,out] | a | (n,n) |
Definition at line 3230 of file lalg_adv.F90.
|
private |
compute the solution to a real system of linear equations A*X = B, where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
[in,out] | a | (n, n) |
[in,out] | b | (n, nrhs) |
[out] | x | (n, nrhs) |
Definition at line 3281 of file lalg_adv.F90.
|
private |
Computes the singular value decomposition of a real M x N matrix a.
[in,out] | a | (m,n) |
[out] | vt | (m,m) (n,n) |
[out] | sg_values | (min(m,n)) |
Definition at line 3359 of file lalg_adv.F90.
|
private |
Computes the inverse of a real M x N matrix, a, using the SVD decomposition.
[in,out] | a | Input (m,n) |
Definition at line 3438 of file lalg_adv.F90.
|
private |
Invert a matrix with the Moore-Penrose pseudo-inverse.
SVD is used to find the U, V and Sigma matrices:
\[ A = U \Sigma V^\dagger \]
Diagonal terms in Sigma <= threshold
are set to zero, and the inverse is constructed as:
\[] A^{-1} \approx V \Sigma U^\dagger \]
[in,out] | a | Input: (m, n) |
Definition at line 3504 of file lalg_adv.F90.
|
private |
Calculate the inverse of a real/complex upper triangular matrix (in unpacked storage). (lower triangular would be a trivial variant of this)
[in,out] | a | (n,n) |
Definition at line 3568 of file lalg_adv.F90.
|
private |
Definition at line 3616 of file lalg_adv.F90.
|
private |
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenproblem in parallel using ScaLAPACK or ELPA on all processors n: dimension of matrix a: input matrix, on exit: contains eigenvectors e: eigenvalues.
[in,out] | a | (n,n) |
[out] | e | (n) |
[in,out] | bof | Bomb on failure. |
Definition at line 3688 of file lalg_adv.F90.
|
private |
An interface to different method to invert a matrix.
The possible methods are: svd, dir, sym, upp For the SVD, an optional argument threshold an be specified For the direct inverse, an optional output determinant can be obtained For the symmetric matrix case, the optional argument uplo must be specified
[in,out] | a | (n,n) |
[out] | det | Determinant of the matrix. Direct inversion only |
[in] | threshold | Threshold for the SVD pseudoinverse |
[in] | uplo | Is the symmetric matrix stored in the upper or lower part? |
Definition at line 3821 of file lalg_adv.F90.
|
private |
This routine calculates a function of a matrix by using an eigenvalue decomposition.
For the hermitian case:
\[ A = V D V^T \implies fun(A) = V fun(D) V^T \]
and in general
\[ A = V D V^-1 \implies fun(A) = V fun(D) V^-1 \]
where \(V\) are the eigenvectors, and \(D\) is a diagonal matrix containing the eigenvalues.
In addition, this function can compute \(fun(factor*A)\) for a complex factor.
This is slow but it is simple to implement, and for the moment it does not affect performance.
[in] | n | dimension of the matrix A |
[in] | factor | complex factor |
[in] | a | matrix A |
[in,out] | fun_a | fun(A) |
[in] | hermitian | is the matrix hermitian? |
[in] | tridiagonal | is the matrix tridiagonal? |
Definition at line 3869 of file lalg_adv.F90.