linalg 1.7.2
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
linalg_eigen.f90
1! linalg_eigen.f90
2
7submodule(linalg) linalg_eigen
8contains
9! ------------------------------------------------------------------------------
10 module subroutine eigen_symm(vecs, a, vals, work, olwork, err)
11 ! Arguments
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
18
19 ! Local Variables
20 character :: jobz
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
28
29 ! Initialization
30 n = size(a, 1)
31 if (vecs) then
32 jobz = 'V'
33 else
34 jobz = 'N'
35 end if
36 if (present(err)) then
37 errmgr => err
38 else
39 errmgr => deferr
40 end if
41
42 ! Input Check
43 flag = 0
44 if (size(a, 2) /= n) then
45 flag = 2
46 else if (size(vals) /= n) then
47 flag = 3
48 end if
49 if (flag /= 0) then
50 ! ERROR: One of the input arrays is not sized correctly
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), &
55 la_array_size_error)
56 return
57 end if
58
59 ! Workspace Query
60 call dsyev(jobz, 'L', n, a, n, vals, temp, -1, flag)
61 lwork = int(temp(1), int32)
62 if (present(olwork)) then
63 olwork = lwork
64 return
65 end if
66
67 ! Local Memory Allocation
68 if (present(work)) then
69 if (size(work) < lwork) then
70 ! ERROR: WORK not sized correctly
71 call errmgr%report_error("eigen_symm", &
72 "Incorrectly sized input array WORK, argument 5.", &
73 la_array_size_error)
74 return
75 end if
76 wptr => work(1:lwork)
77 else
78 allocate(wrk(lwork), stat = istat)
79 if (istat /= 0) then
80 ! ERROR: Out of memory
81 call errmgr%report_error("eigen_symm", &
82 "Insufficient memory available.", &
83 la_out_of_memory_error)
84 return
85 end if
86 wptr => wrk
87 end if
88
89 ! Process
90 call dsyev(jobz, 'L', n, a, n, vals, wptr, lwork, flag)
91 if (flag > 0) then
92 call errmgr%report_error("eigen_symm", &
93 "The algorithm failed to converge.", la_convergence_error)
94 end if
95
96 ! Formatting
97100 format(a, i0, a)
98 end subroutine
99
100! ------------------------------------------------------------------------------
101 module subroutine eigen_asymm(a, vals, vecs, work, olwork, err)
102 ! Arguments
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
109
110 ! Parameters
111 real(real64), parameter :: zero = 0.0d0
112 real(real64), parameter :: two = 2.0d0
113
114 ! Local Variables
115 character :: jobvl, jobvr
116 integer(int32) :: i, j, jp1, n, n1, n2a, n2b, n3a, n3b, istat, flag, &
117 lwork, lwork1
118 real(real64) :: eps
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
127
128 ! Initialization
129 jobvl = 'N'
130 if (present(vecs)) then
131 jobvr = 'V'
132 else
133 jobvr = 'N'
134 end if
135 n = size(a, 1)
136 eps = two * epsilon(eps)
137 if (present(err)) then
138 errmgr => err
139 else
140 errmgr => deferr
141 end if
142
143 ! Input Check
144 flag = 0
145 if (size(a, 2) /= n) then
146 flag = 1
147 else if (size(vals) /= n) then
148 flag = 2
149 else if (present(vecs)) then
150 if (size(vecs, 1) /= n .or. size(vecs, 2) /= n) then
151 flag = 3
152 end if
153 end if
154 if (flag /= 0) then
155 ! ERROR: One of the input arrays is not sized correctly
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), &
160 la_array_size_error)
161 return
162 end if
163
164 ! Workspace Query
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
170 else
171 lwork = lwork1 + 2 * n
172 end if
173 if (present(olwork)) then
174 olwork = lwork
175 return
176 end if
177
178 ! Local Memory Allocation
179 if (present(work)) then
180 if (size(work) < lwork) then
181 ! ERROR: WORK not sized correctly
182 call errmgr%report_error("eigen_asymm", &
183 "Incorrectly sized input array WORK, argument 5.", &
184 la_array_size_error)
185 return
186 end if
187 wptr => work(1:lwork)
188 else
189 allocate(wrk(lwork), stat = istat)
190 if (istat /= 0) then
191 ! ERROR: Out of memory
192 call errmgr%report_error("eigen_asymm", &
193 "Insufficient memory available.", &
194 la_out_of_memory_error)
195 return
196 end if
197 wptr => wrk
198 end if
199
200 ! Locate each array within the workspace array
201 n1 = n
202 n2a = n1 + 1
203 n2b = n2a + n - 1
204 n3a = n2b + 1
205 n3b = n3a + lwork1 - 1
206
207 ! Assign pointers
208 wr => wptr(1:n1)
209 wi => wptr(n2a:n2b)
210 w => wptr(n3a:n3b)
211
212 ! Process
213 if (present(vecs)) then
214 ! Assign a pointer to the eigenvector matrix
215 vr(1:n,1:n) => wptr(n3b+1:lwork)
216
217 ! Compute the eigenvectors and eigenvalues
218 call dgeev(jobvl, jobvr, n, a, n, wr, wi, dummy_mtx, n, vr, n, &
219 w, lwork1, flag)
220
221 ! Check for convergence
222 if (flag > 0) then
223 call errmgr%report_error("eigen_asymm", &
224 "The algorithm failed to converge.", la_convergence_error)
225 return
226 end if
227
228 ! Store the eigenvalues and eigenvectors
229 j = 1
230 do while (j <= n)
231 if (abs(wi(j)) < eps) then
232 ! We've got a real-valued eigenvalue
233 vals(j) = cmplx(wr(j), zero, real64)
234 do i = 1, n
235 vecs(i,j) = cmplx(vr(i,j), zero, real64)
236 end do
237 else
238 ! We've got a complex cojugate pair of eigenvalues
239 jp1 = j + 1
240 vals(j) = cmplx(wr(j), wi(j), real64)
241 vals(jp1) = conjg(vals(j))
242 do i = 1, n
243 vecs(i,j) = cmplx(vr(i,j), vr(i,jp1), real64)
244 vecs(i,jp1) = conjg(vecs(i,j))
245 end do
246
247 ! Increment j and continue the loop
248 j = j + 2
249 cycle
250 end if
251
252 ! Increment j
253 j = j + 1
254 end do
255 else
256 ! Compute just the eigenvalues
257 call dgeev(jobvl, jobvr, n, a, n, wr, wi, dummy_mtx, n, &
258 dummy_mtx, n, w, lwork1, flag)
259
260 ! Check for convergence
261 if (flag > 0) then
262 call errmgr%report_error("eigen_asymm", &
263 "The algorithm failed to converge.", la_convergence_error)
264 return
265 end if
266
267 ! Store the eigenvalues
268 do i = 1, n
269 vals(i) = cmplx(wr(i), wi(i), real64)
270 end do
271 end if
272
273 ! Formatting
274100 format(a, i0, a)
275 end subroutine
276
277! ------------------------------------------------------------------------------
278 module subroutine eigen_gen(a, b, alpha, beta, vecs, work, olwork, err)
279 ! Arguments
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
287
288 ! Parameters
289 real(real64), parameter :: zero = 0.0d0
290 real(real64), parameter :: two = 2.0d0
291
292 ! Local Variables
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
301 real(real64) :: eps
302 class(errors), pointer :: errmgr
303 type(errors), target :: deferr
304 character(len = :), allocatable :: errmsg
305
306 ! Initialization
307 jobvl = 'N'
308 jobvr = 'N'
309 if (present(vecs)) then
310 jobvr = 'V'
311 else
312 jobvr = 'N'
313 end if
314 n = size(a, 1)
315 eps = two * epsilon(eps)
316 if (present(err)) then
317 errmgr => err
318 else
319 errmgr => deferr
320 end if
321
322 ! Input Check
323 flag = 0
324 if (size(a, 2) /= n) then
325 flag = 1
326 else if (size(b, 1) /= n .or. size(b, 2) /= n) then
327 flag = 2
328 else if (size(alpha) /= n) then
329 flag = 3
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
334 end if
335 if (flag /= 0) then
336 ! ERROR: One of the input arrays is not sized correctly
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), &
341 la_array_size_error)
342 return
343 end if
344
345 ! Workspace Query
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
351 lwork = lwork + n
352 end if
353 if (present(vecs)) then
354 lwork = lwork + n * n
355 end if
356 if (present(olwork)) then
357 olwork = lwork
358 return
359 end if
360
361 ! Local Memory Allocation
362 if (present(work)) then
363 if (size(work) < lwork) then
364 ! ERROR: WORK not sized correctly
365 call errmgr%report_error("eigen_gen", &
366 "Incorrectly sized input array WORK, argument 5.", &
367 la_array_size_error)
368 return
369 end if
370 wptr => work(1:lwork)
371 else
372 allocate(wrk(lwork), stat = istat)
373 if (istat /= 0) then
374 ! ERROR: Out of memory
375 call errmgr%report_error("eigen_gen", &
376 "Insufficient memory available.", &
377 la_out_of_memory_error)
378 return
379 end if
380 wptr => wrk
381 end if
382
383 ! Locate each array within the workspace array & assign pointers
384 n1 = n
385 n2a = n1 + 1
386 n2b = n2a + n - 1
387 n3a = n2b + 1
388 n3b = n3a + lwork1 - 1
389 n4b = n3b
390 alphar => wptr(1:n1)
391 alphai => wptr(n2a:n2b)
392 w => wptr(n3a:n3b)
393 if (.not.present(beta)) then
394 n4a = n3b + 1
395 n4b = n4a + n - 1
396 bptr => wptr(n4a:n4b)
397 end if
398
399 ! Process
400 if (present(vecs)) then
401 ! Assign a pointer to the eigenvector matrix
402 vr(1:n,1:n) => wptr(n4b+1:lwork)
403
404 ! Compute the eigenvalues and eigenvectors
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)
408 else
409 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
410 bptr, dummy, n, vr, n, w, lwork1, flag)
411 end if
412
413 ! Check for convergence
414 if (flag > 0) then
415 call errmgr%report_error("eigen_gen", &
416 "The algorithm failed to converge.", la_convergence_error)
417 return
418 end if
419
420 ! Store the eigenvalues and eigenvectors
421 j = 1
422 do while (j <= n)
423 if (abs(alphai(j)) < eps) then
424 ! Real-Valued
425 alpha(j) = cmplx(alphar(j), zero, real64)
426 do i = 1, n
427 vecs(i,j) = cmplx(vr(i,j), zero, real64)
428 end do
429 else
430 ! Complex-Valued
431 jp1 = j + 1
432 alpha(j) = cmplx(alphar(j), alphai(j), real64)
433 alpha(jp1) = cmplx(alphar(jp1), alphai(jp1), real64)
434 do i = 1, n
435 vecs(i,j) = cmplx(vr(i,j), vr(i,jp1), real64)
436 vecs(i,jp1) = conjg(vecs(i,j))
437 end do
438
439 ! Increment j and continue
440 j = j + 2
441 cycle
442 end if
443
444 ! Increment j
445 j = j + 1
446 end do
447 if (.not.present(beta)) alpha = alpha / bptr
448 else
449 ! Compute just the eigenvalues
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)
453 else
454 call dggev(jobvl, jobvr, n, a, n, b, n, alphar, alphai, &
455 bptr, dummy, n, dummy, n, w, lwork1, flag)
456 end if
457
458 ! Check for convergence
459 if (flag > 0) then
460 call errmgr%report_error("eigen_gen", &
461 "The algorithm failed to converge.", la_convergence_error)
462 return
463 end if
464
465 ! Store the eigenvalues
466 do i = 1, n
467 alpha(i) = cmplx(alphar(i), alphai(i), real64)
468 end do
469 if (.not.present(beta)) alpha = alpha / bptr
470 end if
471
472 ! Formatting
473100 format(a, i0, a)
474 end subroutine
475
476! ------------------------------------------------------------------------------
477 module subroutine eigen_cmplx(a, vals, vecs, work, olwork, rwork, err)
478 ! Arguments
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
486
487 ! Local Variables
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
499
500 ! Initialization
501 if (present(err)) then
502 errmgr => err
503 else
504 errmgr => deferr
505 end if
506 jobvl = 'N'
507 if (present(vecs)) then
508 jobvr = 'V'
509 else
510 jobvr = 'N'
511 end if
512 n = size(a, 1)
513 lrwork = 2 * n
514
515 ! Input Check
516 flag = 0
517 if (size(a, 2) /= n) then
518 flag = 1
519 else if (size(vals) /= n) then
520 flag = 2
521 else if (present(vecs)) then
522 if (size(vecs, 1) /= n .or. size(vecs, 2) /= n) then
523 flag = 3
524 end if
525 end if
526 if (flag /= 0) then
527 ! ERROR: One of the input arrays is not sized correctly
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), &
532 la_array_size_error)
533 return
534 end if
535
536 ! Workspace Query
537 call zgeev(jobvl, jobvr, n, a, n, dummy, dummy_mtx, n, dummy_mtx, n, temp, &
538 -1, rdummy, flag)
539 lwork = int(temp(1), int32)
540 if (present(olwork)) then
541 olwork = lwork
542 return
543 end if
544
545 ! Local Memory Allocation
546 if (present(work)) then
547 if (size(work) < lwork) then
548 ! ERROR: WORK not sized correctly
549 call errmgr%report_error("eigen_cmplx", &
550 "Incorrectly sized input array WORK.", &
551 la_array_size_error)
552 return
553 end if
554 wptr => work
555 else
556 allocate(wrk(lwork), stat = flag)
557 if (flag /= 0) then
558 ! ERROR: Out of memory
559 call errmgr%report_error("eigen_cmplx", &
560 "Insufficient memory available.", &
561 la_out_of_memory_error)
562 return
563 end if
564 wptr => wrk
565 end if
566
567 if (present(rwork)) then
568 if (size(work) < lrwork) then
569 ! ERROR: RWORK not sized correctly
570 call errmgr%report_error("eigen_cmplx", &
571 "Incorrectly sized input array RWORK.", &
572 la_array_size_error)
573 return
574 end if
575 rwptr => rwork
576 else
577 allocate(rwrk(lrwork), stat = flag)
578 if (flag /= 0) then
579 ! ERROR: Out of memory
580 call errmgr%report_error("eigen_cmplx", &
581 "Insufficient memory available.", &
582 la_out_of_memory_error)
583 return
584 end if
585 rwptr => rwrk
586 end if
587
588 ! Process
589 if (present(vecs)) then
590 call zgeev(jobvl, jobvr, n, a, n, vals, dummy_mtx, n, vecs, n, &
591 wptr, lwork, rwptr, flag)
592 else
593 call zgeev(jobvl, jobvr, n, a, n, vals, dummy_mtx, n, dummy_mtx, n, &
594 wptr, lwork, rwptr, flag)
595 end if
596
597 if (flag > 0) then
598 call errmgr%report_error("eigen_cmplx", &
599 "The algorithm failed to converge.", &
600 la_convergence_error)
601 return
602 end if
603
604 ! Formatting
605100 format(a, i0, a)
606 end subroutine
607
608! ------------------------------------------------------------------------------
609end submodule
Provides a set of common linear algebra routines.
Definition: linalg.f90:145