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
22 use blas_oct_m
24 use blacs_oct_m
25 use debug_oct_m
26 use global_oct_m
27 use, intrinsic :: iso_fortran_env
28 use lapack_oct_m
30 use math_oct_m
32 use mpi_oct_m
35 use sort_oct_m
36 use utils_oct_m
37#ifdef HAVE_ELPA
38 use elpa
39#endif
40
41 implicit none
42
43 private
44 public :: &
59 zlalg_exp, &
60 zlalg_phi, &
63 lalg_zdni, &
65 lalg_zd2ni, &
69
70 interface lalg_cholesky
71 module procedure dcholesky, zcholesky
72 end interface lalg_cholesky
73
74 interface lalg_geneigensolve
75 module procedure dgeneigensolve, zgeneigensolve
76 end interface lalg_geneigensolve
77
78 interface lalg_eigensolve_nonh
79 module procedure zeigensolve_nonh, deigensolve_nonh
80 end interface lalg_eigensolve_nonh
81
82 interface lalg_eigensolve
83 module procedure deigensolve, zeigensolve
84 end interface lalg_eigensolve
85
88 end interface lalg_eigensolve_tridiagonal
89
92 end interface lalg_eigensolve_parallel
93
96
97 interface lalg_determinant
98 module procedure ddeterminant, zdeterminant
99 end interface lalg_determinant
100
101 interface lalg_inverse
102 module procedure dinverse, zinverse
103 end interface lalg_inverse
104
105 interface lalg_linsyssolve
106 module procedure dlinsyssolve, zlinsyssolve
107 end interface lalg_linsyssolve
108
111 end interface lalg_singular_value_decomp
112
113 interface lalg_pseudo_inverse
115 end interface lalg_pseudo_inverse
117 interface lalg_svd_inverse
118 module procedure dsvd_inverse, zsvd_inverse
119 end interface lalg_svd_inverse
120
123 end interface lalg_lowest_geneigensolve
124
125 interface lalg_lowest_eigensolve
126 module procedure dlowest_eigensolve, zlowest_eigensolve
127 end interface lalg_lowest_eigensolve
128
129 interface lapack_geev
130 module procedure lalg_dgeev, lalg_zgeev
131 end interface lapack_geev
132
133 interface lalg_least_squares
134 module procedure dleast_squares_vec, zleast_squares_vec
135 end interface lalg_least_squares
136
137 interface lalg_matrix_function
139 end interface lalg_matrix_function
140
141 interface lalg_matrix_rank_svd
142 module procedure dmatrix_rank_svd, zmatrix_rank_svd
143 end interface lalg_matrix_rank_svd
144
145contains
146
147 ! ---------------------------------------------------------
149 real(real64) function sfmin()
150 interface
151 real(real64) function dlamch(cmach)
152 import real64
153 implicit none
154 character(1), intent(in) :: cmach
155 end function dlamch
156 end interface
157
158 sfmin = dlamch('S')
159 end function sfmin
160
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
170
171 real(real64), allocatable :: wr(:), wi(:)
172
173 push_sub(lalg_dgeev)
174
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)
179
180 w(1:n) = cmplx(wr(1:n), wi(1:n), real64)
182 safe_deallocate_a(wr)
183 safe_deallocate_a(wi)
184
185 pop_sub(lalg_dgeev)
186 end subroutine lalg_dgeev
187
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
197
198 push_sub(lalg_zgeev)
199
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)
201
202 pop_sub(lalg_zgeev)
203 end subroutine lalg_zgeev
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(:)
231
232 integer :: ii
233
234 push_sub(zlalg_exp)
235
236 safe_allocate(evectors(1:nn, 1:nn))
237
238 if (hermitian) then
239 safe_allocate(evalues(1:nn))
240 safe_allocate(zevalues(1:nn))
241
242 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
243
244 call lalg_eigensolve(nn, evectors, evalues)
245
246 zevalues(1:nn) = exp(pp*evalues(1:nn))
247
248 do ii = 1, nn
249 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
250 end do
251
252 ex(:, :) = matmul(evectors(:, :), ex(:, :))
253
254 safe_deallocate_a(evalues)
255 safe_deallocate_a(zevalues)
256 else
257 safe_allocate(zevalues(1:nn))
258
259 evectors(1:nn, 1:nn) = aa(1:nn, 1:nn)
260
261 call lalg_eigensolve_nonh(nn, evectors, zevalues)
262
263 zevalues(1:nn) = exp(pp*zevalues(1:nn))
264
265 ex(1:nn, 1:nn) = evectors(1:nn, 1:nn)
266
267 call lalg_inverse(nn, evectors, 'dir')
268
269 do ii = 1, nn
270 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
271 end do
272
273 ex(:, :) = matmul(ex(:, :), evectors(:, :))
274
275 safe_deallocate_a(zevalues)
276 end if
277
278 safe_deallocate_a(evectors)
279
280 pop_sub(zlalg_exp)
281 end subroutine zlalg_exp
282
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
306
307 complex(real64), allocatable :: evectors(:, :), zevalues(:)
308 real(real64), allocatable :: evalues(:)
309
310 integer :: ii
311
312 push_sub(zlalg_phi)
313
314 safe_allocate(evectors(1:nn, 1:nn))
315
316 if (hermitian) then
317 safe_allocate(evalues(1:nn))
318 safe_allocate(zevalues(1:nn))
319
320 evectors(:, :) = aa(:, :)
321
322 call lalg_eigensolve(nn, evectors, evalues)
323
324 do ii = 1, nn
325 zevalues(ii) = (exp(pp*evalues(ii)) - m_z1) / (pp*evalues(ii))
326 end do
327
328 do ii = 1, nn
329 ex(1:nn, ii) = zevalues(1:nn)*conjg(evectors(ii, 1:nn))
330 end do
331
332 ex(:, :) = matmul(evectors(:, :), ex(:, :))
333
334 safe_deallocate_a(evalues)
335 safe_deallocate_a(zevalues)
336 else
337 safe_allocate(zevalues(1:nn))
338
339 evectors(:, :) = aa(:, :)
340
341 call lalg_eigensolve_nonh(nn, evectors, zevalues)
342
343 do ii = 1, nn
344 zevalues(ii) = (exp(pp*zevalues(ii)) - m_z1) / (pp*zevalues(ii))
345 end do
346
347 ex(:, :) = evectors(:, :)
348
349 call lalg_inverse(nn, evectors, 'dir')
350
351 do ii = 1, nn
352 evectors(1:nn, ii) = zevalues(1:nn)*evectors(1:nn, ii)
353 end do
354
355 ex(:, :) = matmul(ex(:, :), evectors(:, :))
356
357 safe_deallocate_a(zevalues)
358 end if
359
360 pop_sub(zlalg_phi)
361 end subroutine zlalg_phi
362
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)
367 end function lalg_zdni
368
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)
372 lalg_zduialpha = mmatrix(alpha, gamma) * eigenvec(delta)
373 end function lalg_zduialpha
374
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)
380 end function lalg_zd2ni
381
382 ! ---------------------------------------------------------
387 pure real(real64) function pseudoinverse_default_tolerance(m, n, sg_values) result(tol)
388 integer, intent(in) :: m, n
389 real(real64), intent(in) :: sg_values(:)
390 tol = m_epsilon * m * n * maxval(sg_values)
392
393
399 function lalg_remove_rotation(n, A) result(P)
400 integer, intent(in) :: n
401 real(real64), intent(in) :: a(1:n, 1:n)
402 real(real64) :: p(1:n, 1:n)
403
404 real(real64) :: ata(1:n, 1:n)
405
406 push_sub(lalg_remove_rotation)
407
408 ata = matmul(transpose(a), a)
409 call lalg_matrix_function(n, m_one, ata, p, square_root, .true.)
410
411 pop_sub(lalg_remove_rotation)
412 end function lalg_remove_rotation
413
414#include "undef.F90"
415#include "complex.F90"
416#include "lalg_adv_lapack_inc.F90"
417
418#include "undef.F90"
419#include "real.F90"
420#include "lalg_adv_lapack_inc.F90"
421
422end module lalg_adv_oct_m
423
424!! Local Variables:
425!! mode: f90
426!! coding: utf-8
427!! End:
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
Definition: lalg_adv.F90:192
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.
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:120
complex(real64) function, public lalg_zd2ni(eigenvec, mmatrix, alpha, beta, gamma, delta)
Definition: lalg_adv.F90:471
integer function zmatrix_rank_svd(a, preserve_mat, tol)
Compute the rank of the matrix A using SVD.
Definition: lalg_adv.F90:1674
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:865
complex(real64) function zdeterminant(n, a, preserve_mat)
Invert a real symmetric or complex Hermitian square matrix a.
Definition: lalg_adv.F90:1224
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:746
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:1381
subroutine, public zlalg_exp(nn, pp, aa, ex, hermitian)
Definition: lalg_adv.F90:318
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:3336
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:2146
subroutine zlalg_pseudo_inverse(a, threshold)
Invert a matrix with the Moore-Penrose pseudo-inverse.
Definition: lalg_adv.F90:1608
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:984
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.
Definition: lalg_adv.F90:2964
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:2252
subroutine deigensolve_tridiagonal(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian c...
Definition: lalg_adv.F90:2561
subroutine zinverse(n, a, method, det, threshold, uplo)
An interface to different method to invert a matrix.
Definition: lalg_adv.F90:1861
real(real64) function sfmin()
Auxiliary function.
Definition: lalg_adv.F90:245
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:2886
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:1459
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:2091
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.
Definition: lalg_adv.F90:3407
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:2491
subroutine zeigensolve_tridiagonal(n, a, e, bof, err_code)
Computes all eigenvalues and eigenvectors of a real symmetric tridiagonal matrix. For the Hermitian c...
Definition: lalg_adv.F90:1054
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:1133
complex(real64) function, public lalg_zduialpha(eigenvec, mmatrix, alpha, gamma, delta)
Definition: lalg_adv.F90:465
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:635
subroutine dleast_squares_vec(nn, aa, bb, xx, preserve_mat)
Definition: lalg_adv.F90:3264
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:2371
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:580
complex(real64) function, public lalg_zdni(eigenvec, alpha, beta)
Definition: lalg_adv.F90:459
pure real(real64) function pseudoinverse_default_tolerance(m, n, sg_values)
Computes the default Moore-Penrose pseudoinverse tolerance for zeroing.
Definition: lalg_adv.F90:483
subroutine zleast_squares_vec(nn, aa, bb, xx, preserve_mat)
Definition: lalg_adv.F90:1763
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 ...
Definition: lalg_adv.F90:495
subroutine lalg_zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
Definition: lalg_adv.F90:284
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.
Definition: lalg_adv.F90:1909
integer function dmatrix_rank_svd(a, preserve_mat, tol)
Compute the rank of the matrix A using SVD.
Definition: lalg_adv.F90:3175
subroutine dinverse(n, a, method, det, threshold, uplo)
An interface to different method to invert a matrix.
Definition: lalg_adv.F90:3359
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:1838
real(real64) function ddeterminant(n, a, preserve_mat)
Invert a real symmetric or complex Hermitian square matrix a.
Definition: lalg_adv.F90:2729
subroutine lalg_dgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
Definition: lalg_adv.F90:257
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:2640
subroutine dsvd_inverse(m, n, a, threshold)
Computes the inverse of a real M x N matrix, a, using the SVD decomposition.
Definition: lalg_adv.F90:3043
subroutine, public zlalg_phi(nn, pp, aa, ex, hermitian)
Definition: lalg_adv.F90:396
subroutine dlalg_pseudo_inverse(a, threshold)
Invert a matrix with the Moore-Penrose pseudo-inverse.
Definition: lalg_adv.F90:3109
subroutine zsvd_inverse(m, n, a, threshold)
Computes inverse of a complex MxN matrix, a, using the SVD decomposition.
Definition: lalg_adv.F90:1542
This module contains interfaces for LAPACK routines.
Definition: lapack.F90:120
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module contains interfaces for ScaLAPACK routines Interfaces are from http:
Definition: scalapack.F90:133
This module is intended to contain "only mathematical" functions and procedures.
Definition: sort.F90:119
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
int true(void)