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 = 128) :: errmsg
36 if (
present(err))
then
44 if (
size(a, 2) /= n)
then
46 else if (
size(vals) /= n)
then
51 write(errmsg, 100)
"Input number ", flag, &
52 " is not sized correctly."
53 call errmgr%report_error(
"eigen_symm", trim(errmsg), &
59 call dsyev(jobz,
'L', n, a, n, vals, temp, -1, flag)
60 lwork = int(temp(1), int32)
61 if (
present(olwork))
then
67 if (
present(work))
then
68 if (
size(work) < lwork)
then
70 call errmgr%report_error(
"eigen_symm", &
71 "Incorrectly sized input array WORK, argument 5.", &
77 allocate(wrk(lwork), stat = istat)
80 call errmgr%report_error(
"eigen_symm", &
81 "Insufficient memory available.", &
82 la_out_of_memory_error)
89 call dsyev(jobz,
'L', n, a, n, vals, wptr, lwork, flag)
91 call errmgr%report_error(
"eigen_symm", &
92 "The algorithm failed to converge.", la_convergence_error)
100 module subroutine eigen_asymm(a, vals, vecs, work, olwork, err)
102 real(real64),
intent(inout),
dimension(:,:) :: a
103 complex(real64),
intent(out),
dimension(:) :: vals
104 complex(real64),
intent(out),
optional,
dimension(:,:) :: vecs
105 real(real64),
intent(out),
pointer,
optional,
dimension(:) :: work
106 integer(int32),
intent(out),
optional :: olwork
107 class(errors),
intent(inout),
optional,
target :: err
110 real(real64),
parameter :: zero = 0.0d0
111 real(real64),
parameter :: two = 2.0d0
114 character :: jobvl, jobvr
115 integer(int32) :: i, j, jp1, n, n1, n2a, n2b, n3a, n3b, istat, flag, &
118 real(real64),
dimension(1) :: dummy, temp
119 real(real64),
dimension(1,1) :: dummy_mtx
120 real(real64),
pointer,
dimension(:) :: wr, wi, wptr, w
121 real(real64),
pointer,
dimension(:,:) :: vr
122 real(real64),
allocatable,
target,
dimension(:) :: wrk
123 class(errors),
pointer :: errmgr
124 type(errors),
target :: deferr
125 character(len = 128) :: errmsg
129 if (
present(vecs))
then
135 eps = two * epsilon(eps)
136 if (
present(err))
then
144 if (
size(a, 2) /= n)
then
146 else if (
size(vals) /= n)
then
148 else if (
present(vecs))
then
149 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n)
then
155 write(errmsg, 100)
"Input number ", flag, &
156 " is not sized correctly."
157 call errmgr%report_error(
"eigen_asymm", trim(errmsg), &
163 call dgeev(jobvl, jobvr, n, a, n, dummy, dummy, dummy_mtx, n, &
164 dummy_mtx, n, temp, -1, flag)
165 lwork1 = int(temp(1), int32)
166 if (
present(vecs))
then
167 lwork = lwork1 + 2 * n + n * n
169 lwork = lwork1 + 2 * n
171 if (
present(olwork))
then
177 if (
present(work))
then
178 if (
size(work) < lwork)
then
180 call errmgr%report_error(
"eigen_asymm", &
181 "Incorrectly sized input array WORK, argument 5.", &
185 wptr => work(1:lwork)
187 allocate(wrk(lwork), stat = istat)
190 call errmgr%report_error(
"eigen_asymm", &
191 "Insufficient memory available.", &
192 la_out_of_memory_error)
203 n3b = n3a + lwork1 - 1
211 if (
present(vecs))
then
213 vr(1:n,1:n) => wptr(n3b+1:lwork)
216 call dgeev(jobvl, jobvr, n, a, n, wr, wi, dummy_mtx, n, vr, n, &
221 call errmgr%report_error(
"eigen_asymm", &
222 "The algorithm failed to converge.", la_convergence_error)
229 if (abs(wi(j)) < eps)
then
231 vals(j) = cmplx(wr(j), zero, real64)
233 vecs(i,j) = cmplx(vr(i,j), zero, real64)
238 vals(j) = cmplx(wr(j), wi(j), real64)
239 vals(jp1) = conjg(vals(j))
241 vecs(i,j) = cmplx(vr(i,j), vr(i,jp1), real64)
242 vecs(i,jp1) = conjg(vecs(i,j))
255 call dgeev(jobvl, jobvr, n, a, n, wr, wi, dummy_mtx, n, &
256 dummy_mtx, n, w, lwork1, flag)
260 call errmgr%report_error(
"eigen_asymm", &
261 "The algorithm failed to converge.", la_convergence_error)
267 vals(i) = cmplx(wr(i), wi(i), real64)
276 module subroutine eigen_gen(a, b, alpha, beta, vecs, work, olwork, err)
278 real(real64),
intent(inout),
dimension(:,:) :: a, b
279 complex(real64),
intent(out),
dimension(:) :: alpha
280 real(real64),
intent(out),
optional,
dimension(:) :: beta
281 complex(real64),
intent(out),
optional,
dimension(:,:) :: vecs
282 real(real64),
intent(out),
optional,
pointer,
dimension(:) :: work
283 integer(int32),
intent(out),
optional :: olwork
284 class(errors),
intent(inout),
optional,
target :: err
287 real(real64),
parameter :: zero = 0.0d0
288 real(real64),
parameter :: two = 2.0d0
291 character :: jobvl, jobvr
292 integer(int32) :: i, j, jp1, n, n1, n2a, n2b, n3a, n3b, n4a, n4b, &
293 istat, flag, lwork, lwork1
294 real(real64),
dimension(1) :: temp
295 real(real64),
dimension(1,1) :: dummy
296 real(real64),
pointer,
dimension(:) :: wptr, w, alphar, alphai, bptr
297 real(real64),
pointer,
dimension(:,:) :: vr
298 real(real64),
allocatable,
target,
dimension(:) :: wrk
300 class(errors),
pointer :: errmgr
301 type(errors),
target :: deferr
302 character(len = 128) :: errmsg
307 if (
present(vecs))
then
313 eps = two * epsilon(eps)
314 if (
present(err))
then
322 if (
size(a, 2) /= n)
then
324 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
326 else if (
size(alpha) /= n)
then
328 else if (
present(beta))
then
329 if (
size(beta) /= n) flag = 4
330 else if (
present(vecs))
then
331 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n) flag = 5
335 write(errmsg, 100)
"Input number ", flag, &
336 " is not sized correctly."
337 call errmgr%report_error(
"eigen_gen", trim(errmsg), &
343 call dggev(jobvl, jobvr, n, a, n, b, n, temp, temp, temp, dummy, n, &
344 dummy, n, temp, -1, flag)
345 lwork1 = int(temp(1), int32)
346 lwork = lwork1 + 2 * n
347 if (.not.
present(beta))
then
350 if (
present(vecs))
then
351 lwork = lwork + n * n
353 if (
present(olwork))
then
359 if (
present(work))
then
360 if (
size(work) < lwork)
then
362 call errmgr%report_error(
"eigen_gen", &
363 "Incorrectly sized input array WORK, argument 5.", &
367 wptr => work(1:lwork)
369 allocate(wrk(lwork), stat = istat)
372 call errmgr%report_error(
"eigen_gen", &
373 "Insufficient memory available.", &
374 la_out_of_memory_error)
385 n3b = n3a + lwork1 - 1
388 alphai => wptr(n2a:n2b)
390 if (.not.
present(beta))
then
393 bptr => wptr(n4a:n4b)
397 if (
present(vecs))
then
399 vr(1:n,1:n) => wptr(n4b+1:lwork)
402 if (
present(beta))
then
403 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
404 beta, dummy, n, vr, n, w, lwork1, flag)
406 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
407 bptr, dummy, n, vr, n, w, lwork1, flag)
412 call errmgr%report_error(
"eigen_gen", &
413 "The algorithm failed to converge.", la_convergence_error)
420 if (abs(alphai(j)) < eps)
then
422 alpha(j) = cmplx(alphar(j), zero, real64)
424 vecs(i,j) = cmplx(vr(i,j), zero, real64)
429 alpha(j) = cmplx(alphar(j), alphai(j), real64)
430 alpha(jp1) = cmplx(alphar(jp1), alphai(jp1), real64)
432 vecs(i,j) = cmplx(vr(i,j), vr(i,jp1), real64)
433 vecs(i,jp1) = conjg(vecs(i,j))
444 if (.not.
present(beta)) alpha = alpha / bptr
447 if (
present(beta))
then
448 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
449 beta, dummy, n, dummy, n, w, lwork1, flag)
451 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
452 bptr, dummy, n, dummy, n, w, lwork1, flag)
457 call errmgr%report_error(
"eigen_gen", &
458 "The algorithm failed to converge.", la_convergence_error)
464 alpha(i) = cmplx(alphar(i), alphai(i), real64)
466 if (.not.
present(beta)) alpha = alpha / bptr
474 module subroutine eigen_cmplx(a, vals, vecs, work, olwork, rwork, err)
476 complex(real64),
intent(inout),
dimension(:,:) :: a
477 complex(real64),
intent(out),
dimension(:) :: vals
478 complex(real64),
intent(out),
optional,
dimension(:,:) :: vecs
479 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
480 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
481 integer(int32),
intent(out),
optional :: olwork
482 class(errors),
intent(inout),
optional,
target :: err
485 character :: jobvl, jobvr
486 character(len = 128) :: errmsg
487 integer(int32) :: n, flag, lwork, lrwork
488 real(real64) :: rdummy(1)
489 complex(real64) :: temp(1), dummy(1), dummy_mtx(1,1)
490 complex(real64),
allocatable,
target,
dimension(:) :: wrk
491 complex(real64),
pointer,
dimension(:) :: wptr
492 real(real64),
allocatable,
target,
dimension(:) :: rwrk
493 real(real64),
pointer,
dimension(:) :: rwptr
494 class(errors),
pointer :: errmgr
495 type(errors),
target :: deferr
498 if (
present(err))
then
504 if (
present(vecs))
then
514 if (
size(a, 2) /= n)
then
516 else if (
size(vals) /= n)
then
518 else if (
present(vecs))
then
519 if (
size(vecs, 1) /= n .or.
size(vecs, 2) /= n)
then
525 write(errmsg, 100)
"Input number ", flag, &
526 " is not sized correctly."
527 call errmgr%report_error(
"eigen_cmplx", trim(errmsg), &
533 call zgeev(jobvl, jobvr, n, a, n, dummy, dummy_mtx, n, dummy_mtx, n, temp, &
535 lwork = int(temp(1), int32)
536 if (
present(olwork))
then
542 if (
present(work))
then
543 if (
size(work) < lwork)
then
545 call errmgr%report_error(
"eigen_cmplx", &
546 "Incorrectly sized input array WORK.", &
552 allocate(wrk(lwork), stat = flag)
555 call errmgr%report_error(
"eigen_cmplx", &
556 "Insufficient memory available.", &
557 la_out_of_memory_error)
563 if (
present(rwork))
then
564 if (
size(work) < lrwork)
then
566 call errmgr%report_error(
"eigen_cmplx", &
567 "Incorrectly sized input array RWORK.", &
573 allocate(rwrk(lrwork), stat = flag)
576 call errmgr%report_error(
"eigen_cmplx", &
577 "Insufficient memory available.", &
578 la_out_of_memory_error)
585 if (
present(vecs))
then
586 call zgeev(jobvl, jobvr, n, a, n, vals, dummy_mtx, n, vecs, n, &
587 wptr, lwork, rwptr, flag)
589 call zgeev(jobvl, jobvr, n, a, n, vals, dummy_mtx, n, dummy_mtx, n, &
590 wptr, lwork, rwptr, flag)
594 call errmgr%report_error(
"eigen_cmplx", &
595 "The algorithm failed to converge.", &
596 la_convergence_error)
Provides a set of common linear algebra routines.