26 use,
intrinsic :: iso_fortran_env
148 real(real64) function sfmin()
150 real(real64) function dlamch(cmach)
153 character(1),
intent(in) :: cmach
160 subroutine lalg_dgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
161 character(1),
intent(in) :: jobvl, jobvr
162 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
163 real(real64),
intent(inout) :: a(:, :)
164 complex(real64),
intent(out) :: w(:)
165 real(real64),
intent(out) :: vl(:, :), vr(:, :)
166 real(real64),
intent(out) :: work(:)
167 real(real64),
intent(out) :: rwork(:)
168 integer,
intent(out) :: info
170 real(real64),
allocatable :: wr(:), wi(:)
174 safe_allocate(wr(1:n))
175 safe_allocate(wi(1:n))
177 call dgeev(jobvl, jobvr, n, a(1, 1), lda, wr(1), wi(1), vl(1, 1), ldvl, vr(1, 1), ldvr, work(1), lwork, info)
179 w(1:n) = cmplx(wr(1:n), wi(1:n), real64)
181 safe_deallocate_a(wr)
182 safe_deallocate_a(wi)
187 subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
188 character(1),
intent(in) :: jobvl, jobvr
189 integer,
intent(in) :: n, lda, ldvl, ldvr, lwork
190 complex(real64),
intent(inout) :: a(:, :)
191 complex(real64),
intent(out) :: w(:)
192 complex(real64),
intent(out) :: vl(:, :), vr(:, :)
193 real(real64),
intent(out) :: rwork(:)
194 complex(real64),
intent(out) :: work(:)
195 integer,
intent(out) :: info
199 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 integer,
intent(in) :: nn
223 complex(real64),
intent(in) :: pp
224 complex(real64),
intent(in) :: aa(:, :)
225 complex(real64),
intent(inout) :: ex(:, :)
226 logical,
intent(in) :: hermitian
228 complex(real64),
allocatable :: evectors(:, :), zevalues(:)
229 real(real64),
allocatable :: evalues(:)
235 safe_allocate(evectors(1:nn, 1:nn))
238 safe_allocate(evalues(1:nn))
239 safe_allocate(zevalues(1:nn))
241 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
245 zevalues(1:nn) =
exp(pp*evalues(1:nn))
248 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
251 ex(:, :) = matmul(evectors(:, :), ex(:, :))
253 safe_deallocate_a(evalues)
254 safe_deallocate_a(zevalues)
256 safe_allocate(zevalues(1:nn))
258 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
262 zevalues(1:nn) =
exp(pp*zevalues(1:nn))
264 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
269 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
272 ex(:, :) = matmul(ex(:, :), evectors(:, :))
274 safe_deallocate_a(zevalues)
277 safe_deallocate_a(evectors)
299 subroutine zlalg_phi(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(:, :) = aa(:, :)
324 zevalues(ii) = (
exp(pp*evalues(ii)) - m_z1) / (pp*evalues(ii))
328 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
331 ex(:, :) = matmul(evectors(:, :), ex(:, :))
333 safe_deallocate_a(evalues)
334 safe_deallocate_a(zevalues)
336 safe_allocate(zevalues(1:nn))
338 evectors(:, :) = aa(:, :)
343 zevalues(ii) = (
exp(pp*zevalues(ii)) - m_z1) / (pp*zevalues(ii))
346 ex(:, :) = evectors(:, :)
351 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
354 ex(:, :) = matmul(ex(:, :), evectors(:, :))
356 safe_deallocate_a(zevalues)
375 integer,
intent(in) :: n
376 complex(real64),
contiguous,
intent(in) :: mat(:, :)
377 complex(real64),
contiguous,
intent(out) :: zeigenvec(:, :)
378 complex(real64),
contiguous,
intent(out) :: zeigenval(:)
379 complex(real64),
contiguous,
intent(out) :: zmat(:, :, :)
381 integer :: i, alpha, beta
382 complex(real64),
allocatable :: knaught(:, :)
383 complex(real64),
allocatable :: lambdaminusdm(:, :)
384 complex(real64),
allocatable :: ilambdaminusdm(:, :)
385 complex(real64),
allocatable :: unit(:, :)
389 safe_allocate(unit(1:n, 1:n))
390 safe_allocate(knaught(1:n, 1:n))
391 safe_allocate(lambdaminusdm(1:n, 1:n))
392 safe_allocate(ilambdaminusdm(1:n, 1:n))
405 knaught(alpha, beta) = - zeigenvec(alpha, i)*conjg(zeigenvec(beta, i))
408 knaught = knaught + unit
409 lambdaminusdm = zeigenval(i)*unit - mat
412 zmat(:, :, i) = matmul(ilambdaminusdm, knaught)
415 safe_deallocate_a(unit)
416 safe_deallocate_a(knaught)
417 safe_deallocate_a(lambdaminusdm)
418 safe_deallocate_a(ilambdaminusdm)
427 integer,
intent(in) :: n
428 complex(real64),
contiguous,
intent(in) :: mat(:, :)
429 complex(real64),
contiguous,
intent(out) :: imat(:, :)
432 complex(real64),
allocatable :: u(:, :), vt(:, :), sigma(:, :)
433 real(real64),
allocatable :: sg_values(:)
437 safe_allocate(u(1:n, 1:n))
438 safe_allocate(vt(1:n, 1:n))
439 safe_allocate(sigma(1:n, 1:n))
440 safe_allocate(sg_values(1:n))
447 if (abs(sg_values(i)) <= 1.0e-12_real64 * maxval(abs(sg_values)))
then
450 sigma(i, i) = m_z1 / sg_values(i)
454 vt(:, :) = conjg(transpose(vt(:, :)))
455 u(:, :) = conjg(transpose(u(:, :)))
456 imat = matmul(vt, matmul(sigma, u))
459 vt(:, :) = matmul(mat(:, :), matmul(imat(:, :), mat(:, :))) - mat(:, :)
460 if (maxval(abs(vt)) > 1.0e-10_real64 * maxval(abs(mat)))
then
461 write(*, *) maxval(abs(vt))
464 write(*, *) 1.0e-10_real64 * maxval(abs(mat))
465 write(*, *) maxval(abs(vt)) > 1.0e-10_real64 * maxval(abs(mat))
467 write(message(1),
'(a)')
'Pseudoinverse failed.'
468 call messages_fatal(1)
472 safe_deallocate_a(vt)
473 safe_deallocate_a(sigma)
474 safe_deallocate_a(sg_values)
486 integer,
intent(in) :: n
487 complex(real64),
intent(in) :: mat(:, :)
489 integer :: alpha, beta, gamma, delta
490 complex(real64) :: deltah, zuder_direct, zder_direct, zuder_directplus, zuder_directminus
491 complex(real64),
allocatable :: zeigenvec(:, :), dm(:, :), zeigref_(:, :), zeigenval(:), mmatrix(:, :, :)
492 complex(real64),
allocatable :: zeigplus(:), zeigminus(:), zeig0(:), zeigplusminus(:), zeigminusplus(:)
496 safe_allocate(zeigenvec(1:n, 1:n))
497 safe_allocate(dm(1:n, 1:n))
498 safe_allocate(zeigref_(1:n, 1:n))
499 safe_allocate(zeigenval(1:n))
500 safe_allocate(mmatrix(1:n, 1:n, 1:n))
501 safe_allocate(zeigplus(1:n))
502 safe_allocate(zeigminus(1:n))
503 safe_allocate(zeig0(1:n))
504 safe_allocate(zeigplusminus(1:n))
505 safe_allocate(zeigminusplus(1:n))
513 deltah = (0.000001_real64, 0.000_real64)
517 zder_direct =
lalg_zdni(zeigref_(:, 2), alpha, beta)
518 zuder_direct =
lalg_zduialpha(zeigref_(:, 2), mmatrix(:, :, 2), 2, alpha, beta)
520 zeigenvec(:, :) = dm(:, :)
521 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
523 zeigenvec(:, 1) = zeigenvec(:, 1)/sum(conjg(zeigref_(1:2, 1))*zeigenvec(1:2, 1))
524 zeigenvec(:, 2) = zeigenvec(:, 2)/sum(conjg(zeigref_(1:2, 2))*zeigenvec(1:2, 2))
525 zuder_directplus = zeigenvec(2, 2)
527 zeigenvec(:, :) = dm(:, :)
528 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
530 zeigenvec(:, 1) = zeigenvec(:, 1)/sum(conjg(zeigref_(1:2, 1))*zeigenvec(1:2, 1))
531 zeigenvec(:, 2) = zeigenvec(:, 2)/sum(conjg(zeigref_(1:2, 2))*zeigenvec(1:2, 2))
532 zuder_directminus = zeigenvec(2, 2)
534 write(*,
'(2i1,4f24.12)') alpha, beta, zder_direct, (zeigplus(2) - zeigminus(2))/(m_two * deltah)
535 write(*,
'(2i1,4f24.12)') alpha, beta, &
536 zuder_direct, (zuder_directplus - zuder_directminus) / (m_two * deltah)
540 if (alpha == gamma .and. beta == delta)
then
541 zder_direct =
lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
543 zeigenvec(:, :) = dm(:, :)
546 zeigenvec(:, :) = dm(:, :)
547 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
550 zeigenvec(:, :) = dm(:, :)
551 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
554 write(*,
'(4i1,4f24.12)') alpha, beta, gamma, delta, &
555 zder_direct, (zeigplus(1) + zeigminus(1) - m_two*zeig0(1))/(deltah**2)
557 zder_direct =
lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
559 zeigenvec(:, :) = dm(:, :)
560 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
561 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
564 zeigenvec(:, :) = dm(:, :)
565 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
566 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
569 zeigenvec(:, :) = dm(:, :)
570 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
571 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
574 zeigenvec(:, :) = dm(:, :)
575 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
576 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
579 write(*,
'(4i1,4f24.12)') alpha, beta, gamma, delta, &
580 zder_direct, (zeigplus(1) + zeigminus(1) - zeigplusminus(1) - zeigminusplus(1))/(m_four*deltah**2)
590 safe_deallocate_a(zeigenval)
591 safe_deallocate_a(mmatrix)
592 safe_deallocate_a(zeigenvec)
593 safe_deallocate_a(dm)
594 safe_deallocate_a(zeigref_)
595 safe_deallocate_a(zeigplus)
596 safe_deallocate_a(zeigminus)
597 safe_deallocate_a(zeig0)
598 safe_deallocate_a(zeigplusminus)
599 safe_deallocate_a(zeigminusplus)
603 complex(real64) function lalg_zdni(eigenvec, alpha, beta)
604 integer,
intent(in) :: alpha, beta
605 complex(real64),
intent(in) :: eigenvec(2)
606 lalg_zdni = conjg(eigenvec(alpha)) * eigenvec(beta)
609 complex(real64) function lalg_zduialpha(eigenvec, mmatrix, alpha, gamma, delta)
610 integer,
intent(in) :: alpha, gamma, delta
611 complex(real64),
intent(in) :: eigenvec(2), mmatrix(2, 2)
615 complex(real64) function lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
616 integer,
intent(in) :: alpha, beta, gamma, delta
617 complex(real64),
intent(in) :: eigenvec(2), mmatrix(2, 2)
618 lalg_zd2ni = conjg(mmatrix(alpha, delta) * eigenvec(gamma)) * eigenvec(beta) + &
619 conjg(eigenvec(alpha)) * mmatrix(beta, gamma) * eigenvec(delta)
628 integer,
intent(in) :: m, n
629 real(real64),
intent(in) :: sg_values(:)
630 tol = m_epsilon * m * n * maxval(sg_values)
640 integer,
intent(in) :: n
641 real(real64),
intent(in) :: a(1:n, 1:n)
642 real(real64) :: p(1:n, 1:n)
644 real(real64) :: ata(1:n, 1:n)
648 ata = matmul(transpose(a), a)
655#include "complex.F90"
656#include "lalg_adv_lapack_inc.F90"
660#include "lalg_adv_lapack_inc.F90"
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
lalg_gemm with both the (Hermitian) transpose of A and B.
scales a vector by a constant
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)
pure real(real64) function, public pseudoinverse_default_tolerance(m, n, sg_values)
Computes the default Moore-Penrose pseudoinverse tolerance for zeroing.
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 lalg_zpseudoinverse(n, mat, imat)
Computes the Moore-Penrose pseudoinverse of a complex matrix.
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 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 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)
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 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.