26 use,
intrinsic :: iso_fortran_env
136 real(real64) function sfmin()
138 real(real64) function dlamch(cmach)
141 character(1),
intent(in) :: cmach
148 subroutine lalg_dgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
149 character(1),
intent(in) :: jobvl, jobvr
150 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
151 real(real64),
intent(inout) :: a(:,:)
152 complex(real64),
intent(out) :: w(:)
153 real(real64),
intent(out) :: vl(:,:), vr(:,:)
154 real(real64),
intent(out) :: rwork(:)
155 real(real64),
intent(out) :: work(:)
156 integer,
intent(out) :: info
158 real(real64),
allocatable :: wr(:), wi(:)
162 safe_allocate(wr(1:n))
163 safe_allocate(wi(1:n))
165 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)
167 w(1:n) = cmplx(wr(1:n), wi(1:n), real64)
169 safe_deallocate_a(wr)
170 safe_deallocate_a(wi)
175 subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
176 character(1),
intent(in) :: jobvl, jobvr
177 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
178 complex(real64),
intent(inout) :: a(:,:)
179 complex(real64),
intent(out) :: w(:)
180 complex(real64),
intent(out) :: vl(:,:), vr(:,:)
181 real(real64),
intent(out) :: rwork(:)
182 complex(real64),
intent(out) :: work(:)
183 integer,
intent(out) :: info
187 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)
209 integer,
intent(in) :: n
210 complex(real64),
intent(in) :: factor
211 complex(real64),
intent(in) :: a(:, :)
212 complex(real64),
intent(inout) :: fun_a(:, :)
214 complex(real64) function fun(z)
216 complex(real64),
intent(in) :: z
219 logical,
intent(in) :: hermitian
221 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
222 real(real64),
allocatable :: evalues(:)
228 assert(
size(a, 1) >= n)
230 safe_allocate(evectors(1:n, 1:n))
233 safe_allocate(evalues(1:n))
234 safe_allocate(zevalues(1:n))
236 evectors(1:n, 1:n) = a(1:n, 1:n)
241 zevalues(i) = fun(factor*evalues(i))
245 fun_a(1:n, i) = zevalues(1:n)*conjg(evectors(i, 1:n))
248 fun_a(1:n, 1:n) = matmul(evectors(1:n, 1:n), fun_a(1:n, 1:n))
250 safe_deallocate_a(evalues)
251 safe_deallocate_a(zevalues)
253 safe_allocate(zevalues(1:n))
255 evectors(1:n, 1:n) = a(1:n, 1:n)
260 zevalues(i) = fun(factor*zevalues(i))
263 fun_a(1:n, 1:n) = evectors(1:n, 1:n)
268 evectors(1:n, i) = zevalues(1:n)*evectors(1:n, i)
271 fun_a(1:n, 1:n) = matmul(fun_a(1:n, 1:n), evectors(1:n, 1:n))
273 safe_deallocate_a(zevalues)
276 safe_deallocate_a(evectors)
299 subroutine zlalg_exp(nn, pp, aa, ex, hermitian)
300 integer,
intent(in) :: nn
301 complex(real64),
intent(in) :: pp
302 complex(real64),
intent(in) :: aa(:, :)
303 complex(real64),
intent(inout) :: ex(:, :)
304 logical,
intent(in) :: hermitian
306 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
307 real(real64),
allocatable :: evalues(:)
313 safe_allocate(evectors(1:nn, 1:nn))
316 safe_allocate(evalues(1:nn))
317 safe_allocate(zevalues(1:nn))
319 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
323 zevalues(1:nn) =
exp(pp*evalues(1:nn))
326 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
329 ex(1:nn, 1:nn) = matmul(evectors(1:nn, 1:nn), ex(1:nn, 1:nn))
331 safe_deallocate_a(evalues)
332 safe_deallocate_a(zevalues)
334 safe_allocate(zevalues(1:nn))
336 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
340 zevalues(1:nn) =
exp(pp*zevalues(1:nn))
342 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
347 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
350 ex(1:nn, 1:nn) = matmul(ex(1:nn, 1:nn), evectors(1:nn, 1:nn))
352 safe_deallocate_a(zevalues)
355 safe_deallocate_a(evectors)
377 subroutine zlalg_phi(nn, pp, aa, ex, hermitian)
378 integer,
intent(in) :: nn
379 complex(real64),
intent(in) :: pp
380 complex(real64),
intent(in) :: aa(:, :)
381 complex(real64),
intent(inout) :: ex(:, :)
382 logical,
intent(in) :: hermitian
384 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
385 real(real64),
allocatable :: evalues(:)
391 safe_allocate(evectors(1:nn, 1:nn))
394 safe_allocate(evalues(1:nn))
395 safe_allocate(zevalues(1:nn))
402 zevalues(ii) = (
exp(pp*evalues(ii)) - m_z1) / (pp*evalues(ii))
406 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
409 ex(1:nn, 1:nn) = matmul(evectors(1:nn, 1:nn), ex(1:nn, 1:nn))
411 safe_deallocate_a(evalues)
412 safe_deallocate_a(zevalues)
414 safe_allocate(zevalues(1:nn))
416 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
421 zevalues(ii) = (
exp(pp*zevalues(ii)) - m_z1) / (pp*zevalues(ii))
424 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
429 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
432 ex(1:nn, 1:nn) = matmul(ex(1:nn, 1:nn), evectors(1:nn, 1:nn))
434 safe_deallocate_a(zevalues)
453 integer,
intent(in) :: n
454 complex(real64),
intent(in) :: mat(:, :)
455 complex(real64),
intent(out) :: zeigenvec(:, :)
456 complex(real64),
intent(out) :: zeigenval(:)
457 complex(real64),
intent(out) :: zmat(:, :, :)
459 integer :: i, alpha, beta
460 complex(real64),
allocatable :: knaught(:, :)
461 complex(real64),
allocatable :: lambdaminusdm(:, :)
462 complex(real64),
allocatable :: ilambdaminusdm(:, :)
463 complex(real64),
allocatable :: unit(:, :)
467 safe_allocate(unit(1:n, 1:n))
468 safe_allocate(knaught(1:n, 1:n))
469 safe_allocate(lambdaminusdm(1:n, 1:n))
470 safe_allocate(ilambdaminusdm(1:n, 1:n))
483 knaught(alpha, beta) = - zeigenvec(alpha, i)*conjg(zeigenvec(beta, i))
486 knaught = knaught + unit
487 lambdaminusdm = zeigenval(i)*unit - mat
489 zmat(:, :, i) = matmul(ilambdaminusdm, knaught)
492 safe_deallocate_a(unit)
493 safe_deallocate_a(knaught)
494 safe_deallocate_a(lambdaminusdm)
495 safe_deallocate_a(ilambdaminusdm)
505 integer,
intent(in) :: n
506 complex(real64),
intent(in) :: mat(:, :)
507 complex(real64),
intent(out) :: imat(:, :)
510 complex(real64),
allocatable :: u(:, :), vt(:, :), sigma(:, :)
511 real(real64),
allocatable :: sg_values(:)
515 safe_allocate(u(1:n, 1:n))
516 safe_allocate(vt(1:n, 1:n))
517 safe_allocate(sigma(1:n, 1:n))
518 safe_allocate(sg_values(1:n))
525 if (abs(sg_values(i)) <= 1.0e-12_real64 * maxval(abs(sg_values)))
then
528 sigma(i, i) = m_z1 / sg_values(i)
532 vt = conjg(transpose(vt))
533 u = conjg(transpose(u))
534 imat = matmul(vt, matmul(sigma, u))
537 vt = matmul(mat, matmul(imat, mat)) - mat
538 if (maxval(abs(vt)) > 1.0e-10_real64 * maxval(abs(mat)))
then
539 write(*, *) maxval(abs(vt))
542 write(*, *) 1.0e-10_real64 * maxval(abs(mat))
543 write(*, *) maxval(abs(vt)) > 1.0e-10_real64 * maxval(abs(mat))
545 write(message(1),
'(a)')
'Pseudoinverse failed.'
546 call messages_fatal(1)
550 safe_deallocate_a(vt)
551 safe_deallocate_a(sigma)
552 safe_deallocate_a(sg_values)
564 integer,
intent(in) :: n
565 complex(real64),
intent(in) :: mat(:, :)
567 integer :: alpha, beta, gamma, delta
568 complex(real64) :: deltah, zuder_direct, zder_direct, zuder_directplus, zuder_directminus
569 complex(real64),
allocatable :: zeigenvec(:, :), dm(:, :), zeigref_(:, :), zeigenval(:), mmatrix(:, :, :)
570 complex(real64),
allocatable :: zeigplus(:), zeigminus(:), zeig0(:), zeigplusminus(:), zeigminusplus(:)
574 safe_allocate(zeigenvec(1:n, 1:n))
575 safe_allocate(dm(1:n, 1:n))
576 safe_allocate(zeigref_(1:n, 1:n))
577 safe_allocate(zeigenval(1:n))
578 safe_allocate(mmatrix(1:n, 1:n, 1:n))
579 safe_allocate(zeigplus(1:n))
580 safe_allocate(zeigminus(1:n))
581 safe_allocate(zeig0(1:n))
582 safe_allocate(zeigplusminus(1:n))
583 safe_allocate(zeigminusplus(1:n))
591 deltah = (0.000001_real64, 0.000_real64)
595 zder_direct =
lalg_zdni(zeigref_(:, 2), alpha, beta)
596 zuder_direct =
lalg_zduialpha(zeigref_(:, 2), mmatrix(:, :, 2), 2, alpha, beta)
599 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
601 zeigenvec(:, 1) = zeigenvec(:, 1)/sum(conjg(zeigref_(1:2, 1))*zeigenvec(1:2, 1))
602 zeigenvec(:, 2) = zeigenvec(:, 2)/sum(conjg(zeigref_(1:2, 2))*zeigenvec(1:2, 2))
603 zuder_directplus = zeigenvec(2, 2)
606 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
608 zeigenvec(:, 1) = zeigenvec(:, 1)/sum(conjg(zeigref_(1:2, 1))*zeigenvec(1:2, 1))
609 zeigenvec(:, 2) = zeigenvec(:, 2)/sum(conjg(zeigref_(1:2, 2))*zeigenvec(1:2, 2))
610 zuder_directminus = zeigenvec(2, 2)
612 write(*,
'(2i1,4f24.12)') alpha, beta, zder_direct, (zeigplus(2) - zeigminus(2))/(m_two * deltah)
613 write(*,
'(2i1,4f24.12)') alpha, beta, &
614 zuder_direct, (zuder_directplus - zuder_directminus) / (m_two * deltah)
618 if (alpha == gamma .and. beta == delta)
then
619 zder_direct =
lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
625 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
629 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
632 write(*,
'(4i1,4f24.12)') alpha, beta, gamma, delta, &
633 zder_direct, (zeigplus(1) + zeigminus(1) - m_two*zeig0(1))/(deltah**2)
635 zder_direct =
lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
638 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
639 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
643 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
644 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
648 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
649 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
653 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
654 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
657 write(*,
'(4i1,4f24.12)') alpha, beta, gamma, delta, &
658 zder_direct, (zeigplus(1) + zeigminus(1) - zeigplusminus(1) - zeigminusplus(1))/(m_four*deltah**2)
668 safe_deallocate_a(zeigenval)
669 safe_deallocate_a(mmatrix)
670 safe_deallocate_a(zeigenvec)
671 safe_deallocate_a(dm)
672 safe_deallocate_a(zeigref_)
673 safe_deallocate_a(zeigplus)
674 safe_deallocate_a(zeigminus)
675 safe_deallocate_a(zeig0)
676 safe_deallocate_a(zeigplusminus)
677 safe_deallocate_a(zeigminusplus)
681 complex(real64) function lalg_zdni(eigenvec, alpha, beta)
682 integer,
intent(in) :: alpha, beta
683 complex(real64) :: eigenvec(2)
684 lalg_zdni = conjg(eigenvec(alpha)) * eigenvec(beta)
687 complex(real64) function lalg_zduialpha(eigenvec, mmatrix, alpha, gamma, delta)
688 integer,
intent(in) :: alpha, gamma, delta
689 complex(real64) :: eigenvec(2), mmatrix(2, 2)
693 complex(real64) function lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
694 integer,
intent(in) :: alpha, beta, gamma, delta
695 complex(real64) :: eigenvec(2), mmatrix(2, 2)
696 lalg_zd2ni = conjg(mmatrix(alpha, delta) * eigenvec(gamma)) * eigenvec(beta) + &
697 conjg(eigenvec(alpha)) * mmatrix(beta, gamma) * eigenvec(delta)
706 integer,
intent(in) :: m, n
707 real(real64),
intent(in) :: sg_values(:)
708 tol = m_epsilon * m * n * maxval(sg_values)
713#include "complex.F90"
714#include "lalg_adv_lapack_inc.F90"
718#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.
subroutine zmatrix_norm2(m, n, a, norm_l2, preserve_mat)
Norm of a 2D matrix.
complex(real64) function, public lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
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 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 NxN 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 zinverse(n, a, method, det, threshold, uplo)
An interface to different method to invert a matrix.
subroutine, public lalg_zpseudoinverse(n, mat, imat)
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 deigensolve(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A.
subroutine, public lalg_zeigenderivatives(n, mat, zeigenvec, zeigenval, zmat)
subroutine dmatrix_norm2(m, n, a, norm_l2, preserve_mat)
Norm of a 2D matrix.
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,...
subroutine, public lalg_check_zeigenderivatives(n, mat)
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)
subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
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, public zlalg_matrix_function(n, factor, a, fun_a, fun, hermitian)
This routine calculates a function of a matrix by using an eigenvalue decomposition.
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 inverse of a real NxN matrix a(:,:) using the SVD decomposition
subroutine, public zlalg_phi(nn, pp, aa, ex, hermitian)
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 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.