7submodule(
linalg) linalg_eigen
10 module subroutine eigen_symm(vecs, a, vals, work, olwork, err)
12 logical,
intent(in) :: vecs
13 real(real64),
intent(inout),
dimension(:,:) :: a
14 real(real64),
intent(out),
dimension(:) :: vals
15 real(real64),
intent(out),
pointer,
optional,
dimension(:) :: work
16 integer(int32),
intent(out),
optional :: olwork
17 class(errors),
intent(inout),
optional,
target :: err
21 integer(int32) :: n, istat, flag, lwork
22 real(real64),
pointer,
dimension(:) :: wptr
23 real(real64),
allocatable,
target,
dimension(:) :: wrk
24 real(real64),
dimension(1) :: temp
25 class(errors),
pointer :: errmgr
26 type(errors),
target :: deferr
27 character(len = :),
allocatable :: errmsg
36 if (
present(err))
then
44 if (
size(a, 2) /= n)
then
46 else if (
size(vals) /= n)
then
51 allocate(
character(len = 256) :: errmsg)
52 write(errmsg, 100)
"Input number ", flag, &
53 " is not sized correctly."
54 call errmgr%report_error(
"eigen_symm", trim(errmsg), &
60 call dsyev(jobz,
'L', n, a, n, vals, temp, -1, flag)
61 lwork = int(temp(1), int32)
62 if (
present(olwork))
then
68 if (
present(work))
then
69 if (
size(work) < lwork)
then
71 call errmgr%report_error(
"eigen_symm", &
72 "Incorrectly sized input array WORK, argument 5.", &
78 allocate(wrk(lwork), stat = istat)
81 call errmgr%report_error(
"eigen_symm", &
82 "Insufficient memory available.", &
83 la_out_of_memory_error)
90 call dsyev(jobz,
'L', n, a, n, vals, wptr, lwork, flag)
92 call errmgr%report_error(
"eigen_symm", &
93 "The algorithm failed to converge.", la_convergence_error)
101 module subroutine eigen_asymm(a, vals, vecs, work, olwork, err)
103 real(real64),
intent(inout),
dimension(:,:) :: a
104 complex(real64),
intent(out),
dimension(:) :: vals
105 complex(real64),
intent(out),
optional,
dimension(:,:) :: vecs
106 real(real64),
intent(out),
pointer,
optional,
dimension(:) :: work
107 integer(int32),
intent(out),
optional :: olwork
108 class(errors),
intent(inout),
optional,
target :: err
111 real(real64),
parameter :: zero = 0.0d0
112 real(real64),
parameter :: two = 2.0d0
115 character :: jobvl, jobvr
116 integer(int32) :: i, j, jp1, n, n1, n2a, n2b, n3a, n3b, istat, flag, &
119 real(real64),
dimension(1) :: dummy, temp
120 real(real64),
dimension(1,1) :: dummy_mtx
121 real(real64),
pointer,
dimension(:) :: wr, wi, wptr, w
122 real(real64),
pointer,
dimension(:,:) :: vr
123 real(real64),
allocatable,
target,
dimension(:) :: wrk
124 class(errors),
pointer :: errmgr
125 type(errors),
target :: deferr
126 character(len = :),
allocatable :: errmsg
130 if (
present(vecs))
then
136 eps = two * epsilon(eps)
137 if (
present(err))
then
145 if (
size(a, 2) /= n)
then
147 else if (
size(vals) /= n)
then
149 else if (
present(vecs))
then
150 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n)
then
156 allocate(
character(len = 256) :: errmsg)
157 write(errmsg, 100)
"Input number ", flag, &
158 " is not sized correctly."
159 call errmgr%report_error(
"eigen_asymm", trim(errmsg), &
165 call dgeev(jobvl, jobvr, n, a, n, dummy, dummy, dummy_mtx, n, &
166 dummy_mtx, n, temp, -1, flag)
167 lwork1 = int(temp(1), int32)
168 if (
present(vecs))
then
169 lwork = lwork1 + 2 * n + n * n
171 lwork = lwork1 + 2 * n
173 if (
present(olwork))
then
179 if (
present(work))
then
180 if (
size(work) < lwork)
then
182 call errmgr%report_error(
"eigen_asymm", &
183 "Incorrectly sized input array WORK, argument 5.", &
187 wptr => work(1:lwork)
189 allocate(wrk(lwork), stat = istat)
192 call errmgr%report_error(
"eigen_asymm", &
193 "Insufficient memory available.", &
194 la_out_of_memory_error)
205 n3b = n3a + lwork1 - 1
213 if (
present(vecs))
then
215 vr(1:n,1:n) => wptr(n3b+1:lwork)
218 call dgeev(jobvl, jobvr, n, a, n, wr, wi, dummy_mtx, n, vr, n, &
223 call errmgr%report_error(
"eigen_asymm", &
224 "The algorithm failed to converge.", la_convergence_error)
231 if (abs(wi(j)) < eps)
then
233 vals(j) = cmplx(wr(j), zero, real64)
235 vecs(i,j) = cmplx(vr(i,j), zero, real64)
240 vals(j) = cmplx(wr(j), wi(j), real64)
241 vals(jp1) = conjg(vals(j))
243 vecs(i,j) = cmplx(vr(i,j), vr(i,jp1), real64)
244 vecs(i,jp1) = conjg(vecs(i,j))
257 call dgeev(jobvl, jobvr, n, a, n, wr, wi, dummy_mtx, n, &
258 dummy_mtx, n, w, lwork1, flag)
262 call errmgr%report_error(
"eigen_asymm", &
263 "The algorithm failed to converge.", la_convergence_error)
269 vals(i) = cmplx(wr(i), wi(i), real64)
278 module subroutine eigen_gen(a, b, alpha, beta, vecs, work, olwork, err)
280 real(real64),
intent(inout),
dimension(:,:) :: a, b
281 complex(real64),
intent(out),
dimension(:) :: alpha
282 real(real64),
intent(out),
optional,
dimension(:) :: beta
283 complex(real64),
intent(out),
optional,
dimension(:,:) :: vecs
284 real(real64),
intent(out),
optional,
pointer,
dimension(:) :: work
285 integer(int32),
intent(out),
optional :: olwork
286 class(errors),
intent(inout),
optional,
target :: err
289 real(real64),
parameter :: zero = 0.0d0
290 real(real64),
parameter :: two = 2.0d0
293 character :: jobvl, jobvr
294 integer(int32) :: i, j, jp1, n, n1, n2a, n2b, n3a, n3b, n4a, n4b, &
295 istat, flag, lwork, lwork1
296 real(real64),
dimension(1) :: temp
297 real(real64),
dimension(1,1) :: dummy
298 real(real64),
pointer,
dimension(:) :: wptr, w, alphar, alphai, bptr
299 real(real64),
pointer,
dimension(:,:) :: vr
300 real(real64),
allocatable,
target,
dimension(:) :: wrk
302 class(errors),
pointer :: errmgr
303 type(errors),
target :: deferr
304 character(len = :),
allocatable :: errmsg
309 if (
present(vecs))
then
315 eps = two * epsilon(eps)
316 if (
present(err))
then
324 if (
size(a, 2) /= n)
then
326 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
328 else if (
size(alpha) /= n)
then
330 else if (
present(beta))
then
331 if (
size(beta) /= n) flag = 4
332 else if (
present(vecs))
then
333 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n) flag = 5
337 allocate(
character(len = 256) :: errmsg)
338 write(errmsg, 100)
"Input number ", flag, &
339 " is not sized correctly."
340 call errmgr%report_error(
"eigen_gen", trim(errmsg), &
346 call dggev(jobvl, jobvr, n, a, n, b, n, temp, temp, temp, dummy, n, &
347 dummy, n, temp, -1, flag)
348 lwork1 = int(temp(1), int32)
349 lwork = lwork1 + 2 * n
350 if (.not.
present(beta))
then
353 if (
present(vecs))
then
354 lwork = lwork + n * n
356 if (
present(olwork))
then
362 if (
present(work))
then
363 if (
size(work) < lwork)
then
365 call errmgr%report_error(
"eigen_gen", &
366 "Incorrectly sized input array WORK, argument 5.", &
370 wptr => work(1:lwork)
372 allocate(wrk(lwork), stat = istat)
375 call errmgr%report_error(
"eigen_gen", &
376 "Insufficient memory available.", &
377 la_out_of_memory_error)
388 n3b = n3a + lwork1 - 1
391 alphai => wptr(n2a:n2b)
393 if (.not.
present(beta))
then
396 bptr => wptr(n4a:n4b)
400 if (
present(vecs))
then
402 vr(1:n,1:n) => wptr(n4b+1:lwork)
405 if (
present(beta))
then
406 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
407 beta, dummy, n, vr, n, w, lwork1, flag)
409 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
410 bptr, dummy, n, vr, n, w, lwork1, flag)
415 call errmgr%report_error(
"eigen_gen", &
416 "The algorithm failed to converge.", la_convergence_error)
423 if (abs(alphai(j)) < eps)
then
425 alpha(j) = cmplx(alphar(j), zero, real64)
427 vecs(i,j) = cmplx(vr(i,j), zero, real64)
432 alpha(j) = cmplx(alphar(j), alphai(j), real64)
433 alpha(jp1) = cmplx(alphar(jp1), alphai(jp1), real64)
435 vecs(i,j) = cmplx(vr(i,j), vr(i,jp1), real64)
436 vecs(i,jp1) = conjg(vecs(i,j))
447 if (.not.
present(beta)) alpha = alpha / bptr
450 if (
present(beta))
then
451 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
452 beta, dummy, n, dummy, n, w, lwork1, flag)
454 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
455 bptr, dummy, n, dummy, n, w, lwork1, flag)
460 call errmgr%report_error(
"eigen_gen", &
461 "The algorithm failed to converge.", la_convergence_error)
467 alpha(i) = cmplx(alphar(i), alphai(i), real64)
469 if (.not.
present(beta)) alpha = alpha / bptr
477 module subroutine eigen_cmplx(a, vals, vecs, work, olwork, rwork, err)
479 complex(real64),
intent(inout),
dimension(:,:) :: a
480 complex(real64),
intent(out),
dimension(:) :: vals
481 complex(real64),
intent(out),
optional,
dimension(:,:) :: vecs
482 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
483 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
484 integer(int32),
intent(out),
optional :: olwork
485 class(errors),
intent(inout),
optional,
target :: err
488 character :: jobvl, jobvr
489 character(len = :),
allocatable :: errmsg
490 integer(int32) :: n, flag, lwork, lrwork
491 real(real64) :: rdummy(1)
492 complex(real64) :: temp(1), dummy(1), dummy_mtx(1,1)
493 complex(real64),
allocatable,
target,
dimension(:) :: wrk
494 complex(real64),
pointer,
dimension(:) :: wptr
495 real(real64),
allocatable,
target,
dimension(:) :: rwrk
496 real(real64),
pointer,
dimension(:) :: rwptr
497 class(errors),
pointer :: errmgr
498 type(errors),
target :: deferr
501 if (
present(err))
then
507 if (
present(vecs))
then
517 if (
size(a, 2) /= n)
then
519 else if (
size(vals) /= n)
then
521 else if (
present(vecs))
then
522 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n)
then
528 allocate(
character(len = 256) :: errmsg)
529 write(errmsg, 100)
"Input number ", flag, &
530 " is not sized correctly."
531 call errmgr%report_error(
"eigen_cmplx", trim(errmsg), &
537 call zgeev(jobvl, jobvr, n, a, n, dummy, dummy_mtx, n, dummy_mtx, n, temp, &
539 lwork = int(temp(1), int32)
540 if (
present(olwork))
then
546 if (
present(work))
then
547 if (
size(work) < lwork)
then
549 call errmgr%report_error(
"eigen_cmplx", &
550 "Incorrectly sized input array WORK.", &
556 allocate(wrk(lwork), stat = flag)
559 call errmgr%report_error(
"eigen_cmplx", &
560 "Insufficient memory available.", &
561 la_out_of_memory_error)
567 if (
present(rwork))
then
568 if (
size(work) < lrwork)
then
570 call errmgr%report_error(
"eigen_cmplx", &
571 "Incorrectly sized input array RWORK.", &
577 allocate(rwrk(lrwork), stat = flag)
580 call errmgr%report_error(
"eigen_cmplx", &
581 "Insufficient memory available.", &
582 la_out_of_memory_error)
589 if (
present(vecs))
then
590 call zgeev(jobvl, jobvr, n, a, n, vals, dummy_mtx, n, vecs, n, &
591 wptr, lwork, rwptr, flag)
593 call zgeev(jobvl, jobvr, n, a, n, vals, dummy_mtx, n, dummy_mtx, n, &
594 wptr, lwork, rwptr, flag)
598 call errmgr%report_error(
"eigen_cmplx", &
599 "The algorithm failed to converge.", &
600 la_convergence_error)
Provides a set of common linear algebra routines.