Octopus
lalg_adv.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module lalg_adv_oct_m
23 use blacs_oct_m
24 use debug_oct_m
25 use global_oct_m
26 use, intrinsic :: iso_fortran_env
27 use lapack_oct_m
29 use mpi_oct_m
32 use sort_oct_m
33 use utils_oct_m
34#ifdef HAVE_ELPA
35 use elpa
36#endif
37
38 implicit none
39
40 private
41 public :: &
54 zlalg_exp, &
55 zlalg_phi, &
60 lalg_zdni, &
62 lalg_zd2ni, &
65
66
67
68 interface lalg_cholesky
69 module procedure dcholesky, zcholesky
70 end interface lalg_cholesky
71
72 interface lalg_geneigensolve
73 module procedure dgeneigensolve, zgeneigensolve
74 end interface lalg_geneigensolve
75
76 interface lalg_eigensolve_nonh
77 module procedure zeigensolve_nonh, deigensolve_nonh
78 end interface lalg_eigensolve_nonh
79
80 interface lalg_eigensolve
81 module procedure deigensolve, zeigensolve
82 end interface lalg_eigensolve
83
86 end interface lalg_eigensolve_parallel
87
90
91 interface lalg_determinant
92 module procedure ddeterminant, zdeterminant
93 end interface lalg_determinant
94
95 interface lalg_inverse
96 module procedure dinverse, zinverse
97 end interface lalg_inverse
98
99 interface lalg_linsyssolve
100 module procedure dlinsyssolve, zlinsyssolve
101 end interface lalg_linsyssolve
102
105 end interface lalg_singular_value_decomp
106
107 interface lalg_svd_inverse
108 module procedure dsvd_inverse, zsvd_inverse
109 end interface lalg_svd_inverse
110
111
115
116 interface lalg_lowest_eigensolve
117 module procedure dlowest_eigensolve, zlowest_eigensolve
118 end interface lalg_lowest_eigensolve
119
120 interface lapack_geev
121 module procedure lalg_dgeev, lalg_zgeev
122 end interface lapack_geev
123
124 interface lalg_least_squares
125 module procedure dleast_squares_vec, zleast_squares_vec
126 end interface lalg_least_squares
127
128 interface lalg_matrix_norm2
129 module procedure dmatrix_norm2, zmatrix_norm2
130 end interface lalg_matrix_norm2
131
132contains
133
134 ! ---------------------------------------------------------
136 real(real64) function sfmin()
137 interface
138 real(real64) function dlamch(cmach)
139 import real64
140 implicit none
141 character(1), intent(in) :: cmach
142 end function dlamch
143 end interface
144
145 sfmin = dlamch('S')
146 end function sfmin
147
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
157
158 real(real64), allocatable :: wr(:), wi(:)
159
160 push_sub(lalg_dgeev)
162 safe_allocate(wr(1:n))
163 safe_allocate(wi(1:n))
164
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)
166
167 w(1:n) = cmplx(wr(1:n), wi(1:n), real64)
168
169 safe_deallocate_a(wr)
170 safe_deallocate_a(wi)
171
172 pop_sub(lalg_dgeev)
173 end subroutine lalg_dgeev
174
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
185 push_sub(lalg_zgeev)
186
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)
189 pop_sub(lalg_zgeev)
190 end subroutine lalg_zgeev
191
208 subroutine zlalg_matrix_function(n, factor, a, fun_a, fun, hermitian)
209 integer, intent(in) :: n
210 complex(real64), intent(in) :: factor
211 complex(real64), intent(in) :: a(:, :)
212 complex(real64), intent(inout) :: fun_a(:, :)
213 interface
214 complex(real64) function fun(z)
215 import real64
216 complex(real64), intent(in) :: z
217 end
218 end interface
219 logical, intent(in) :: hermitian
220
221 complex(real64), allocatable :: evectors(:, :), zevalues(:)
222 real(real64), allocatable :: evalues(:)
223
224 integer :: i
225
226 push_sub(zlalg_matrix_function)
227
228 assert(size(a, 1) >= n)
230 safe_allocate(evectors(1:n, 1:n))
231
232 if (hermitian) then
233 safe_allocate(evalues(1:n))
234 safe_allocate(zevalues(1:n))
235
236 evectors(1:n, 1:n) = a(1:n, 1:n)
237
238 call lalg_eigensolve(n, evectors, evalues)
239
240 do i = 1, n
241 zevalues(i) = fun(factor*evalues(i))
242 end do
243
244 do i = 1, n
245 fun_a(1:n, i) = zevalues(1:n)*conjg(evectors(i, 1:n))
246 end do
247
248 fun_a(1:n, 1:n) = matmul(evectors(1:n, 1:n), fun_a(1:n, 1:n))
249
250 safe_deallocate_a(evalues)
251 safe_deallocate_a(zevalues)
252 else
253 safe_allocate(zevalues(1:n))
254
255 evectors(1:n, 1:n) = a(1:n, 1:n)
256
257 call lalg_eigensolve_nonh(n, evectors, zevalues)
258
259 do i = 1, n
260 zevalues(i) = fun(factor*zevalues(i))
261 end do
262
263 fun_a(1:n, 1:n) = evectors(1:n, 1:n)
264
265 call lalg_inverse(n, evectors, 'dir')
266
267 do i = 1, n
268 evectors(1:n, i) = zevalues(1:n)*evectors(1:n, i)
269 end do
270
271 fun_a(1:n, 1:n) = matmul(fun_a(1:n, 1:n), evectors(1:n, 1:n))
272
273 safe_deallocate_a(zevalues)
274 end if
275
276 safe_deallocate_a(evectors)
277
278 pop_sub(zlalg_matrix_function)
279 end subroutine zlalg_matrix_function
280
281
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
305
306 complex(real64), allocatable :: evectors(:, :), zevalues(:)
307 real(real64), allocatable :: evalues(:)
308
309 integer :: ii
310
311 push_sub(zlalg_exp)
312
313 safe_allocate(evectors(1:nn, 1:nn))
314
315 if (hermitian) then
316 safe_allocate(evalues(1:nn))
317 safe_allocate(zevalues(1:nn))
318
319 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
320
321 call lalg_eigensolve(nn, evectors, evalues)
322
323 zevalues(1:nn) = exp(pp*evalues(1:nn))
324
325 do ii = 1, nn
326 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
327 end do
328
329 ex(1:nn, 1:nn) = matmul(evectors(1:nn, 1:nn), ex(1:nn, 1:nn))
330
331 safe_deallocate_a(evalues)
332 safe_deallocate_a(zevalues)
333 else
334 safe_allocate(zevalues(1:nn))
335
336 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
337
338 call lalg_eigensolve_nonh(nn, evectors, zevalues)
339
340 zevalues(1:nn) = exp(pp*zevalues(1:nn))
341
342 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
343
344 call lalg_inverse(nn, evectors, 'dir')
345
346 do ii = 1, nn
347 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
348 end do
349
350 ex(1:nn, 1:nn) = matmul(ex(1:nn, 1:nn), evectors(1:nn, 1:nn))
351
352 safe_deallocate_a(zevalues)
353 end if
354
355 safe_deallocate_a(evectors)
356
357 pop_sub(zlalg_exp)
358 end subroutine zlalg_exp
359
360
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
383
384 complex(real64), allocatable :: evectors(:, :), zevalues(:)
385 real(real64), allocatable :: evalues(:)
386
387 integer :: ii
388
389 push_sub(zlalg_phi)
390
391 safe_allocate(evectors(1:nn, 1:nn))
393 if (hermitian) then
394 safe_allocate(evalues(1:nn))
395 safe_allocate(zevalues(1:nn))
396
397 evectors = aa
398
399 call lalg_eigensolve(nn, evectors, evalues)
400
401 do ii = 1, nn
402 zevalues(ii) = (exp(pp*evalues(ii)) - m_z1) / (pp*evalues(ii))
403 end do
404
405 do ii = 1, nn
406 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
407 end do
408
409 ex(1:nn, 1:nn) = matmul(evectors(1:nn, 1:nn), ex(1:nn, 1:nn))
410
411 safe_deallocate_a(evalues)
412 safe_deallocate_a(zevalues)
413 else
414 safe_allocate(zevalues(1:nn))
415
416 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
417
418 call lalg_eigensolve_nonh(nn, evectors, zevalues)
419
420 do ii = 1, nn
421 zevalues(ii) = (exp(pp*zevalues(ii)) - m_z1) / (pp*zevalues(ii))
422 end do
423
424 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
425
426 call lalg_inverse(nn, evectors, 'dir')
427
428 do ii = 1, nn
429 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
430 end do
431
432 ex(1:nn, 1:nn) = matmul(ex(1:nn, 1:nn), evectors(1:nn, 1:nn))
433
434 safe_deallocate_a(zevalues)
435 end if
436
437 pop_sub(zlalg_phi)
438 end subroutine zlalg_phi
439
440
452 subroutine lalg_zeigenderivatives(n, mat, zeigenvec, zeigenval, zmat)
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(:, :, :)
458
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(:, :)
464
465 push_sub(lalg_zeigenderivatives)
466
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))
471
472 zeigenvec = mat
473 call lalg_eigensolve_nonh(n, zeigenvec, zeigenval)
474
475 unit = m_z0
476 do i = 1, n
477 unit(i, i) = m_z1
478 end do
479
480 do i = 1, n
481 do alpha = 1, n
482 do beta = 1, n
483 knaught(alpha, beta) = - zeigenvec(alpha, i)*conjg(zeigenvec(beta, i))
484 end do
485 end do
486 knaught = knaught + unit
487 lambdaminusdm = zeigenval(i)*unit - mat
488 call lalg_zpseudoinverse(n, lambdaminusdm, ilambdaminusdm)
489 zmat(:, :, i) = matmul(ilambdaminusdm, knaught)
490 end do
491
492 safe_deallocate_a(unit)
493 safe_deallocate_a(knaught)
494 safe_deallocate_a(lambdaminusdm)
495 safe_deallocate_a(ilambdaminusdm)
497 end subroutine lalg_zeigenderivatives
498
499
504 subroutine lalg_zpseudoinverse(n, mat, imat)
505 integer, intent(in) :: n
506 complex(real64), intent(in) :: mat(:, :)
507 complex(real64), intent(out) :: imat(:, :)
508
509 integer :: i
510 complex(real64), allocatable :: u(:, :), vt(:, :), sigma(:, :)
511 real(real64), allocatable :: sg_values(:)
512
513 push_sub(lalg_zpseudoinverse)
514
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))
519
520 imat = mat
521 call lalg_singular_value_decomp(n, n, imat, u, vt, sg_values)
522
523 sigma = m_z0
524 do i = 1, n
525 if (abs(sg_values(i)) <= 1.0e-12_real64 * maxval(abs(sg_values))) then
526 sigma(i, i) = m_z0
527 else
528 sigma(i, i) = m_z1 / sg_values(i)
529 end if
530 end do
531
532 vt = conjg(transpose(vt))
533 u = conjg(transpose(u))
534 imat = matmul(vt, matmul(sigma, u))
535
536 ! Check if we truly have a pseudoinverse
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))
540 write(*, *) vt
541 write(*, *)
542 write(*, *) 1.0e-10_real64 * maxval(abs(mat))
543 write(*, *) maxval(abs(vt)) > 1.0e-10_real64 * maxval(abs(mat))
544 write(*, *) mat
545 write(message(1), '(a)') 'Pseudoinverse failed.'
546 call messages_fatal(1)
547 end if
548
549 safe_deallocate_a(u)
550 safe_deallocate_a(vt)
551 safe_deallocate_a(sigma)
552 safe_deallocate_a(sg_values)
553 pop_sub(lalg_zpseudoinverse)
554 end subroutine lalg_zpseudoinverse
555
556
563 subroutine lalg_check_zeigenderivatives(n, mat)
564 integer, intent(in) :: n
565 complex(real64), intent(in) :: mat(:, :)
566
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(:)
571
573
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))
584
585 assert(n == 2)
586
587 dm = mat
588 call lalg_zeigenderivatives(2, dm, zeigref_, zeigenval, mmatrix)
589
590
591 deltah = (0.000001_real64, 0.000_real64)
592 !deltah = M_z1 * maxval(abs(dm)) * 0.001_real64
593 do alpha = 1, n
594 do beta = 1, n
595 zder_direct = lalg_zdni(zeigref_(:, 2), alpha, beta)
596 zuder_direct = lalg_zduialpha(zeigref_(:, 2), mmatrix(:, :, 2), 2, alpha, beta)
598 zeigenvec = dm
599 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
600 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigplus)
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)
604
605 zeigenvec = dm
606 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
607 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigminus)
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)
611
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)
615
616 do gamma = 1, n
617 do delta = 1, n
618 if (alpha == gamma .and. beta == delta) then
619 zder_direct = lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
620
621 zeigenvec = dm
622 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeig0)
623
624 zeigenvec = dm
625 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
626 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigplus)
627
628 zeigenvec = dm
629 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
630 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigminus)
631
632 write(*, '(4i1,4f24.12)') alpha, beta, gamma, delta, &
633 zder_direct, (zeigplus(1) + zeigminus(1) - m_two*zeig0(1))/(deltah**2)
634 else
635 zder_direct = lalg_zd2ni(zeigref_(:, 1), mmatrix(:, :, 1), alpha, beta, gamma, delta)
636
637 zeigenvec = dm
638 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
639 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
640 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigplus)
641
642 zeigenvec = dm
643 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
644 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
645 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigminus)
646
647 zeigenvec = dm
648 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) + deltah
649 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) - deltah
650 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigplusminus)
651
652 zeigenvec = dm
653 zeigenvec(alpha, beta) = zeigenvec(alpha, beta) - deltah
654 zeigenvec(gamma, delta) = zeigenvec(gamma, delta) + deltah
655 call lalg_eigensolve_nonh(2, zeigenvec(:, :), zeigminusplus)
657 write(*, '(4i1,4f24.12)') alpha, beta, gamma, delta, &
658 zder_direct, (zeigplus(1) + zeigminus(1) - zeigplusminus(1) - zeigminusplus(1))/(m_four*deltah**2)
659
660
661 end if
662 end do
663 end do
664
665 end do
666 end do
667
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)
679 end subroutine lalg_check_zeigenderivatives
680
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)
685 end function lalg_zdni
686
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)
690 lalg_zduialpha = mmatrix(alpha, gamma) * eigenvec(delta)
691 end function lalg_zduialpha
692
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)
698 end function lalg_zd2ni
699
700 ! ---------------------------------------------------------
705 pure real(real64) function pseudoinverse_default_tolerance(m, n, sg_values) result(tol)
706 integer, intent(in) :: m, n
707 real(real64), intent(in) :: sg_values(:)
708 tol = m_epsilon * m * n * maxval(sg_values)
710
711
712#include "undef.F90"
713#include "complex.F90"
714#include "lalg_adv_lapack_inc.F90"
715
716#include "undef.F90"
717#include "real.F90"
718#include "lalg_adv_lapack_inc.F90"
719
720end module lalg_adv_oct_m
721
722!! Local Variables:
723!! mode: f90
724!! coding: utf-8
725!! End:
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
Definition: lalg_adv.F90:184
double exp(double __x) __attribute__((__nothrow__
This module contains interfaces for BLACS routines Interfaces are from http:
Definition: blacs.F90:27
This module provides the BLACS processor grid.
subroutine zmatrix_norm2(m, n, a, norm_l2, preserve_mat)
Norm of a 2D matrix.
Definition: lalg_adv.F90:1547
complex(real64) function, public lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
Definition: lalg_adv.F90:787
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...
Definition: lalg_adv.F90:1150
complex(real64) function zdeterminant(n, a, preserve_mat)
Invert a real symmetric or complex Hermitian square matrix a.
Definition: lalg_adv.F90:1432
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...
Definition: lalg_adv.F90:1031
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...
Definition: lalg_adv.F90:1632
subroutine, public zlalg_exp(nn, pp, aa, ex, hermitian)
Definition: lalg_adv.F90:393
subroutine deigensolve_parallel(n, a, e, bof, err_code)
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenprobl...
Definition: lalg_adv.F90:3298
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...
Definition: lalg_adv.F90:2265
subroutine zeigensolve(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A.
Definition: lalg_adv.F90:1269
subroutine dsingular_value_decomp(m, n, a, u, vt, sg_values, preserve_mat)
computes the singular value decomposition of a real NxN matrix a(:,:)
Definition: lalg_adv.F90:3046
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...
Definition: lalg_adv.F90:2367
subroutine zinverse(n, a, method, det, threshold, uplo)
An interface to different method to invert a matrix.
Definition: lalg_adv.F90:2104
subroutine, public lalg_zpseudoinverse(n, mat, imat)
Definition: lalg_adv.F90:598
real(real64) function sfmin()
Auxiliary function.
Definition: lalg_adv.F90:230
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 ...
Definition: lalg_adv.F90:2967
subroutine zsingular_value_decomp(m, n, a, u, vt, sg_values, preserve_mat)
computes the singular value decomposition of a complex MxN matrix a(:,:)
Definition: lalg_adv.F90:1711
subroutine dcholesky(n, a, bof, err_code)
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a,...
Definition: lalg_adv.F90:2212
subroutine deigensolve(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric or hermitian square matrix A.
Definition: lalg_adv.F90:2606
subroutine, public lalg_zeigenderivatives(n, mat, zeigenvec, zeigenval, zmat)
Definition: lalg_adv.F90:546
subroutine dmatrix_norm2(m, n, a, norm_l2, preserve_mat)
Norm of a 2D matrix.
Definition: lalg_adv.F90:2882
subroutine zlowest_eigensolve(k, n, a, e, v, preserve_mat)
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem,...
Definition: lalg_adv.F90:1341
complex(real64) function, public lalg_zduialpha(eigenvec, mmatrix, alpha, gamma, delta)
Definition: lalg_adv.F90:781
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...
Definition: lalg_adv.F90:929
subroutine dleast_squares_vec(nn, aa, bb, xx, preserve_mat)
Definition: lalg_adv.F90:3226
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...
Definition: lalg_adv.F90:2486
subroutine zcholesky(n, a, bof, err_code)
Compute the Cholesky decomposition of real symmetric or complex Hermitian positive definite matrix a,...
Definition: lalg_adv.F90:876
subroutine, public lalg_check_zeigenderivatives(n, mat)
Definition: lalg_adv.F90:657
complex(real64) function, public lalg_zdni(eigenvec, alpha, beta)
Definition: lalg_adv.F90:775
pure real(real64) function pseudoinverse_default_tolerance(m, n, sg_values)
Computes the default Moore-Penrose pseudoinverse tolerance for zeroing.
Definition: lalg_adv.F90:799
subroutine zleast_squares_vec(nn, aa, bb, xx, preserve_mat)
Definition: lalg_adv.F90:1896
subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
Definition: lalg_adv.F90:269
subroutine dinverse(n, a, method, det, threshold, uplo)
An interface to different method to invert a matrix.
Definition: lalg_adv.F90:3431
subroutine zeigensolve_parallel(n, a, e, bof, err_code)
Computes all the eigenvalues and the eigenvectors of a real symmetric or complex Hermitian eigenprobl...
Definition: lalg_adv.F90:1971
real(real64) function ddeterminant(n, a, preserve_mat)
Invert a real symmetric or complex Hermitian square matrix a.
Definition: lalg_adv.F90:2767
subroutine lalg_dgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
Definition: lalg_adv.F90:242
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.
Definition: lalg_adv.F90:302
subroutine dlowest_eigensolve(k, n, a, e, v, preserve_mat)
Computes the k lowest eigenvalues and the eigenvectors of a standard symmetric-definite eigenproblem,...
Definition: lalg_adv.F90:2678
subroutine dsvd_inverse(m, n, a, threshold)
computes inverse of a real NxN matrix a(:,:) using the SVD decomposition
Definition: lalg_adv.F90:3124
subroutine, public zlalg_phi(nn, pp, aa, ex, hermitian)
Definition: lalg_adv.F90:471
subroutine zsvd_inverse(m, n, a, threshold)
computes inverse of a complex MxN matrix a(:,:) using the SVD decomposition
Definition: lalg_adv.F90:1793
This module contains interfaces for LAPACK routines.
Definition: lapack.F90:118
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
Definition: scalapack.F90:131
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:117
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118