27 use,
intrinsic :: iso_fortran_env
149 real(real64) function sfmin()
151 real(real64) function dlamch(cmach)
154 character(1),
intent(in) :: cmach
161 subroutine lalg_dgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
162 character(1),
intent(in) :: jobvl, jobvr
163 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
164 real(real64),
intent(inout) :: a(:, :)
165 complex(real64),
intent(out) :: w(:)
166 real(real64),
intent(out) :: vl(:, :), vr(:, :)
167 real(real64),
intent(out) :: rwork(:)
168 real(real64),
intent(out) :: work(:)
169 integer,
intent(out) :: info
171 real(real64),
allocatable :: wr(:), wi(:)
175 safe_allocate(wr(1:n))
176 safe_allocate(wi(1:n))
178 call dgeev(jobvl, jobvr, n, a(1, 1), lda, wr(1), wi(1), vl(1, 1), ldvl, vr(1, 1), ldvr, work(1), lwork, rwork(1), info)
180 w(1:n) = cmplx(wr(1:n), wi(1:n), real64)
182 safe_deallocate_a(wr)
183 safe_deallocate_a(wi)
188 subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
189 character(1),
intent(in) :: jobvl, jobvr
190 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
191 complex(real64),
intent(inout) :: a(:, :)
192 complex(real64),
intent(out) :: w(:)
193 complex(real64),
intent(out) :: vl(:, :), vr(:, :)
194 real(real64),
intent(out) :: rwork(:)
195 complex(real64),
intent(out) :: work(:)
196 integer,
intent(out) :: info
200 call zgeev(jobvl, jobvr, n, a(1, 1), lda, w(1), vl(1, 1), ldvl, vr(1, 1), ldvr, work(1), lwork, rwork(1), info)
222 subroutine zlalg_exp(nn, pp, aa, ex, hermitian)
223 integer,
intent(in) :: nn
224 complex(real64),
intent(in) :: pp
225 complex(real64),
intent(in) :: aa(:, :)
226 complex(real64),
intent(inout) :: ex(:, :)
227 logical,
intent(in) :: hermitian
229 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
230 real(real64),
allocatable :: evalues(:)
236 safe_allocate(evectors(1:nn, 1:nn))
239 safe_allocate(evalues(1:nn))
240 safe_allocate(zevalues(1:nn))
242 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
246 zevalues(1:nn) =
exp(pp*evalues(1:nn))
249 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
252 ex(:, :) = matmul(evectors(:, :), ex(:, :))
254 safe_deallocate_a(evalues)
255 safe_deallocate_a(zevalues)
257 safe_allocate(zevalues(1:nn))
259 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
263 zevalues(1:nn) =
exp(pp*zevalues(1:nn))
265 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
270 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
273 ex(:, :) = matmul(ex(:, :), evectors(:, :))
275 safe_deallocate_a(zevalues)
278 safe_deallocate_a(evectors)
300 subroutine zlalg_phi(nn, pp, aa, ex, hermitian)
301 integer,
intent(in) :: nn
302 complex(real64),
intent(in) :: pp
303 complex(real64),
intent(in) :: aa(:, :)
304 complex(real64),
intent(inout) :: ex(:, :)
305 logical,
intent(in) :: hermitian
307 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
308 real(real64),
allocatable :: evalues(:)
314 safe_allocate(evectors(1:nn, 1:nn))
317 safe_allocate(evalues(1:nn))
318 safe_allocate(zevalues(1:nn))
320 evectors(:, :) = aa(:, :)
325 zevalues(ii) = (
exp(pp*evalues(ii)) - m_z1) / (pp*evalues(ii))
329 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
332 ex(:, :) = matmul(evectors(:, :), ex(:, :))
334 safe_deallocate_a(evalues)
335 safe_deallocate_a(zevalues)
337 safe_allocate(zevalues(1:nn))
339 evectors(:, :) = aa(:, :)
344 zevalues(ii) = (
exp(pp*zevalues(ii)) - m_z1) / (pp*zevalues(ii))
347 ex(:, :) = evectors(:, :)
352 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
355 ex(:, :) = matmul(ex(:, :), evectors(:, :))
357 safe_deallocate_a(zevalues)
363 complex(real64) function lalg_zdni(eigenvec, alpha, beta)
364 integer,
intent(in) :: alpha, beta
365 complex(real64) :: eigenvec(2)
366 lalg_zdni = conjg(eigenvec(alpha)) * eigenvec(beta)
369 complex(real64) function lalg_zduialpha(eigenvec, mmatrix, alpha, gamma, delta)
370 integer,
intent(in) :: alpha, gamma, delta
371 complex(real64) :: eigenvec(2), mmatrix(2, 2)
375 complex(real64) function lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
376 integer,
intent(in) :: alpha, beta, gamma, delta
377 complex(real64) :: eigenvec(2), mmatrix(2, 2)
378 lalg_zd2ni = conjg(mmatrix(alpha, delta) * eigenvec(gamma)) * eigenvec(beta) + &
379 conjg(eigenvec(alpha)) * mmatrix(beta, gamma) * eigenvec(delta)
388 integer,
intent(in) :: m, n
389 real(real64),
intent(in) :: sg_values(:)
390 tol = m_epsilon * m * n * maxval(sg_values)
400 integer,
intent(in) :: n
401 real(real64),
intent(in) :: a(1:n, 1:n)
402 real(real64) :: p(1:n, 1:n)
404 real(real64) :: ata(1:n, 1:n)
408 ata = matmul(transpose(a), a)
415#include "complex.F90"
416#include "lalg_adv_lapack_inc.F90"
420#include "lalg_adv_lapack_inc.F90"
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
double exp(double __x) __attribute__((__nothrow__
This module contains interfaces for BLACS routines Interfaces are from http:
This module provides the BLACS processor grid.
This module contains interfaces for BLAS routines You should not use these routines directly....
complex(real64) function, public lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
integer function zmatrix_rank_svd(a, preserve_mat, tol)
Compute the rank of the matrix A using SVD.
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 gener...
complex(real64) function zdeterminant(n, a, preserve_mat)
Invert a real symmetric or complex Hermitian square matrix a.
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) e...
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...
subroutine, public zlalg_exp(nn, pp, aa, ex, hermitian)
subroutine deigensolve_parallel(n, a, e, bof, err_code)
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenprobl...
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 generalize...
subroutine zlalg_pseudo_inverse(a, threshold)
Invert a matrix with the Moore-Penrose pseudo-inverse.
subroutine zeigensolve(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A.
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.
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) e...
subroutine deigensolve_tridiagonal(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian c...
subroutine zinverse(n, a, method, det, threshold, uplo)
An interface to different method to invert a matrix.
real(real64) function sfmin()
Auxiliary function.
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 ...
subroutine zsingular_value_decomp(m, n, a, u, vt, sg_values, preserve_mat)
Computes the singular value decomposition of a complex MxN matrix a.
subroutine dcholesky(n, a, bof, err_code)
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a,...
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.
subroutine deigensolve(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A.
subroutine zeigensolve_tridiagonal(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian c...
subroutine zlowest_eigensolve(k, n, a, e, v, preserve_mat)
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem,...
complex(real64) function, public lalg_zduialpha(eigenvec, mmatrix, alpha, gamma, delta)
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 generalize...
subroutine dleast_squares_vec(nn, aa, bb, xx, preserve_mat)
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 gener...
subroutine zcholesky(n, a, bof, err_code)
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a,...
complex(real64) function, public lalg_zdni(eigenvec, alpha, beta)
pure real(real64) function pseudoinverse_default_tolerance(m, n, sg_values)
Computes the default Moore-Penrose pseudoinverse tolerance for zeroing.
subroutine zleast_squares_vec(nn, aa, bb, xx, preserve_mat)
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 ...
subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
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.
integer function dmatrix_rank_svd(a, preserve_mat, tol)
Compute the rank of the matrix A using SVD.
subroutine dinverse(n, a, method, det, threshold, uplo)
An interface to different method to invert a matrix.
subroutine zeigensolve_parallel(n, a, e, bof, err_code)
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenprobl...
real(real64) function ddeterminant(n, a, preserve_mat)
Invert a real symmetric or complex Hermitian square matrix a.
subroutine lalg_dgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
subroutine dlowest_eigensolve(k, n, a, e, v, preserve_mat)
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem,...
subroutine dsvd_inverse(m, n, a, threshold)
Computes the inverse of a real M x N matrix, a, using the SVD decomposition.
subroutine, public zlalg_phi(nn, pp, aa, ex, hermitian)
subroutine dlalg_pseudo_inverse(a, threshold)
Invert a matrix with the Moore-Penrose pseudo-inverse.
subroutine zsvd_inverse(m, n, a, threshold)
Computes inverse of a complex MxN matrix, a, using the SVD decomposition.
This module contains interfaces for LAPACK routines.
This module is intended to contain "only mathematical" functions and procedures.
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
This module is intended to contain "only mathematical" functions and procedures.
This module is intended to contain simple general-purpose utility functions and procedures.