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