26 use,
intrinsic :: iso_fortran_env
143 real(real64) function sfmin()
145 real(real64) function dlamch(cmach)
148 character(1),
intent(in) :: cmach
155 subroutine lalg_dgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
156 character(1),
intent(in) :: jobvl, jobvr
157 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
158 real(real64),
intent(inout) :: a(:,:)
159 complex(real64),
intent(out) :: w(:)
160 real(real64),
intent(out) :: vl(:,:), vr(:,:)
161 real(real64),
intent(out) :: rwork(:)
162 real(real64),
intent(out) :: work(:)
163 integer,
intent(out) :: info
165 real(real64),
allocatable :: wr(:), wi(:)
169 safe_allocate(wr(1:n))
170 safe_allocate(wi(1:n))
172 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)
174 w(1:n) = cmplx(wr(1:n), wi(1:n), real64)
176 safe_deallocate_a(wr)
177 safe_deallocate_a(wi)
182 subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
183 character(1),
intent(in) :: jobvl, jobvr
184 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
185 complex(real64),
intent(inout) :: a(:,:)
186 complex(real64),
intent(out) :: w(:)
187 complex(real64),
intent(out) :: vl(:,:), vr(:,:)
188 real(real64),
intent(out) :: rwork(:)
189 complex(real64),
intent(out) :: work(:)
190 integer,
intent(out) :: info
194 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)
217 integer,
intent(in) :: nn
218 complex(real64),
intent(in) :: pp
219 complex(real64),
intent(in) :: aa(:, :)
220 complex(real64),
intent(inout) :: ex(:, :)
221 logical,
intent(in) :: hermitian
223 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
224 real(real64),
allocatable :: evalues(:)
230 safe_allocate(evectors(1:nn, 1:nn))
233 safe_allocate(evalues(1:nn))
234 safe_allocate(zevalues(1:nn))
236 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
240 zevalues(1:nn) =
exp(pp*evalues(1:nn))
243 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
246 ex(1:nn, 1:nn) = matmul(evectors(1:nn, 1:nn), ex(1:nn, 1:nn))
248 safe_deallocate_a(evalues)
249 safe_deallocate_a(zevalues)
251 safe_allocate(zevalues(1:nn))
253 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
257 zevalues(1:nn) =
exp(pp*zevalues(1:nn))
259 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
264 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
267 ex(1:nn, 1:nn) = matmul(ex(1:nn, 1:nn), evectors(1:nn, 1:nn))
269 safe_deallocate_a(zevalues)
272 safe_deallocate_a(evectors)
294 subroutine zlalg_phi(nn, pp, aa, ex, hermitian)
295 integer,
intent(in) :: nn
296 complex(real64),
intent(in) :: pp
297 complex(real64),
intent(in) :: aa(:, :)
298 complex(real64),
intent(inout) :: ex(:, :)
299 logical,
intent(in) :: hermitian
301 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
302 real(real64),
allocatable :: evalues(:)
308 safe_allocate(evectors(1:nn, 1:nn))
311 safe_allocate(evalues(1:nn))
312 safe_allocate(zevalues(1:nn))
319 zevalues(ii) = (
exp(pp*evalues(ii)) - m_z1) / (pp*evalues(ii))
323 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
326 ex(1:nn, 1:nn) = matmul(evectors(1:nn, 1:nn), ex(1:nn, 1:nn))
328 safe_deallocate_a(evalues)
329 safe_deallocate_a(zevalues)
331 safe_allocate(zevalues(1:nn))
333 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
338 zevalues(ii) = (
exp(pp*zevalues(ii)) - m_z1) / (pp*zevalues(ii))
341 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
346 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
349 ex(1:nn, 1:nn) = matmul(ex(1:nn, 1:nn), evectors(1:nn, 1:nn))
351 safe_deallocate_a(zevalues)
370 integer,
intent(in) :: n
371 complex(real64),
intent(in) :: mat(:, :)
372 complex(real64),
intent(out) :: zeigenvec(:, :)
373 complex(real64),
intent(out) :: zeigenval(:)
374 complex(real64),
intent(out) :: zmat(:, :, :)
376 integer :: i, alpha, beta
377 complex(real64),
allocatable :: knaught(:, :)
378 complex(real64),
allocatable :: lambdaminusdm(:, :)
379 complex(real64),
allocatable :: ilambdaminusdm(:, :)
380 complex(real64),
allocatable :: unit(:, :)
384 safe_allocate(unit(1:n, 1:n))
385 safe_allocate(knaught(1:n, 1:n))
386 safe_allocate(lambdaminusdm(1:n, 1:n))
387 safe_allocate(ilambdaminusdm(1:n, 1:n))
400 knaught(alpha, beta) = - zeigenvec(alpha, i)*conjg(zeigenvec(beta, i))
403 knaught = knaught + unit
404 lambdaminusdm = zeigenval(i)*unit - mat
406 zmat(:, :, i) = matmul(ilambdaminusdm, knaught)
409 safe_deallocate_a(unit)
410 safe_deallocate_a(knaught)
411 safe_deallocate_a(lambdaminusdm)
412 safe_deallocate_a(ilambdaminusdm)
422 integer,
intent(in) :: n
423 complex(real64),
intent(in) :: mat(:, :)
424 complex(real64),
intent(out) :: imat(:, :)
427 complex(real64),
allocatable :: u(:, :), vt(:, :), sigma(:, :)
428 real(real64),
allocatable :: sg_values(:)
432 safe_allocate(u(1:n, 1:n))
433 safe_allocate(vt(1:n, 1:n))
434 safe_allocate(sigma(1:n, 1:n))
435 safe_allocate(sg_values(1:n))
442 if (abs(sg_values(i)) <= 1.0e-12_real64 * maxval(abs(sg_values)))
then
445 sigma(i, i) = m_z1 / sg_values(i)
449 vt = conjg(transpose(vt))
450 u = conjg(transpose(u))
451 imat = matmul(vt, matmul(sigma, u))
454 vt = matmul(mat, matmul(imat, mat)) - mat
455 if (maxval(abs(vt)) > 1.0e-10_real64 * maxval(abs(mat)))
then
456 write(*, *) maxval(abs(vt))
459 write(*, *) 1.0e-10_real64 * maxval(abs(mat))
460 write(*, *) maxval(abs(vt)) > 1.0e-10_real64 * maxval(abs(mat))
462 write(message(1),
'(a)')
'Pseudoinverse failed.'
463 call messages_fatal(1)
467 safe_deallocate_a(vt)
468 safe_deallocate_a(sigma)
469 safe_deallocate_a(sg_values)
481 integer,
intent(in) :: n
482 complex(real64),
intent(in) :: mat(:, :)
484 integer :: alpha, beta, gamma, delta
485 complex(real64) :: deltah, zuder_direct, zder_direct, zuder_directplus, zuder_directminus
486 complex(real64),
allocatable :: zeigenvec(:, :), dm(:, :), zeigref_(:, :), zeigenval(:), mmatrix(:, :, :)
487 complex(real64),
allocatable :: zeigplus(:), zeigminus(:), zeig0(:), zeigplusminus(:), zeigminusplus(:)
491 safe_allocate(zeigenvec(1:n, 1:n))
492 safe_allocate(dm(1:n, 1:n))
493 safe_allocate(zeigref_(1:n, 1:n))
494 safe_allocate(zeigenval(1:n))
495 safe_allocate(mmatrix(1:n, 1:n, 1:n))
496 safe_allocate(zeigplus(1:n))
497 safe_allocate(zeigminus(1:n))
498 safe_allocate(zeig0(1:n))
499 safe_allocate(zeigplusminus(1:n))
500 safe_allocate(zeigminusplus(1:n))
508 deltah = (0.000001_real64, 0.000_real64)
512 zder_direct =
lalg_zdni(zeigref_(:, 2), alpha, beta)
513 zuder_direct =
lalg_zduialpha(zeigref_(:, 2), mmatrix(:, :, 2), 2, alpha, beta)
516 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
518 zeigenvec(:, 1) = zeigenvec(:, 1)/sum(conjg(zeigref_(1:2, 1))*zeigenvec(1:2, 1))
519 zeigenvec(:, 2) = zeigenvec(:, 2)/sum(conjg(zeigref_(1:2, 2))*zeigenvec(1:2, 2))
520 zuder_directplus = zeigenvec(2, 2)
523 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
525 zeigenvec(:, 1) = zeigenvec(:, 1)/sum(conjg(zeigref_(1:2, 1))*zeigenvec(1:2, 1))
526 zeigenvec(:, 2) = zeigenvec(:, 2)/sum(conjg(zeigref_(1:2, 2))*zeigenvec(1:2, 2))
527 zuder_directminus = zeigenvec(2, 2)
529 write(*,
'(2i1,4f24.12)') alpha, beta, zder_direct, (zeigplus(2) - zeigminus(2))/(m_two * deltah)
530 write(*,
'(2i1,4f24.12)') alpha, beta, &
531 zuder_direct, (zuder_directplus - zuder_directminus) / (m_two * deltah)
535 if (alpha == gamma .and. beta == delta)
then
536 zder_direct =
lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
542 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
546 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
549 write(*,
'(4i1,4f24.12)') alpha, beta, gamma, delta, &
550 zder_direct, (zeigplus(1) + zeigminus(1) - m_two*zeig0(1))/(deltah**2)
552 zder_direct =
lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
555 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
556 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
560 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
561 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
565 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
566 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
570 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
571 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
574 write(*,
'(4i1,4f24.12)') alpha, beta, gamma, delta, &
575 zder_direct, (zeigplus(1) + zeigminus(1) - zeigplusminus(1) - zeigminusplus(1))/(m_four*deltah**2)
585 safe_deallocate_a(zeigenval)
586 safe_deallocate_a(mmatrix)
587 safe_deallocate_a(zeigenvec)
588 safe_deallocate_a(dm)
589 safe_deallocate_a(zeigref_)
590 safe_deallocate_a(zeigplus)
591 safe_deallocate_a(zeigminus)
592 safe_deallocate_a(zeig0)
593 safe_deallocate_a(zeigplusminus)
594 safe_deallocate_a(zeigminusplus)
598 complex(real64) function lalg_zdni(eigenvec, alpha, beta)
599 integer,
intent(in) :: alpha, beta
600 complex(real64) :: eigenvec(2)
601 lalg_zdni = conjg(eigenvec(alpha)) * eigenvec(beta)
604 complex(real64) function lalg_zduialpha(eigenvec, mmatrix, alpha, gamma, delta)
605 integer,
intent(in) :: alpha, gamma, delta
606 complex(real64) :: eigenvec(2), mmatrix(2, 2)
610 complex(real64) function lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
611 integer,
intent(in) :: alpha, beta, gamma, delta
612 complex(real64) :: eigenvec(2), mmatrix(2, 2)
613 lalg_zd2ni = conjg(mmatrix(alpha, delta) * eigenvec(gamma)) * eigenvec(beta) + &
614 conjg(eigenvec(alpha)) * mmatrix(beta, gamma) * eigenvec(delta)
623 integer,
intent(in) :: m, n
624 real(real64),
intent(in) :: sg_values(:)
625 tol = m_epsilon * m * n * maxval(sg_values)
635 integer,
intent(in) :: n
636 real(real64),
intent(in) :: a(1:n, 1:n)
637 real(real64) :: p(1:n, 1:n)
639 real(real64) :: ata(1:n, 1:n)
643 ata = matmul(transpose(a), a)
650#include "complex.F90"
651#include "lalg_adv_lapack_inc.F90"
655#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)
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 dlalg_matrix_function(n, factor, a, fun_a, fun, hermitian)
This routine calculates a function of a matrix by using an eigenvalue decomposition.
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 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.