linalg 1.7.4
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
linalg_basic.f90
1! linalg_basic.f90
2
3submodule(linalg) linalg_basic
4 use blas
5 use lapack
6 implicit none
7contains
8! ******************************************************************************
9! MATRIX MULTIPLICATION ROUTINES
10! ------------------------------------------------------------------------------
11 module subroutine mtx_mult_mtx(transa, transb, alpha, a, b, beta, c, err)
12 ! Arguments
13 logical, intent(in) :: transa, transb
14 real(real64), intent(in) :: alpha, beta
15 real(real64), intent(in), dimension(:,:) :: a, b
16 real(real64), intent(inout), dimension(:,:) :: c
17 class(errors), intent(inout), optional, target :: err
18
19 ! Parameters
20 real(real64), parameter :: zero = 0.0d0
21 real(real64), parameter :: one = 1.0d0
22
23 ! Local Variables
24 character :: ta, tb
25 integer(int32) :: m, n, k, lda, ldb, flag
26 class(errors), pointer :: errmgr
27 type(errors), target :: deferr
28 character(len = :), allocatable :: errmsg
29
30 ! Initialization
31 m = size(c, 1)
32 n = size(c, 2)
33 if (transa) then ! K = # of columns in op(A) (# of rows in op(B))
34 k = size(a, 1)
35 ta = 'T'
36 lda = k
37 else
38 k = size(a, 2)
39 ta = 'N'
40 lda = m
41 end if
42 if (transb) then
43 tb = 'T'
44 ldb = n
45 else
46 tb = 'N'
47 ldb = k
48 end if
49 if (present(err)) then
50 errmgr => err
51 else
52 errmgr => deferr
53 end if
54
55 ! Input Check
56 flag = 0
57 if (transa) then
58 if (size(a, 2) /= m) flag = 4
59 else
60 if (size(a, 1) /= m) flag = 4
61 end if
62 if (transb) then
63 if (size(b, 2) /= k .or. size(b, 1) /= n) flag = 5
64 else
65 if (size(b, 1) /= k .or. size(b, 2) /= n) flag = 5
66 end if
67 if (flag /= 0) then
68 ! ERROR: Matrix dimensions mismatch
69 allocate(character(len = 256) :: errmsg)
70 write(errmsg, 100) &
71 "Matrix dimension mismatch. Input number ", flag, &
72 " was not sized correctly."
73 call errmgr%report_error("mtx_mult_mtx", errmsg, &
74 la_array_size_error)
75 return
76 end if
77
78 ! Call DGEMM
79 call dgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
80
81 ! Formatting
82100 format(a, i0, a)
83 end subroutine
84
85! ------------------------------------------------------------------------------
86 module subroutine mtx_mult_vec(trans, alpha, a, b, beta, c, err)
87 ! Arguments
88 logical, intent(in) :: trans
89 real(real64), intent(in) :: alpha, beta
90 real(real64), intent(in), dimension(:,:) :: a
91 real(real64), intent(in), dimension(:) :: b
92 real(real64), intent(inout), dimension(:) :: c
93 class(errors), intent(inout), optional, target :: err
94
95 ! Local Variables
96 character :: t
97 integer(int32) :: m, n, flag
98 class(errors), pointer :: errmgr
99 type(errors), target :: deferr
100 character(len = :), allocatable :: errmsg
101
102 ! Initialization
103 m = size(a, 1)
104 n = size(a, 2)
105 t = 'N'
106 if (trans) t = 'T'
107 if (present(err)) then
108 errmgr => err
109 else
110 errmgr => deferr
111 end if
112
113 ! Input Check
114 flag = 0
115 if (trans) then
116 if (size(b) /= m) then
117 flag = 4
118 else if (size(c) /= n) then
119 flag = 6
120 end if
121 else
122 if (size(b) /= n) then
123 flag = 4
124 else if (size(c) /= m) then
125 flag = 6
126 end if
127 end if
128 if (flag /= 0) then
129 ! ERROR: Matrix dimensions mismatch
130 allocate(character(len = 256) :: errmsg)
131 write(errmsg, 100) &
132 "Matrix dimension mismatch. Input number ", flag, &
133 " was not sized correctly."
134 call errmgr%report_error("mtx_mult_vec", errmsg, &
135 la_array_size_error)
136 return
137 end if
138
139 ! Call DGEMV
140 call dgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
141
142 ! Formatting
143100 format(a, i0, a)
144 end subroutine
145
146! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
147! COMPLEX VALUED VERSIONS !
148! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
149 module subroutine cmtx_mult_mtx(opa, opb, alpha, a, b, beta, c, err)
150 ! Arguments
151 integer(int32), intent(in) :: opa, opb
152 complex(real64), intent(in) :: alpha, beta
153 complex(real64), intent(in), dimension(:,:) :: a, b
154 complex(real64), intent(inout), dimension(:,:) :: c
155 class(errors), intent(inout), optional, target :: err
156
157 ! Parameters
158 real(real64), parameter :: zero = 0.0d0
159 real(real64), parameter :: one = 1.0d0
160
161 ! Local Variables
162 character :: ta, tb
163 integer(int32) :: m, n, k, lda, ldb, flag
164 class(errors), pointer :: errmgr
165 type(errors), target :: deferr
166 character(len = :), allocatable :: errmsg
167
168 ! Initialization
169 m = size(c, 1)
170 n = size(c, 2)
171 if (opa == la_transpose) then ! K = # of columns in op(A) (# of rows in op(B))
172 k = size(a, 1)
173 ta = 'T'
174 lda = k
175 else if (opa == la_hermitian_transpose) then
176 k = size(a, 1)
177 ta = 'C'
178 lda = k
179 else
180 k = size(a, 2)
181 ta = 'N'
182 lda = m
183 end if
184 if (opb == la_transpose) then
185 tb = 'T'
186 ldb = n
187 else if (opb == la_hermitian_transpose) then
188 tb = 'C'
189 ldb = n
190 else
191 tb = 'N'
192 ldb = k
193 end if
194 if (present(err)) then
195 errmgr => err
196 else
197 errmgr => deferr
198 end if
199
200 ! Input Check
201 flag = 0
202 if (opa == la_transpose .or. opa == la_hermitian_transpose) then
203 if (size(a, 2) /= m) flag = 4
204 else
205 if (size(a, 1) /= m) flag = 4
206 end if
207 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
208 if (size(b, 2) /= k .or. size(b, 1) /= n) flag = 5
209 else
210 if (size(b, 1) /= k .or. size(b, 2) /= n) flag = 5
211 end if
212 if (flag /= 0) then
213 ! ERROR: Matrix dimensions mismatch
214 allocate(character(len = 256) :: errmsg)
215 write(errmsg, 100) &
216 "Matrix dimension mismatch. Input number ", flag, &
217 " was not sized correctly."
218 call errmgr%report_error("cmtx_mult_mtx", errmsg, &
219 la_array_size_error)
220 return
221 end if
222
223 ! Call ZGEMM
224 call zgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
225
226 ! Formatting
227100 format(a, i0, a)
228 end subroutine
229
230! ------------------------------------------------------------------------------
231 module subroutine cmtx_mult_vec(opa, alpha, a, b, beta, c, err)
232 ! Arguments
233 integer(int32), intent(in) :: opa
234 complex(real64), intent(in) :: alpha, beta
235 complex(real64), intent(in), dimension(:,:) :: a
236 complex(real64), intent(in), dimension(:) :: b
237 complex(real64), intent(inout), dimension(:) :: c
238 class(errors), intent(inout), optional, target :: err
239
240 ! Local Variables
241 character :: t
242 integer(int32) :: m, n, flag
243 class(errors), pointer :: errmgr
244 type(errors), target :: deferr
245 character(len = :), allocatable :: errmsg
246
247 ! Initialization
248 m = size(a, 1)
249 n = size(a, 2)
250 if (opa == la_transpose) then
251 t = 'T'
252 else if (opa == la_hermitian_transpose) then
253 t = 'C'
254 else
255 t = 'N'
256 end if
257 if (present(err)) then
258 errmgr => err
259 else
260 errmgr => deferr
261 end if
262
263 ! Input Check
264 flag = 0
265 if (opa == la_transpose .or. opa == la_hermitian_transpose) then
266 if (size(b) /= m) then
267 flag = 4
268 else if (size(c) /= n) then
269 flag = 6
270 end if
271 else
272 if (size(b) /= n) then
273 flag = 4
274 else if (size(c) /= m) then
275 flag = 6
276 end if
277 end if
278 if (flag /= 0) then
279 ! ERROR: Matrix dimensions mismatch
280 allocate(character(len = 256) :: errmsg)
281 write(errmsg, 100) &
282 "Matrix dimension mismatch. Input number ", flag, &
283 " was not sized correctly."
284 call errmgr%report_error("cmtx_mult_vec", errmsg, &
285 la_array_size_error)
286 return
287 end if
288
289 ! Call ZGEMV
290 call zgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
291
292 ! Formatting
293100 format(a, i0, a)
294 end subroutine
295
296! ******************************************************************************
297! RANK 1 UPDATE
298! ------------------------------------------------------------------------------
299 module subroutine rank1_update_dbl(alpha, x, y, a, err)
300 ! Arguments
301 real(real64), intent(in) :: alpha
302 real(real64), intent(in), dimension(:) :: x, y
303 real(real64), intent(inout), dimension(:,:) :: a
304 class(errors), intent(inout), optional, target :: err
305
306 ! Parameters
307 real(real64), parameter :: zero = 0.0d0
308
309 ! Local Variables
310 integer(int32) :: j, m, n
311 real(real64) :: temp
312 class(errors), pointer :: errmgr
313 type(errors), target :: deferr
314
315 ! Initialization
316 m = size(x)
317 n = size(y)
318 if (present(err)) then
319 errmgr => err
320 else
321 errmgr => deferr
322 end if
323
324 ! Input Check
325 if (size(a, 1) /= m .or. size(a, 2) /= n) then
326 ! ERROR: Matrix dimension array
327 call errmgr%report_error("rank1_update_dbl", &
328 "Matrix dimension mismatch.", la_array_size_error)
329 return
330 end if
331
332 ! Process
333 do j = 1, n
334 if (y(j) /= zero) then
335 temp = alpha * y(j)
336 a(:,j) = a(:,j) + temp * x
337 end if
338 end do
339 end subroutine
340
341! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
342! COMPLEX VALUED VERSIONS !
343! xxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxxx !
344 module subroutine rank1_update_cmplx(alpha, x, y, a, err)
345 ! Arguments
346 complex(real64), intent(in) :: alpha
347 complex(real64), intent(in), dimension(:) :: x, y
348 complex(real64), intent(inout), dimension(:,:) :: a
349 class(errors), intent(inout), optional, target :: err
350
351 ! Parameters
352 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
353
354 ! Local Variables
355 integer(int32) :: j, m, n
356 complex(real64) :: temp
357 class(errors), pointer :: errmgr
358 type(errors), target :: deferr
359
360 ! Initialization
361 m = size(x)
362 n = size(y)
363 if (present(err)) then
364 errmgr => err
365 else
366 errmgr => deferr
367 end if
368
369 ! Input Check
370 if (size(a, 1) /= m .or. size(a, 2) /= n) then
371 ! ERROR: Matrix dimension array
372 call errmgr%report_error("rank1_update_cmplx", &
373 "Matrix dimension mismatch.", la_array_size_error)
374 return
375 end if
376
377 ! Process
378 do j = 1, n
379 if (y(j) /= zero) then
380 temp = alpha * conjg(y(j))
381 a(:,j) = a(:,j) + temp * x
382 end if
383 end do
384 end subroutine
385
386! ******************************************************************************
387! DIAGONAL MATRIX MULTIPLICATION ROUTINES
388! ------------------------------------------------------------------------------
389 module subroutine diag_mtx_mult_mtx(lside, trans, alpha, a, b, beta, c, err)
390 ! Arguments
391 logical, intent(in) :: lside, trans
392 real(real64) :: alpha, beta
393 real(real64), intent(in), dimension(:) :: a
394 real(real64), intent(in), dimension(:,:) :: b
395 real(real64), intent(inout), dimension(:,:) :: c
396 class(errors), intent(inout), optional, target :: err
397
398 ! Parameters
399 real(real64), parameter :: zero = 0.0d0
400 real(real64), parameter :: one = 1.0d0
401
402 ! Local Variables
403 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
404 real(real64) :: temp
405 class(errors), pointer :: errmgr
406 type(errors), target :: deferr
407 character(len = :), allocatable :: errmsg
408
409 ! Initialization
410 m = size(c, 1)
411 n = size(c, 2)
412 k = size(a)
413 nrowb = size(b, 1)
414 ncolb = size(b, 2)
415 if (present(err)) then
416 errmgr => err
417 else
418 errmgr => deferr
419 end if
420
421 ! Input Check
422 flag = 0
423 if (lside) then
424 if (k > m) then
425 flag = 4
426 else
427 if (trans) then
428 ! Compute C = alpha * A * B**T + beta * C
429 if (nrowb /= n .or. ncolb < k) flag = 5
430 else
431 ! Compute C = alpha * A * B + beta * C
432 if (nrowb < k .or. ncolb /= n) flag = 5
433 end if
434 end if
435 else
436 if (k > n) then
437 flag = 4
438 else
439 if (trans) then
440 ! Compute C = alpha * B**T * A + beta * C
441 if (ncolb /= m .or. nrowb < k) flag = 5
442 else
443 ! Compute C = alpha * B * A + beta * C
444 if (nrowb /= m .or. ncolb < k) flag = 5
445 end if
446 end if
447 end if
448 if (flag /= 0) then
449 ! ERROR: One of the input arrays is not sized correctly
450 allocate(character(len = 256) :: errmsg)
451 write(errmsg, 100) "Input number ", flag, &
452 " is not sized correctly."
453 call errmgr%report_error("diag_mtx_mult_mtx", trim(errmsg), &
454 la_array_size_error)
455 return
456 end if
457
458 ! Deal with ALPHA == 0
459 if (alpha == 0) then
460 if (beta == zero) then
461 c = zero
462 else if (beta /= one) then
463 c = beta * c
464 end if
465 return
466 end if
467
468 ! Process
469 if (lside) then
470 if (trans) then
471 ! Compute C = alpha * A * B**T + beta * C
472 do i = 1, k
473 if (beta == zero) then
474 c(i,:) = zero
475 else if (beta /= one) then
476 c(i,:) = beta * c(i,:)
477 end if
478 temp = alpha * a(i)
479 c(i,:) = c(i,:) + temp * b(:,i)
480 end do
481 else
482 ! Compute C = alpha * A * B + beta * C
483 do i = 1, k
484 if (beta == zero) then
485 c(i,:) = zero
486 else if (beta /= one) then
487 c(i,:) = beta * c(i,:)
488 end if
489 temp = alpha * a(i)
490 c(i,:) = c(i,:) + temp * b(i,:)
491 end do
492 end if
493
494 ! Handle extra rows
495 if (m > k) then
496 if (beta == zero) then
497 c(k+1:m,:) = zero
498 else
499 c(k+1:m,:) = beta * c(k+1:m,:)
500 end if
501 end if
502 else
503 if (trans) then
504 ! Compute C = alpha * B**T * A + beta * C
505 do i = 1, k
506 if (beta == zero) then
507 c(:,i) = zero
508 else if (beta /= one) then
509 c(:,i) = beta * c(:,i)
510 end if
511 temp = alpha * a(i)
512 c(:,i) = c(:,i) + temp * b(i,:)
513 end do
514 else
515 ! Compute C = alpha * B * A + beta * C
516 do i = 1, k
517 if (beta == zero) then
518 c(:,i) = zero
519 else if (beta /= one) then
520 c(:,i) = beta * c(:,i)
521 end if
522 temp = alpha * a(i)
523 c(:,i) = c(:,i) + temp * b(:,i)
524 end do
525 end if
526
527 ! Handle extra columns
528 if (n > k) then
529 if (beta == zero) then
530 c(:,k+1:m) = zero
531 else if (beta /= one) then
532 c(:,k+1:m) = beta * c(:,k+1:m)
533 end if
534 end if
535 end if
536
537 ! Formatting
538100 format(a, i0, a)
539 end subroutine
540
541! ------------------------------------------------------------------------------
542 module subroutine diag_mtx_mult_mtx2(lside, alpha, a, b, err)
543 ! Arguments
544 logical, intent(in) :: lside
545 real(real64), intent(in) :: alpha
546 real(real64), intent(in), dimension(:) :: a
547 real(real64), intent(inout), dimension(:,:) :: b
548 class(errors), intent(inout), optional, target :: err
549
550 ! Parameters
551 real(real64), parameter :: zero = 0.0d0
552 real(real64), parameter :: one = 1.0d0
553
554 ! Local Variables
555 integer(int32) :: i, m, n, k
556 real(real64) :: temp
557 class(errors), pointer :: errmgr
558 type(errors), target :: deferr
559
560 ! Initialization
561 m = size(b, 1)
562 n = size(b, 2)
563 k = size(a)
564 if (present(err)) then
565 errmgr => err
566 else
567 errmgr => deferr
568 end if
569
570 ! Input Check
571 if ((lside .and. k > m) .or. (.not.lside .and. k > n)) then
572 ! ERROR: One of the input arrays is not sized correctly
573 call errmgr%report_error("diag_mtx_mult_mtx2", &
574 "Input number 3 is not sized correctly.", &
575 la_array_size_error)
576 return
577 end if
578
579 ! Process
580 if (lside) then
581 ! Compute B = alpha * A * B
582 do i = 1, k
583 temp = alpha * a(i)
584 b(i,:) = temp * b(i,:)
585 end do
586 if (m > k) b(k+1:m,:) = zero
587 else
588 ! Compute B = alpha * B * A
589 do i = 1, k
590 temp = alpha * a(i)
591 b(:,i) = temp * b(:,i)
592 end do
593 if (n > k) b(:,k+1:n) = zero
594 end if
595 end subroutine
596
597! ------------------------------------------------------------------------------
598 module subroutine diag_mtx_mult_mtx3(lside, trans, alpha, a, b, beta, c, err)
599 ! Arguments
600 logical, intent(in) :: lside, trans
601 real(real64) :: alpha, beta
602 complex(real64), intent(in), dimension(:) :: a
603 real(real64), intent(in), dimension(:,:) :: b
604 complex(real64), intent(inout), dimension(:,:) :: c
605 class(errors), intent(inout), optional, target :: err
606
607 ! Parameters
608 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
609 complex(real64), parameter :: one = (1.0d0, 0.0d0)
610
611 ! Local Variables
612 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
613 complex(real64) :: temp
614 class(errors), pointer :: errmgr
615 type(errors), target :: deferr
616 character(len = :), allocatable :: errmsg
617
618 ! Initialization
619 m = size(c, 1)
620 n = size(c, 2)
621 k = size(a)
622 nrowb = size(b, 1)
623 ncolb = size(b, 2)
624 if (present(err)) then
625 errmgr => err
626 else
627 errmgr => deferr
628 end if
629
630 ! Input Check
631 flag = 0
632 if (lside) then
633 if (k > m) then
634 flag = 4
635 else
636 if (trans) then
637 ! Compute C = alpha * A * B**T + beta * C
638 if (nrowb /= n .or. ncolb < k) flag = 5
639 else
640 ! Compute C = alpha * A * B + beta * C
641 if (nrowb < k .or. ncolb /= n) flag = 5
642 end if
643 end if
644 else
645 if (k > n) then
646 flag = 4
647 else
648 if (trans) then
649 ! Compute C = alpha * B**T * A + beta * C
650 if (ncolb /= m .or. nrowb < k) flag = 5
651 else
652 ! Compute C = alpha * B * A + beta * C
653 if (nrowb /= m .or. ncolb < k) flag = 5
654 end if
655 end if
656 end if
657 if (flag /= 0) then
658 ! ERROR: One of the input arrays is not sized correctly
659 allocate(character(len = 256) :: errmsg)
660 write(errmsg, 100) "Input number ", flag, &
661 " is not sized correctly."
662 call errmgr%report_error("diag_mtx_mult_mtx3", trim(errmsg), &
663 la_array_size_error)
664 return
665 end if
666
667 ! Deal with ALPHA == 0
668 if (alpha == 0) then
669 if (beta == zero) then
670 c = zero
671 else if (beta /= one) then
672 c = beta * c
673 end if
674 return
675 end if
676
677 ! Process
678 if (lside) then
679 if (trans) then
680 ! Compute C = alpha * A * B**T + beta * C
681 do i = 1, k
682 if (beta == zero) then
683 c(i,:) = zero
684 else if (beta /= one) then
685 c(i,:) = beta * c(i,:)
686 end if
687 temp = alpha * a(i)
688 c(i,:) = c(i,:) + temp * b(:,i)
689 end do
690 else
691 ! Compute C = alpha * A * B + beta * C
692 do i = 1, k
693 if (beta == zero) then
694 c(i,:) = zero
695 else if (beta /= one) then
696 c(i,:) = beta * c(i,:)
697 end if
698 temp = alpha * a(i)
699 c(i,:) = c(i,:) + temp * b(i,:)
700 end do
701 end if
702
703 ! Handle extra rows
704 if (m > k) then
705 if (beta == zero) then
706 c(k+1:m,:) = zero
707 else
708 c(k+1:m,:) = beta * c(k+1:m,:)
709 end if
710 end if
711 else
712 if (trans) then
713 ! Compute C = alpha * B**T * A + beta * C
714 do i = 1, k
715 if (beta == zero) then
716 c(:,i) = zero
717 else if (beta /= one) then
718 c(:,i) = beta * c(:,i)
719 end if
720 temp = alpha * a(i)
721 c(:,i) = c(:,i) + temp * b(i,:)
722 end do
723 else
724 ! Compute C = alpha * B * A + beta * C
725 do i = 1, k
726 if (beta == zero) then
727 c(:,i) = zero
728 else if (beta /= one) then
729 c(:,i) = beta * c(:,i)
730 end if
731 temp = alpha * a(i)
732 c(:,i) = c(:,i) + temp * b(:,i)
733 end do
734 end if
735
736 ! Handle extra columns
737 if (n > k) then
738 if (beta == zero) then
739 c(:,k+1:m) = zero
740 else if (beta /= one) then
741 c(:,k+1:m) = beta * c(:,k+1:m)
742 end if
743 end if
744 end if
745
746 ! Formatting
747100 format(a, i0, a)
748 end subroutine
749
750! ------------------------------------------------------------------------------
751 module subroutine diag_mtx_mult_mtx4(lside, opb, alpha, a, b, beta, c, err)
752 ! Arguments
753 logical, intent(in) :: lside
754 integer(int32), intent(in) :: opb
755 real(real64) :: alpha, beta
756 complex(real64), intent(in), dimension(:) :: a
757 complex(real64), intent(in), dimension(:,:) :: b
758 complex(real64), intent(inout), dimension(:,:) :: c
759 class(errors), intent(inout), optional, target :: err
760
761 ! Parameters
762 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
763 complex(real64), parameter :: one = (1.0d0, 0.0d0)
764
765 ! Local Variables
766 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
767 complex(real64) :: temp
768 class(errors), pointer :: errmgr
769 type(errors), target :: deferr
770 character(len = :), allocatable :: errmsg
771
772 ! Initialization
773 m = size(c, 1)
774 n = size(c, 2)
775 k = size(a)
776 nrowb = size(b, 1)
777 ncolb = size(b, 2)
778 if (present(err)) then
779 errmgr => err
780 else
781 errmgr => deferr
782 end if
783
784 ! Input Check
785 flag = 0
786 if (lside) then
787 if (k > m) then
788 flag = 4
789 else
790 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
791 ! Compute C = alpha * A * B**T + beta * C
792 if (nrowb /= n .or. ncolb < k) flag = 5
793 else
794 ! Compute C = alpha * A * B + beta * C
795 if (nrowb < k .or. ncolb /= n) flag = 5
796 end if
797 end if
798 else
799 if (k > n) then
800 flag = 4
801 else
802 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
803 ! Compute C = alpha * B**T * A + beta * C
804 if (ncolb /= m .or. nrowb < k) flag = 5
805 else
806 ! Compute C = alpha * B * A + beta * C
807 if (nrowb /= m .or. ncolb < k) flag = 5
808 end if
809 end if
810 end if
811 if (flag /= 0) then
812 ! ERROR: One of the input arrays is not sized correctly
813 allocate(character(len = 256) :: errmsg)
814 write(errmsg, 100) "Input number ", flag, &
815 " is not sized correctly."
816 call errmgr%report_error("diag_mtx_mult_mtx4", trim(errmsg), &
817 la_array_size_error)
818 return
819 end if
820
821 ! Deal with ALPHA == 0
822 if (alpha == 0) then
823 if (beta == zero) then
824 c = zero
825 else if (beta /= one) then
826 c = beta * c
827 end if
828 return
829 end if
830
831 ! Process
832 if (lside) then
833 if (opb == la_transpose) then
834 ! Compute C = alpha * A * B**T + beta * C
835 do i = 1, k
836 if (beta == zero) then
837 c(i,:) = zero
838 else if (beta /= one) then
839 c(i,:) = beta * c(i,:)
840 end if
841 temp = alpha * a(i)
842 c(i,:) = c(i,:) + temp * b(:,i)
843 end do
844 else if (opb == la_hermitian_transpose) then
845 ! Compute C = alpha * A * B**H + beta * C
846 do i = 1, k
847 if (beta == zero) then
848 c(i,:) = zero
849 else if (beta /= one) then
850 c(i,:) = beta * c(i,:)
851 end if
852 temp = alpha * a(i)
853 c(i,:) = c(i,:) + temp * conjg(b(:,i))
854 end do
855 else
856 ! Compute C = alpha * A * B + beta * C
857 do i = 1, k
858 if (beta == zero) then
859 c(i,:) = zero
860 else if (beta /= one) then
861 c(i,:) = beta * c(i,:)
862 end if
863 temp = alpha * a(i)
864 c(i,:) = c(i,:) + temp * b(i,:)
865 end do
866 end if
867
868 ! Handle extra rows
869 if (m > k) then
870 if (beta == zero) then
871 c(k+1:m,:) = zero
872 else
873 c(k+1:m,:) = beta * c(k+1:m,:)
874 end if
875 end if
876 else
877 if (opb == la_transpose) then
878 ! Compute C = alpha * B**T * A + beta * C
879 do i = 1, k
880 if (beta == zero) then
881 c(:,i) = zero
882 else if (beta /= one) then
883 c(:,i) = beta * c(:,i)
884 end if
885 temp = alpha * a(i)
886 c(:,i) = c(:,i) + temp * b(i,:)
887 end do
888 else if (opb == la_hermitian_transpose) then
889 ! Compute C = alpha * B**H * A + beta * C
890 do i = 1, k
891 if (beta == zero) then
892 c(:,i) = zero
893 else if (beta /= one) then
894 c(:,i) = beta * c(:,i)
895 end if
896 temp = alpha * a(i)
897 c(:,i) = c(:,i) + temp * conjg(b(i,:))
898 end do
899 else
900 ! Compute C = alpha * B * A + beta * C
901 do i = 1, k
902 if (beta == zero) then
903 c(:,i) = zero
904 else if (beta /= one) then
905 c(:,i) = beta * c(:,i)
906 end if
907 temp = alpha * a(i)
908 c(:,i) = c(:,i) + temp * b(:,i)
909 end do
910 end if
911
912 ! Handle extra columns
913 if (n > k) then
914 if (beta == zero) then
915 c(:,k+1:m) = zero
916 else if (beta /= one) then
917 c(:,k+1:m) = beta * c(:,k+1:m)
918 end if
919 end if
920 end if
921
922 ! Formatting
923100 format(a, i0, a)
924 end subroutine
925
926! ------------------------------------------------------------------------------
927 module subroutine diag_mtx_mult_mtx_cmplx(lside, opb, alpha, a, b, beta, c, err)
928 ! Arguments
929 logical, intent(in) :: lside
930 integer(int32), intent(in) :: opb
931 complex(real64) :: alpha, beta
932 complex(real64), intent(in), dimension(:) :: a
933 complex(real64), intent(in), dimension(:,:) :: b
934 complex(real64), intent(inout), dimension(:,:) :: c
935 class(errors), intent(inout), optional, target :: err
936
937 ! Parameters
938 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
939 complex(real64), parameter :: one = (1.0d0, 0.0d0)
940
941 ! Local Variables
942 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
943 complex(real64) :: temp
944 class(errors), pointer :: errmgr
945 type(errors), target :: deferr
946 character(len = :), allocatable :: errmsg
947
948 ! Initialization
949 m = size(c, 1)
950 n = size(c, 2)
951 k = size(a)
952 nrowb = size(b, 1)
953 ncolb = size(b, 2)
954 if (present(err)) then
955 errmgr => err
956 else
957 errmgr => deferr
958 end if
959
960 ! Input Check
961 flag = 0
962 if (lside) then
963 if (k > m) then
964 flag = 4
965 else
966 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
967 ! Compute C = alpha * A * B**T + beta * C
968 if (nrowb /= n .or. ncolb < k) flag = 5
969 else
970 ! Compute C = alpha * A * B + beta * C
971 if (nrowb < k .or. ncolb /= n) flag = 5
972 end if
973 end if
974 else
975 if (k > n) then
976 flag = 4
977 else
978 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
979 ! Compute C = alpha * B**T * A + beta * C
980 if (ncolb /= m .or. nrowb < k) flag = 5
981 else
982 ! Compute C = alpha * B * A + beta * C
983 if (nrowb /= m .or. ncolb < k) flag = 5
984 end if
985 end if
986 end if
987 if (flag /= 0) then
988 ! ERROR: One of the input arrays is not sized correctly
989 allocate(character(len = 256) :: errmsg)
990 write(errmsg, 100) "Input number ", flag, &
991 " is not sized correctly."
992 call errmgr%report_error("diag_mtx_mult_mtx_cmplx", trim(errmsg), &
993 la_array_size_error)
994 return
995 end if
996
997 ! Deal with ALPHA == 0
998 if (alpha == 0) then
999 if (beta == zero) then
1000 c = zero
1001 else if (beta /= one) then
1002 c = beta * c
1003 end if
1004 return
1005 end if
1006
1007 ! Process
1008 if (lside) then
1009 if (opb == la_transpose) then
1010 ! Compute C = alpha * A * B**T + beta * C
1011 do i = 1, k
1012 if (beta == zero) then
1013 c(i,:) = zero
1014 else if (beta /= one) then
1015 c(i,:) = beta * c(i,:)
1016 end if
1017 temp = alpha * a(i)
1018 c(i,:) = c(i,:) + temp * b(:,i)
1019 end do
1020 else if (opb == la_hermitian_transpose) then
1021 ! Compute C = alpha * A * B**H + beta * C
1022 do i = 1, k
1023 if (beta == zero) then
1024 c(i,:) = zero
1025 else if (beta /= one) then
1026 c(i,:) = beta * c(i,:)
1027 end if
1028 temp = alpha * a(i)
1029 c(i,:) = c(i,:) + temp * conjg(b(:,i))
1030 end do
1031 else
1032 ! Compute C = alpha * A * B + beta * C
1033 do i = 1, k
1034 if (beta == zero) then
1035 c(i,:) = zero
1036 else if (beta /= one) then
1037 c(i,:) = beta * c(i,:)
1038 end if
1039 temp = alpha * a(i)
1040 c(i,:) = c(i,:) + temp * b(i,:)
1041 end do
1042 end if
1043
1044 ! Handle extra rows
1045 if (m > k) then
1046 if (beta == zero) then
1047 c(k+1:m,:) = zero
1048 else
1049 c(k+1:m,:) = beta * c(k+1:m,:)
1050 end if
1051 end if
1052 else
1053 if (opb == la_transpose) then
1054 ! Compute C = alpha * B**T * A + beta * C
1055 do i = 1, k
1056 if (beta == zero) then
1057 c(:,i) = zero
1058 else if (beta /= one) then
1059 c(:,i) = beta * c(:,i)
1060 end if
1061 temp = alpha * a(i)
1062 c(:,i) = c(:,i) + temp * b(i,:)
1063 end do
1064 else if (opb == la_hermitian_transpose) then
1065 ! Compute C = alpha * B**H * A + beta * C
1066 do i = 1, k
1067 if (beta == zero) then
1068 c(:,i) = zero
1069 else if (beta /= one) then
1070 c(:,i) = beta * c(:,i)
1071 end if
1072 temp = alpha * a(i)
1073 c(:,i) = c(:,i) + temp * conjg(b(i,:))
1074 end do
1075 else
1076 ! Compute C = alpha * B * A + beta * C
1077 do i = 1, k
1078 if (beta == zero) then
1079 c(:,i) = zero
1080 else if (beta /= one) then
1081 c(:,i) = beta * c(:,i)
1082 end if
1083 temp = alpha * a(i)
1084 c(:,i) = c(:,i) + temp * b(:,i)
1085 end do
1086 end if
1087
1088 ! Handle extra columns
1089 if (n > k) then
1090 if (beta == zero) then
1091 c(:,k+1:m) = zero
1092 else if (beta /= one) then
1093 c(:,k+1:m) = beta * c(:,k+1:m)
1094 end if
1095 end if
1096 end if
1097
1098 ! Formatting
1099100 format(a, i0, a)
1100 end subroutine
1101
1102! ------------------------------------------------------------------------------
1103 module subroutine diag_mtx_mult_mtx2_cmplx(lside, alpha, a, b, err)
1104 ! Arguments
1105 logical, intent(in) :: lside
1106 complex(real64), intent(in) :: alpha
1107 complex(real64), intent(in), dimension(:) :: a
1108 complex(real64), intent(inout), dimension(:,:) :: b
1109 class(errors), intent(inout), optional, target :: err
1110
1111 ! Parameters
1112 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1113 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1114
1115 ! Local Variables
1116 integer(int32) :: i, m, n, k
1117 complex(real64) :: temp
1118 class(errors), pointer :: errmgr
1119 type(errors), target :: deferr
1120
1121 ! Initialization
1122 m = size(b, 1)
1123 n = size(b, 2)
1124 k = size(a)
1125 if (present(err)) then
1126 errmgr => err
1127 else
1128 errmgr => deferr
1129 end if
1130
1131 ! Input Check
1132 if ((lside .and. k > m) .or. (.not.lside .and. k > n)) then
1133 ! ERROR: One of the input arrays is not sized correctly
1134 call errmgr%report_error("diag_mtx_mult_mtx2_cmplx", &
1135 "Input number 3 is not sized correctly.", &
1136 la_array_size_error)
1137 return
1138 end if
1139
1140 ! Process
1141 if (lside) then
1142 ! Compute B = alpha * A * B
1143 do i = 1, k
1144 temp = alpha * a(i)
1145 b(i,:) = temp * b(i,:)
1146 end do
1147 if (m > k) b(k+1:m,:) = zero
1148 else
1149 ! Compute B = alpha * B * A
1150 do i = 1, k
1151 temp = alpha * a(i)
1152 b(:,i) = temp * b(:,i)
1153 end do
1154 if (n > k) b(:,k+1:n) = zero
1155 end if
1156 end subroutine
1157
1158! ------------------------------------------------------------------------------
1159 module subroutine diag_mtx_mult_mtx_mix(lside, opb, alpha, a, b, beta, c, err)
1160 ! Arguments
1161 logical, intent(in) :: lside
1162 integer(int32), intent(in) :: opb
1163 complex(real64) :: alpha, beta
1164 real(real64), intent(in), dimension(:) :: a
1165 complex(real64), intent(in), dimension(:,:) :: b
1166 complex(real64), intent(inout), dimension(:,:) :: c
1167 class(errors), intent(inout), optional, target :: err
1168
1169 ! Parameters
1170 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1171 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1172
1173 ! Local Variables
1174 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
1175 complex(real64) :: temp
1176 class(errors), pointer :: errmgr
1177 type(errors), target :: deferr
1178 character(len = :), allocatable :: errmsg
1179
1180 ! Initialization
1181 m = size(c, 1)
1182 n = size(c, 2)
1183 k = size(a)
1184 nrowb = size(b, 1)
1185 ncolb = size(b, 2)
1186 if (present(err)) then
1187 errmgr => err
1188 else
1189 errmgr => deferr
1190 end if
1191
1192 ! Input Check
1193 flag = 0
1194 if (lside) then
1195 if (k > m) then
1196 flag = 4
1197 else
1198 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
1199 ! Compute C = alpha * A * B**T + beta * C
1200 if (nrowb /= n .or. ncolb < k) flag = 5
1201 else
1202 ! Compute C = alpha * A * B + beta * C
1203 if (nrowb < k .or. ncolb /= n) flag = 5
1204 end if
1205 end if
1206 else
1207 if (k > n) then
1208 flag = 4
1209 else
1210 if (opb == la_transpose .or. opb == la_hermitian_transpose) then
1211 ! Compute C = alpha * B**T * A + beta * C
1212 if (ncolb /= m .or. nrowb < k) flag = 5
1213 else
1214 ! Compute C = alpha * B * A + beta * C
1215 if (nrowb /= m .or. ncolb < k) flag = 5
1216 end if
1217 end if
1218 end if
1219 if (flag /= 0) then
1220 ! ERROR: One of the input arrays is not sized correctly
1221 allocate(character(len = 256) :: errmsg)
1222 write(errmsg, 100) "Input number ", flag, &
1223 " is not sized correctly."
1224 call errmgr%report_error("diag_mtx_mult_mtx_mix", trim(errmsg), &
1225 la_array_size_error)
1226 return
1227 end if
1228
1229 ! Deal with ALPHA == 0
1230 if (alpha == 0) then
1231 if (beta == zero) then
1232 c = zero
1233 else if (beta /= one) then
1234 c = beta * c
1235 end if
1236 return
1237 end if
1238
1239 ! Process
1240 if (lside) then
1241 if (opb == la_transpose) then
1242 ! Compute C = alpha * A * B**T + beta * C
1243 do i = 1, k
1244 if (beta == zero) then
1245 c(i,:) = zero
1246 else if (beta /= one) then
1247 c(i,:) = beta * c(i,:)
1248 end if
1249 temp = alpha * a(i)
1250 c(i,:) = c(i,:) + temp * b(:,i)
1251 end do
1252 else if (opb == la_hermitian_transpose) then
1253 ! Compute C = alpha * A * B**H + beta * C
1254 do i = 1, k
1255 if (beta == zero) then
1256 c(i,:) = zero
1257 else if (beta /= one) then
1258 c(i,:) = beta * c(i,:)
1259 end if
1260 temp = alpha * a(i)
1261 c(i,:) = c(i,:) + temp * conjg(b(:,i))
1262 end do
1263 else
1264 ! Compute C = alpha * A * B + beta * C
1265 do i = 1, k
1266 if (beta == zero) then
1267 c(i,:) = zero
1268 else if (beta /= one) then
1269 c(i,:) = beta * c(i,:)
1270 end if
1271 temp = alpha * a(i)
1272 c(i,:) = c(i,:) + temp * b(i,:)
1273 end do
1274 end if
1275
1276 ! Handle extra rows
1277 if (m > k) then
1278 if (beta == zero) then
1279 c(k+1:m,:) = zero
1280 else
1281 c(k+1:m,:) = beta * c(k+1:m,:)
1282 end if
1283 end if
1284 else
1285 if (opb == la_transpose) then
1286 ! Compute C = alpha * B**T * A + beta * C
1287 do i = 1, k
1288 if (beta == zero) then
1289 c(:,i) = zero
1290 else if (beta /= one) then
1291 c(:,i) = beta * c(:,i)
1292 end if
1293 temp = alpha * a(i)
1294 c(:,i) = c(:,i) + temp * b(i,:)
1295 end do
1296 else if (opb == la_hermitian_transpose) then
1297 ! Compute C = alpha * B**H * A + beta * C
1298 do i = 1, k
1299 if (beta == zero) then
1300 c(:,i) = zero
1301 else if (beta /= one) then
1302 c(:,i) = beta * c(:,i)
1303 end if
1304 temp = alpha * a(i)
1305 c(:,i) = c(:,i) + temp * conjg(b(i,:))
1306 end do
1307 else
1308 ! Compute C = alpha * B * A + beta * C
1309 do i = 1, k
1310 if (beta == zero) then
1311 c(:,i) = zero
1312 else if (beta /= one) then
1313 c(:,i) = beta * c(:,i)
1314 end if
1315 temp = alpha * a(i)
1316 c(:,i) = c(:,i) + temp * b(:,i)
1317 end do
1318 end if
1319
1320 ! Handle extra columns
1321 if (n > k) then
1322 if (beta == zero) then
1323 c(:,k+1:m) = zero
1324 else if (beta /= one) then
1325 c(:,k+1:m) = beta * c(:,k+1:m)
1326 end if
1327 end if
1328 end if
1329
1330 ! Formatting
1331100 format(a, i0, a)
1332 end subroutine
1333
1334! ------------------------------------------------------------------------------
1335 module subroutine diag_mtx_mult_mtx2_mix(lside, alpha, a, b, err)
1336 ! Arguments
1337 logical, intent(in) :: lside
1338 complex(real64), intent(in) :: alpha
1339 real(real64), intent(in), dimension(:) :: a
1340 complex(real64), intent(inout), dimension(:,:) :: b
1341 class(errors), intent(inout), optional, target :: err
1342
1343 ! Parameters
1344 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1345 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1346
1347 ! Local Variables
1348 integer(int32) :: i, m, n, k
1349 complex(real64) :: temp
1350 class(errors), pointer :: errmgr
1351 type(errors), target :: deferr
1352
1353 ! Initialization
1354 m = size(b, 1)
1355 n = size(b, 2)
1356 k = size(a)
1357 if (present(err)) then
1358 errmgr => err
1359 else
1360 errmgr => deferr
1361 end if
1362
1363 ! Input Check
1364 if ((lside .and. k > m) .or. (.not.lside .and. k > n)) then
1365 ! ERROR: One of the input arrays is not sized correctly
1366 call errmgr%report_error("diag_mtx_mult_mtx2_cmplx", &
1367 "Input number 3 is not sized correctly.", &
1368 la_array_size_error)
1369 return
1370 end if
1371
1372 ! Process
1373 if (lside) then
1374 ! Compute B = alpha * A * B
1375 do i = 1, k
1376 temp = alpha * a(i)
1377 b(i,:) = temp * b(i,:)
1378 end do
1379 if (m > k) b(k+1:m,:) = zero
1380 else
1381 ! Compute B = alpha * B * A
1382 do i = 1, k
1383 temp = alpha * a(i)
1384 b(:,i) = temp * b(:,i)
1385 end do
1386 if (n > k) b(:,k+1:n) = zero
1387 end if
1388 end subroutine
1389
1390! ******************************************************************************
1391! BASIC OPERATION ROUTINES
1392! ------------------------------------------------------------------------------
1393 pure module function trace_dbl(x) result(y)
1394 ! Arguments
1395 real(real64), intent(in), dimension(:,:) :: x
1396 real(real64) :: y
1397
1398 ! Parameters
1399 real(real64), parameter :: zero = 0.0d0
1400
1401 ! Local Variables
1402 integer(int32) :: i, m, n, mn
1403
1404 ! Initialization
1405 y = zero
1406 m = size(x, 1)
1407 n = size(x, 2)
1408 mn = min(m, n)
1409
1410 ! Process
1411 do i = 1, mn
1412 y = y + x(i,i)
1413 end do
1414 end function
1415
1416! ------------------------------------------------------------------------------
1417 pure module function trace_cmplx(x) result(y)
1418 ! Arguments
1419 complex(real64), intent(in), dimension(:,:) :: x
1420 complex(real64) :: y
1421
1422 ! Parameters
1423 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1424
1425 ! Local Variables
1426 integer(int32) :: i, m, n, mn
1427
1428 ! Initialization
1429 y = zero
1430 m = size(x, 1)
1431 n = size(x, 2)
1432 mn = min(m, n)
1433
1434 ! Process
1435 do i = 1, mn
1436 y = y + x(i,i)
1437 end do
1438 end function
1439
1440! ------------------------------------------------------------------------------
1441 module function mtx_rank_dbl(a, tol, work, olwork, err) result(rnk)
1442 ! Arguments
1443 real(real64), intent(inout), dimension(:,:) :: a
1444 real(real64), intent(in), optional :: tol
1445 real(real64), intent(out), target, optional, dimension(:) :: work
1446 integer(int32), intent(out), optional :: olwork
1447 class(errors), intent(inout), optional, target :: err
1448 integer(int32) :: rnk
1449
1450 ! Local Variables
1451 integer(int32) :: i, m, n, mn, istat, lwork, flag
1452 real(real64), pointer, dimension(:) :: wptr, s, w
1453 real(real64), allocatable, target, dimension(:) :: wrk
1454 real(real64) :: t, tref, smlnum
1455 real(real64), dimension(1) :: dummy, temp
1456 class(errors), pointer :: errmgr
1457 type(errors), target :: deferr
1458 character(len = :), allocatable :: errmsg
1459
1460 ! Initialization
1461 m = size(a, 1)
1462 n = size(a, 2)
1463 mn = min(m, n)
1464 smlnum = dlamch('s')
1465 rnk = 0
1466 if (present(err)) then
1467 errmgr => err
1468 else
1469 errmgr => deferr
1470 end if
1471
1472 ! Workspace Query
1473 !call svd(a, a(1:mn,1), olwork = lwork)
1474 call dgesvd('N', 'N', m, n, a, m, dummy, dummy, m, dummy, n, temp, &
1475 -1, flag)
1476 lwork = int(temp(1), int32) + mn
1477 if (present(olwork)) then
1478 olwork = lwork
1479 return
1480 end if
1481
1482 ! Local Memory Allocation
1483 if (present(work)) then
1484 if (size(work) < lwork) then
1485 ! ERROR: WORK not sized correctly
1486 call errmgr%report_error("mtx_rank", &
1487 "Incorrectly sized input array WORK, argument 5.", &
1488 la_array_size_error)
1489 return
1490 end if
1491 wptr => work(1:lwork)
1492 else
1493 allocate(wrk(lwork), stat = istat)
1494 if (istat /= 0) then
1495 ! ERROR: Out of memory
1496 call errmgr%report_error("mtx_rank", &
1497 "Insufficient memory available.", &
1498 la_out_of_memory_error)
1499 return
1500 end if
1501 wptr => wrk
1502 end if
1503 s => wptr(1:mn)
1504 w => wptr(mn+1:lwork)
1505
1506 ! Compute the singular values of A
1507 call dgesvd('N', 'N', m, n, a, m, s, dummy, m, dummy, n, w, &
1508 lwork - mn, flag)
1509 if (flag > 0) then
1510 allocate(character(len = 256) :: errmsg)
1511 write(errmsg, 100) flag, " superdiagonals could not " // &
1512 "converge to zero as part of the QR iteration process."
1513 call errmgr%report_warning("mtx_rank", errmsg, la_convergence_error)
1514 end if
1515
1516 ! Determine the threshold tolerance for the singular values such that
1517 ! singular values less than the threshold result in zero when inverted.
1518 tref = max(m, n) * epsilon(t) * s(1)
1519 if (present(tol)) then
1520 t = tol
1521 else
1522 t = tref
1523 end if
1524 if (t < smlnum) then
1525 ! ! The supplied tolerance is too small, simply fall back to the
1526 ! ! default, but issue a warning to the user
1527 ! t = tref
1528 ! call report_warning("mtx_rank", "The supplied tolerance was " // &
1529 ! "smaller than a value that would result in an overflow " // &
1530 ! "condition, or is negative; therefore, the tolerance has " // &
1531 ! "been reset to its default value.")
1532 end if
1533
1534 ! Count the singular values that are larger than the tolerance value
1535 do i = 1, mn
1536 if (s(i) < t) exit
1537 rnk = rnk + 1
1538 end do
1539
1540 ! Formatting
1541100 format(i0, a)
1542 end function
1543
1544! ------------------------------------------------------------------------------
1545 module function mtx_rank_cmplx(a, tol, work, olwork, rwork, err) result(rnk)
1546 ! Arguments
1547 complex(real64), intent(inout), dimension(:,:) :: a
1548 real(real64), intent(in), optional :: tol
1549 complex(real64), intent(out), target, optional, dimension(:) :: work
1550 integer(int32), intent(out), optional :: olwork
1551 real(real64), intent(out), target, optional, dimension(:) :: rwork
1552 class(errors), intent(inout), optional, target :: err
1553 integer(int32) :: rnk
1554
1555 ! External Function Interfaces
1556 interface
1557 function dlamch(cmach) result(x)
1558 use, intrinsic :: iso_fortran_env, only : real64
1559 character, intent(in) :: cmach
1560 real(real64) :: x
1561 end function
1562 end interface
1563
1564 ! Local Variables
1565 integer(int32) :: i, m, n, mn, istat, lwork, flag, lrwork
1566 real(real64), pointer, dimension(:) :: s, rwptr, rw
1567 real(real64), allocatable, target, dimension(:) :: rwrk
1568 complex(real64), allocatable, target, dimension(:) :: wrk
1569 complex(real64), pointer, dimension(:) :: wptr
1570 real(real64) :: t, tref, smlnum
1571 real(real64), dimension(1) :: dummy
1572 complex(real64), dimension(1) :: cdummy, temp
1573 class(errors), pointer :: errmgr
1574 type(errors), target :: deferr
1575 character(len = :), allocatable :: errmsg
1576
1577 ! Initialization
1578 m = size(a, 1)
1579 n = size(a, 2)
1580 mn = min(m, n)
1581 lrwork = 6 * mn
1582 smlnum = dlamch('s')
1583 rnk = 0
1584 if (present(err)) then
1585 errmgr => err
1586 else
1587 errmgr => deferr
1588 end if
1589
1590 ! Workspace Query
1591 call zgesvd('N', 'N', m, n, a, m, dummy, cdummy, m, cdummy, n, temp, &
1592 -1, dummy, flag)
1593 lwork = int(temp(1), int32)
1594 if (present(olwork)) then
1595 olwork = lwork
1596 return
1597 end if
1598
1599 ! Local Memory Allocation
1600 if (present(work)) then
1601 if (size(work) < lwork) then
1602 ! ERROR: WORK not sized correctly
1603 call errmgr%report_error("mtx_rank_cmplx", &
1604 "Incorrectly sized input array WORK, argument 5.", &
1605 la_array_size_error)
1606 return
1607 end if
1608 wptr => work(1:lwork)
1609 else
1610 allocate(wrk(lwork), stat = istat)
1611 if (istat /= 0) then
1612 ! ERROR: Out of memory
1613 call errmgr%report_error("mtx_rank_cmplx", &
1614 "Insufficient memory available.", &
1615 la_out_of_memory_error)
1616 return
1617 end if
1618 wptr => wrk
1619 end if
1620
1621 if (present(rwork)) then
1622 if (size(rwork) < lrwork) then
1623 ! ERROR: RWORK not sized correctly
1624 call errmgr%report_error("mtx_rank_cmplx", &
1625 "Incorrectly sized input array RWORK.", &
1626 la_array_size_error)
1627 return
1628 end if
1629 rwptr => rwork(1:lrwork)
1630 else
1631 allocate(rwrk(lrwork), stat = istat)
1632 if (istat /= 0) then
1633 end if
1634 rwptr => rwrk
1635 end if
1636 s => rwptr(1:mn)
1637 rw => rwptr(mn+1:lrwork)
1638
1639 ! Compute the singular values of A
1640 call zgesvd('N', 'N', m, n, a, m, s, cdummy, m, cdummy, n, wptr, &
1641 lwork - mn, rw, flag)
1642 if (flag > 0) then
1643 allocate(character(len = 256) :: errmsg)
1644 write(errmsg, 100) flag, " superdiagonals could not " // &
1645 "converge to zero as part of the QR iteration process."
1646 call errmgr%report_warning("mtx_rank_cmplx", errmsg, la_convergence_error)
1647 end if
1648
1649 ! Determine the threshold tolerance for the singular values such that
1650 ! singular values less than the threshold result in zero when inverted.
1651 tref = max(m, n) * epsilon(t) * s(1)
1652 if (present(tol)) then
1653 t = tol
1654 else
1655 t = tref
1656 end if
1657 if (t < smlnum) then
1658 ! ! The supplied tolerance is too small, simply fall back to the
1659 ! ! default, but issue a warning to the user
1660 ! t = tref
1661 ! call report_warning("mtx_rank", "The supplied tolerance was " // &
1662 ! "smaller than a value that would result in an overflow " // &
1663 ! "condition, or is negative; therefore, the tolerance has " // &
1664 ! "been reset to its default value.")
1665 end if
1666
1667 ! Count the singular values that are larger than the tolerance value
1668 do i = 1, mn
1669 if (s(i) < t) exit
1670 rnk = rnk + 1
1671 end do
1672
1673 ! Formatting
1674100 format(i0, a)
1675 end function
1676
1677! ------------------------------------------------------------------------------
1678 module function det_dbl(a, iwork, err) result(x)
1679 ! Arguments
1680 real(real64), intent(inout), dimension(:,:) :: a
1681 integer(int32), intent(out), target, optional, dimension(:) :: iwork
1682 class(errors), intent(inout), optional, target :: err
1683 real(real64) :: x
1684
1685 ! Parameters
1686 real(real64), parameter :: zero = 0.0d0
1687 real(real64), parameter :: one = 1.0d0
1688 real(real64), parameter :: ten = 1.0d1
1689 real(real64), parameter :: p1 = 1.0d-1
1690
1691 ! Local Variables
1692 integer(int32) :: i, ep, n, istat, flag
1693 integer(int32), pointer, dimension(:) :: ipvt
1694 integer(int32), allocatable, target, dimension(:) :: iwrk
1695 real(real64) :: temp
1696 class(errors), pointer :: errmgr
1697 type(errors), target :: deferr
1698
1699 ! Initialization
1700 n = size(a, 1)
1701 x = zero
1702 if (present(err)) then
1703 errmgr => err
1704 else
1705 errmgr => deferr
1706 end if
1707
1708 ! Input Check
1709 if (size(a, 2) /= n) then
1710 call errmgr%report_error("det", &
1711 "The supplied matrix must be square.", la_array_size_error)
1712 return
1713 end if
1714
1715 ! Local Memory Allocation
1716 if (present(iwork)) then
1717 if (size(iwork) < n) then
1718 ! ERROR: WORK not sized correctly
1719 call errmgr%report_error("det", &
1720 "Incorrectly sized input array IWORK, argument 2.", &
1721 la_array_size_error)
1722 return
1723 end if
1724 ipvt => iwork(1:n)
1725 else
1726 allocate(iwrk(n), stat = istat)
1727 if (istat /= 0) then
1728 ! ERROR: Out of memory
1729 call errmgr%report_error("det", &
1730 "Insufficient memory available.", &
1731 la_out_of_memory_error)
1732 return
1733 end if
1734 ipvt => iwrk
1735 end if
1736
1737 ! Compute the LU factorization of A
1738 call dgetrf(n, n, a, n, ipvt, flag)
1739 if (flag > 0) then
1740 ! A singular matrix has a determinant of zero
1741 x = zero
1742 return
1743 end if
1744
1745 ! Compute the product of the diagonal of A
1746 temp = one
1747 ep = 0
1748 do i = 1, n
1749 if (ipvt(i) /= i) temp = -temp
1750
1751 temp = a(i,i) * temp
1752 if (temp == zero) then
1753 x = zero
1754 exit
1755 end if
1756
1757 do while (abs(temp) < one)
1758 temp = ten * temp
1759 ep = ep - 1
1760 end do
1761
1762 do while (abs(temp) > ten)
1763 temp = p1 * temp
1764 ep = ep + 1
1765 end do
1766 end do
1767 x = temp * ten**ep
1768 end function
1769
1770! ------------------------------------------------------------------------------
1771 module function det_cmplx(a, iwork, err) result(x)
1772 ! Arguments
1773 complex(real64), intent(inout), dimension(:,:) :: a
1774 integer(int32), intent(out), target, optional, dimension(:) :: iwork
1775 class(errors), intent(inout), optional, target :: err
1776 complex(real64) :: x
1777
1778 ! Parameters
1779 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
1780 complex(real64), parameter :: one = (1.0d0, 0.0d0)
1781 complex(real64), parameter :: ten = (1.0d1, 0.0d0)
1782 complex(real64), parameter :: p1 = (1.0d-1, 0.0d0)
1783 real(real64), parameter :: real_one = 1.0d0
1784 real(real64), parameter :: real_ten = 1.0d1
1785
1786 ! Local Variables
1787 integer(int32) :: i, ep, n, istat, flag
1788 integer(int32), pointer, dimension(:) :: ipvt
1789 integer(int32), allocatable, target, dimension(:) :: iwrk
1790 complex(real64) :: temp
1791 class(errors), pointer :: errmgr
1792 type(errors), target :: deferr
1793
1794 ! Initialization
1795 n = size(a, 1)
1796 x = zero
1797 if (present(err)) then
1798 errmgr => err
1799 else
1800 errmgr => deferr
1801 end if
1802
1803 ! Input Check
1804 if (size(a, 2) /= n) then
1805 call errmgr%report_error("det_cmplx", &
1806 "The supplied matrix must be square.", la_array_size_error)
1807 return
1808 end if
1809
1810 ! Local Memory Allocation
1811 if (present(iwork)) then
1812 if (size(iwork) < n) then
1813 ! ERROR: WORK not sized correctly
1814 call errmgr%report_error("det_cmplx", &
1815 "Incorrectly sized input array IWORK, argument 2.", &
1816 la_array_size_error)
1817 return
1818 end if
1819 ipvt => iwork(1:n)
1820 else
1821 allocate(iwrk(n), stat = istat)
1822 if (istat /= 0) then
1823 ! ERROR: Out of memory
1824 call errmgr%report_error("det_cmplx", &
1825 "Insufficient memory available.", &
1826 la_out_of_memory_error)
1827 return
1828 end if
1829 ipvt => iwrk
1830 end if
1831
1832 ! Compute the LU factorization of A
1833 call zgetrf(n, n, a, n, ipvt, flag)
1834 if (flag > 0) then
1835 ! A singular matrix has a determinant of zero
1836 x = zero
1837 return
1838 end if
1839
1840 ! Compute the product of the diagonal of A
1841 temp = one
1842 ep = 0
1843 do i = 1, n
1844 if (ipvt(i) /= i) temp = -temp
1845
1846 temp = a(i,i) * temp
1847 if (temp == zero) then
1848 x = zero
1849 exit
1850 end if
1851
1852 do while (abs(temp) < real_one)
1853 temp = ten * temp
1854 ep = ep - 1
1855 end do
1856
1857 do while (abs(temp) > real_ten)
1858 temp = p1 * temp
1859 ep = ep + 1
1860 end do
1861 end do
1862 x = temp * ten**ep
1863 end function
1864
1865! ******************************************************************************
1866! ARRAY SWAPPING ROUTINE
1867! ------------------------------------------------------------------------------
1868 module subroutine swap_dbl(x, y, err)
1869 ! Arguments
1870 real(real64), intent(inout), dimension(:) :: x, y
1871 class(errors), intent(inout), optional, target :: err
1872
1873 ! Local Variables
1874 integer(int32) :: i, n
1875 real(real64) :: temp
1876 class(errors), pointer :: errmgr
1877 type(errors), target :: deferr
1878
1879 ! Initialization
1880 n = size(x)
1881 if (present(err)) then
1882 errmgr => err
1883 else
1884 errmgr => deferr
1885 end if
1886
1887 ! Input Check
1888 if (size(y) /= n) then
1889 call errmgr%report_error("swap", &
1890 "The input arrays are not the same size.", &
1891 la_array_size_error)
1892 return
1893 end if
1894
1895 ! Process
1896 do i = 1, n
1897 temp = x(i)
1898 x(i) = y(i)
1899 y(i) = temp
1900 end do
1901 end subroutine
1902
1903! ------------------------------------------------------------------------------
1904 module subroutine swap_cmplx(x, y, err)
1905 ! Arguments
1906 complex(real64), intent(inout), dimension(:) :: x, y
1907 class(errors), intent(inout), optional, target :: err
1908
1909 ! Local Variables
1910 integer(int32) :: i, n
1911 complex(real64) :: temp
1912 class(errors), pointer :: errmgr
1913 type(errors), target :: deferr
1914
1915 ! Initialization
1916 n = size(x)
1917 if (present(err)) then
1918 errmgr => err
1919 else
1920 errmgr => deferr
1921 end if
1922
1923 ! Input Check
1924 if (size(y) /= n) then
1925 call errmgr%report_error("swap_cmplx", &
1926 "The input arrays are not the same size.", &
1927 la_array_size_error)
1928 return
1929 end if
1930
1931 ! Process
1932 do i = 1, n
1933 temp = x(i)
1934 x(i) = y(i)
1935 y(i) = temp
1936 end do
1937 end subroutine
1938
1939! ******************************************************************************
1940! ARRAY MULTIPLICIATION ROUTINES
1941! ------------------------------------------------------------------------------
1942 module subroutine recip_mult_array_dbl(a, x)
1943 ! Arguments
1944 real(real64), intent(in) :: a
1945 real(real64), intent(inout), dimension(:) :: x
1946
1947 ! External Function Interfaces
1948 interface
1949 function dlamch(cmach) result(x)
1950 use, intrinsic :: iso_fortran_env, only : real64
1951 character, intent(in) :: cmach
1952 real(real64) :: x
1953 end function
1954 end interface
1955
1956 ! Parameters
1957 real(real64), parameter :: zero = 0.0d0
1958 real(real64), parameter :: one = 1.0d0
1959 real(real64), parameter :: twotho = 2.0d3
1960
1961 ! Local Variables
1962 logical :: done
1963 real(real64) :: bignum, cden, cden1, cnum, cnum1, mul, smlnum
1964
1965 ! Initialization
1966 smlnum = dlamch('s')
1967 bignum = one / smlnum
1968 if (log10(bignum) > twotho) then
1969 smlnum = sqrt(smlnum)
1970 bignum = sqrt(bignum)
1971 end if
1972
1973 ! Initialize the denominator to A, and the numerator to ONE
1974 cden = a
1975 cnum = one
1976
1977 ! Process
1978 do
1979 cden1 = cden * smlnum
1980 cnum1 = cnum / bignum
1981 if (abs(cden1) > abs(cnum) .and. cnum /= zero) then
1982 mul = smlnum
1983 done = .false.
1984 cden = cden1
1985 else if (abs(cnum1) > abs(cden)) then
1986 mul = bignum
1987 done = .false.
1988 cnum = cnum1
1989 else
1990 mul = cnum / cden
1991 done = .true.
1992 end if
1993
1994 ! Scale the vector X by MUL
1995 x = mul * x
1996
1997 ! Exit if done
1998 if (done) exit
1999 end do
2000 end subroutine
2001
2002! ******************************************************************************
2003! TRIANGULAR MATRIX MULTIPLICATION ROUTINES
2004! ------------------------------------------------------------------------------
2005 module subroutine tri_mtx_mult_dbl(upper, alpha, a, beta, b, err)
2006 ! Arguments
2007 logical, intent(in) :: upper
2008 real(real64), intent(in) :: alpha, beta
2009 real(real64), intent(in), dimension(:,:) :: a
2010 real(real64), intent(inout), dimension(:,:) :: b
2011 class(errors), intent(inout), optional, target :: err
2012
2013 ! Parameters
2014 real(real64), parameter :: zero = 0.0d0
2015
2016 ! Local Variables
2017 integer(int32) :: i, j, k, n, d1, d2, flag
2018 real(real64) :: temp
2019 class(errors), pointer :: errmgr
2020 type(errors), target :: deferr
2021 character(len = :), allocatable :: errmsg
2022
2023 ! Initialization
2024 n = size(a, 1)
2025 d1 = n
2026 d2 = n
2027 if (present(err)) then
2028 errmgr => err
2029 else
2030 errmgr => deferr
2031 end if
2032
2033 ! Input Check
2034 flag = 0
2035 if (size(a, 2) /= n) then
2036 flag = 3
2037 d2 = size(a, 2)
2038 else if (size(b, 1) /= n .or. size(b, 2) /= n) then
2039 flag = 5
2040 d1 = size(b, 1)
2041 d2 = size(b, 2)
2042 end if
2043 if (flag /= 0) then
2044 ! ERROR: Incorrectly sized matrix
2045 allocate(character(len = 256) :: errmsg)
2046 write(errmsg, 100) "The matrix at input ", flag, &
2047 " was not sized appropriately. A matrix of ", n, "-by-", n, &
2048 "was expected, but a matrix of ", d1, "-by-", d2, " was found."
2049 call errmgr%report_error("tri_mtx_mult_dbl", trim(errmsg), &
2050 la_array_size_error)
2051 return
2052 end if
2053
2054 ! Process
2055 if (upper) then
2056 ! Form: B = alpha * A**T * A + beta * B
2057 if (beta == zero) then
2058 do j = 1, n
2059 do i = 1, j
2060 temp = zero
2061 do k = 1, j
2062 temp = temp + a(k,i) * a(k,j)
2063 end do
2064 temp = alpha * temp
2065 b(i,j) = temp
2066 if (i /= j) b(j,i) = temp
2067 end do
2068 end do
2069 else
2070 do j = 1, n
2071 do i = 1, j
2072 temp = zero
2073 do k = 1, j
2074 temp = temp + a(k,i) * a(k,j)
2075 end do
2076 temp = alpha * temp
2077 b(i,j) = temp + beta * b(i,j)
2078 if (i /= j) b(j,i) = temp + beta * b(j,i)
2079 end do
2080 end do
2081 end if
2082 else
2083 ! Form: B = alpha * A * A**T + beta * B
2084 if (beta == zero) then
2085 do j = 1, n
2086 do i = j, n
2087 temp = zero
2088 do k = 1, j
2089 temp = temp + a(i,k) * a(j,k)
2090 end do
2091 temp = alpha * temp
2092 b(i,j) = temp
2093 if (i /= j) b(j,i) = temp
2094 end do
2095 end do
2096 else
2097 do j = 1, n
2098 do i = j, n
2099 temp = zero
2100 do k = 1, j
2101 temp = temp + a(i,k) * a(j,k)
2102 end do
2103 temp = alpha * temp
2104 b(i,j) = temp + beta * b(i,j)
2105 if (i /= j) b(j,i) = temp + beta * b(j,i)
2106 end do
2107 end do
2108 end if
2109 end if
2110
2111 ! Formatting
2112100 format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2113 end subroutine
2114
2115! ------------------------------------------------------------------------------
2116 module subroutine tri_mtx_mult_cmplx(upper, alpha, a, beta, b, err)
2117 ! Arguments
2118 logical, intent(in) :: upper
2119 complex(real64), intent(in) :: alpha, beta
2120 complex(real64), intent(in), dimension(:,:) :: a
2121 complex(real64), intent(inout), dimension(:,:) :: b
2122 class(errors), intent(inout), optional, target :: err
2123
2124 ! Parameters
2125 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
2126
2127 ! Local Variables
2128 integer(int32) :: i, j, k, n, d1, d2, flag
2129 complex(real64) :: temp
2130 class(errors), pointer :: errmgr
2131 type(errors), target :: deferr
2132 character(len = :), allocatable :: errmsg
2133
2134 ! Initialization
2135 n = size(a, 1)
2136 d1 = n
2137 d2 = n
2138 if (present(err)) then
2139 errmgr => err
2140 else
2141 errmgr => deferr
2142 end if
2143
2144 ! Input Check
2145 flag = 0
2146 if (size(a, 2) /= n) then
2147 flag = 3
2148 d2 = size(a, 2)
2149 else if (size(b, 1) /= n .or. size(b, 2) /= n) then
2150 flag = 5
2151 d1 = size(b, 1)
2152 d2 = size(b, 2)
2153 end if
2154 if (flag /= 0) then
2155 ! ERROR: Incorrectly sized matrix
2156 allocate(character(len = 256) :: errmsg)
2157 write(errmsg, 100) "The matrix at input ", flag, &
2158 " was not sized appropriately. A matrix of ", n, "-by-", n, &
2159 "was expected, but a matrix of ", d1, "-by-", d2, " was found."
2160 call errmgr%report_error("tri_mtx_mult_cmplx", trim(errmsg), &
2161 la_array_size_error)
2162 return
2163 end if
2164
2165 ! Process
2166 if (upper) then
2167 ! Form: B = alpha * A**T * A + beta * B
2168 if (beta == zero) then
2169 do j = 1, n
2170 do i = 1, j
2171 temp = zero
2172 do k = 1, j
2173 temp = temp + a(k,i) * a(k,j)
2174 end do
2175 temp = alpha * temp
2176 b(i,j) = temp
2177 if (i /= j) b(j,i) = temp
2178 end do
2179 end do
2180 else
2181 do j = 1, n
2182 do i = 1, j
2183 temp = zero
2184 do k = 1, j
2185 temp = temp + a(k,i) * a(k,j)
2186 end do
2187 temp = alpha * temp
2188 b(i,j) = temp + beta * b(i,j)
2189 if (i /= j) b(j,i) = temp + beta * b(j,i)
2190 end do
2191 end do
2192 end if
2193 else
2194 ! Form: B = alpha * A * A**T + beta * B
2195 if (beta == zero) then
2196 do j = 1, n
2197 do i = j, n
2198 temp = zero
2199 do k = 1, j
2200 temp = temp + a(i,k) * a(j,k)
2201 end do
2202 temp = alpha * temp
2203 b(i,j) = temp
2204 if (i /= j) b(j,i) = temp
2205 end do
2206 end do
2207 else
2208 do j = 1, n
2209 do i = j, n
2210 temp = zero
2211 do k = 1, j
2212 temp = temp + a(i,k) * a(j,k)
2213 end do
2214 temp = alpha * temp
2215 b(i,j) = temp + beta * b(i,j)
2216 if (i /= j) b(j,i) = temp + beta * b(j,i)
2217 end do
2218 end do
2219 end if
2220 end if
2221
2222 ! Formatting
2223100 format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2224 end subroutine
2225
2226! ******************************************************************************
2227! BANDED MATRIX MULTIPLICATION ROUTINES
2228! ------------------------------------------------------------------------------
2229 module subroutine band_mtx_vec_mult_dbl(trans, kl, ku, alpha, a, x, beta, &
2230 y, err)
2231 ! Arguments
2232 logical, intent(in) :: trans
2233 integer(int32), intent(in) :: kl, ku
2234 real(real64), intent(in) :: alpha, beta
2235 real(real64), intent(in), dimension(:,:) :: a
2236 real(real64), intent(in), dimension(:) :: x
2237 real(real64), intent(inout), dimension(:) :: y
2238 class(errors), intent(inout), optional, target :: err
2239
2240 ! Local Variables
2241 integer(int32) :: m, n
2242 class(errors), pointer :: errmgr
2243 type(errors), target :: deferr
2244
2245 ! Initialization
2246 if (present(err)) then
2247 errmgr => err
2248 else
2249 errmgr => deferr
2250 end if
2251 if (trans) then
2252 m = size(x)
2253 n = size(y)
2254 else
2255 m = size(y)
2256 n = size(x)
2257 end if
2258
2259 ! Input Checking
2260 if (kl < 0) go to 10
2261 if (ku < 0) go to 20
2262 if (size(a, 1) /= kl + ku + 1) go to 30
2263 if (size(a, 2) /= n) go to 30
2264
2265 ! Process
2266 if (trans) then
2267 call dgbmv("T", m, n, kl, ku, alpha, a, size(a, 1), x, 1, beta, y, 1)
2268 else
2269 call dgbmv("N", m, n, kl, ku, alpha, a, size(a, 1), x, 1, beta, y, 1)
2270 end if
2271
2272 ! End
2273 return
2274
2275 ! KL < 0
227610 continue
2277 call errmgr%report_error("band_mtx_vec_mult_dbl", &
2278 "The number of subdiagonals must be at least 0.", &
2279 la_invalid_input_error)
2280 return
2281
2282 ! KU < 0
228320 continue
2284 call errmgr%report_error("band_mtx_vec_mult_dbl", &
2285 "The number of superdiagonals must be at least 0.", &
2286 la_invalid_input_error)
2287 return
2288
2289 ! A is incorrectly sized
229030 continue
2291 call errmgr%report_error("band_mtx_vec_mult_dbl", &
2292 "The size of matrix A is not compatible with the other vectors.", &
2293 la_array_size_error)
2294 return
2295 end subroutine
2296
2297! ------------------------------------------------------------------------------
2298 module subroutine band_mtx_vec_mult_cmplx(trans, kl, ku, alpha, a, x, &
2299 beta, y, err)
2300 ! Arguments
2301 integer(int32), intent(in) :: trans
2302 integer(int32), intent(in) :: kl, ku
2303 complex(real64), intent(in) :: alpha, beta
2304 complex(real64), intent(in), dimension(:,:) :: a
2305 complex(real64), intent(in), dimension(:) :: x
2306 complex(real64), intent(inout), dimension(:) :: y
2307 class(errors), intent(inout), optional, target :: err
2308
2309 ! Local Variables
2310 character :: op
2311 logical :: trns
2312 integer(int32) :: m, n
2313 class(errors), pointer :: errmgr
2314 type(errors), target :: deferr
2315
2316 ! Initialization
2317 if (present(err)) then
2318 errmgr => err
2319 else
2320 errmgr => deferr
2321 end if
2322 if (trans == la_transpose) then
2323 op = "T"
2324 trns = .true.
2325 else if (trans == la_hermitian_transpose) then
2326 op = "C"
2327 trns = .true.
2328 else
2329 op = "N"
2330 trns = .false.
2331 end if
2332 if (trns) then
2333 m = size(x)
2334 n = size(y)
2335 else
2336 m = size(y)
2337 n = size(x)
2338 end if
2339
2340 ! Input Checking
2341 if (kl < 0) go to 10
2342 if (ku < 0) go to 20
2343 if (size(a, 1) /= kl + ku + 1) go to 30
2344 if (size(a, 2) /= n) go to 30
2345
2346 ! Process
2347 call zgbmv(op, m, n, kl, ku, alpha, a, size(a, 1), x, 1, beta, y, 1)
2348
2349 ! End
2350 return
2351
2352 ! KL < 0
235310 continue
2354 call errmgr%report_error("band_mtx_vec_mult_cmplx", &
2355 "The number of subdiagonals must be at least 0.", &
2356 la_invalid_input_error)
2357 return
2358
2359 ! KU < 0
236020 continue
2361 call errmgr%report_error("band_mtx_vec_mult_cmplx", &
2362 "The number of superdiagonals must be at least 0.", &
2363 la_invalid_input_error)
2364 return
2365
2366 ! A is incorrectly sized
236730 continue
2368 call errmgr%report_error("band_mtx_vec_mult_cmplx", &
2369 "The size of matrix A is not compatible with the other vectors.", &
2370 la_array_size_error)
2371 return
2372 end subroutine
2373
2374! ------------------------------------------------------------------------------
2375 module subroutine band_to_full_mtx_dbl(kl, ku, b, f, err)
2376 ! Arguments
2377 integer(int32), intent(in) :: kl, ku
2378 real(real64), intent(in), dimension(:,:) :: b
2379 real(real64), intent(out), dimension(:,:) :: f
2380 class(errors), intent(inout), optional, target :: err
2381
2382 ! Parameters
2383 real(real64), parameter :: zero = 0.0d0
2384
2385 ! Local Variables
2386 class(errors), pointer :: errmgr
2387 type(errors), target :: deferr
2388 integer(int32) :: i, j, k, m, n, i1, i2
2389
2390 ! Initialization
2391 if (present(err)) then
2392 errmgr => err
2393 else
2394 errmgr => deferr
2395 end if
2396 m = size(f, 1)
2397 n = size(f, 2)
2398
2399 ! Input Check
2400 if (kl < 0) go to 10
2401 if (ku < 0) go to 20
2402 if (size(b, 2) /= n) go to 30
2403 if (size(b, 1) /= kl + ku + 1) go to 40
2404
2405 ! Process
2406 do j = 1, n
2407 k = ku + 1 - j
2408 i1 = max(1, j - ku)
2409 i2 = min(m, j + kl)
2410 do i = 1, i1 - 1
2411 f(i,j) = zero
2412 end do
2413 do i = i1, i2
2414 f(i,j) = b(k+i,j)
2415 end do
2416 do i = i2 + 1, m
2417 f(i,j) = zero
2418 end do
2419 end do
2420
2421 ! End
2422 return
2423
2424 ! KL < 0
242510 continue
2426 call errmgr%report_error("band_to_full_mtx_dbl", &
2427 "The number of subdiagonals must be at least 0.", &
2428 la_invalid_input_error)
2429 return
2430
2431 ! KU < 0
243220 continue
2433 call errmgr%report_error("band_to_full_mtx_dbl", &
2434 "The number of superdiagonals must be at least 0.", &
2435 la_invalid_input_error)
2436 return
2437
2438 ! A is incorrectly sized
243930 continue
2440 call errmgr%report_error("band_to_full_mtx_dbl", &
2441 "The number of columns in the banded matrix does not match " // &
2442 "the number of columns in the full matrix.", &
2443 la_array_size_error)
2444 return
2445
244640 continue
2447 call errmgr%report_error("band_to_full_mtx_dbl", &
2448 "The number of rows in the banded matrix does not align with " // &
2449 "the number of sub and super-diagonals specified.", &
2450 la_array_size_error)
2451 return
2452 end subroutine
2453
2454! ------------------------------------------------------------------------------
2455 module subroutine band_to_full_mtx_cmplx(kl, ku, b, f, err)
2456 ! Arguments
2457 integer(int32), intent(in) :: kl, ku
2458 complex(real64), intent(in), dimension(:,:) :: b
2459 complex(real64), intent(out), dimension(:,:) :: f
2460 class(errors), intent(inout), optional, target :: err
2461
2462 ! Parameters
2463 complex(real64), parameter :: zero = (0.0d0, 0.0d0)
2464
2465 ! Local Variables
2466 class(errors), pointer :: errmgr
2467 type(errors), target :: deferr
2468 integer(int32) :: i, j, k, m, n, i1, i2
2469
2470 ! Initialization
2471 if (present(err)) then
2472 errmgr => err
2473 else
2474 errmgr => deferr
2475 end if
2476 m = size(f, 1)
2477 n = size(f, 2)
2478
2479 ! Input Check
2480 if (kl < 0) go to 10
2481 if (ku < 0) go to 20
2482 if (size(b, 2) /= n) go to 30
2483 if (size(b, 1) /= kl + ku + 1) go to 40
2484
2485 ! Process
2486 do j = 1, n
2487 k = ku + 1 - j
2488 i1 = max(1, j - ku)
2489 i2 = min(m, j + kl)
2490 do i = 1, i1 - 1
2491 f(i,j) = zero
2492 end do
2493 do i = i1, i2
2494 f(i,j) = b(k+i,j)
2495 end do
2496 do i = i2 + 1, m
2497 f(i,j) = zero
2498 end do
2499 end do
2500
2501 ! End
2502 return
2503
2504 ! KL < 0
250510 continue
2506 call errmgr%report_error("band_to_full_mtx_cmplx", &
2507 "The number of subdiagonals must be at least 0.", &
2508 la_invalid_input_error)
2509 return
2510
2511 ! KU < 0
251220 continue
2513 call errmgr%report_error("band_to_full_mtx_cmplx", &
2514 "The number of superdiagonals must be at least 0.", &
2515 la_invalid_input_error)
2516 return
2517
2518 ! A is incorrectly sized
251930 continue
2520 call errmgr%report_error("band_to_full_mtx_cmplx", &
2521 "The number of columns in the banded matrix does not match " // &
2522 "the number of columns in the full matrix.", &
2523 la_array_size_error)
2524 return
2525
252640 continue
2527 call errmgr%report_error("band_to_full_mtx_cmplx", &
2528 "The number of rows in the banded matrix does not align with " // &
2529 "the number of sub and super-diagonals specified.", &
2530 la_array_size_error)
2531 return
2532 end subroutine
2533
2534! ------------------------------------------------------------------------------
2535 module subroutine band_diag_mtx_mult_dbl(left, m, kl, ku, alpha, a, b, err)
2536 ! Arguments
2537 logical, intent(in) :: left
2538 integer(int32), intent(in) :: m, kl, ku
2539 real(real64), intent(in) :: alpha
2540 real(real64), intent(inout), dimension(:,:) :: a
2541 real(real64), intent(in), dimension(:) :: b
2542 class(errors), intent(inout), optional, target :: err
2543
2544 ! Parameters
2545 real(real64), parameter :: one = 1.0d0
2546
2547 ! Local Variables
2548 class(errors), pointer :: errmgr
2549 type(errors), target :: deferr
2550 integer(int32) :: i, i1, i2, j, k, n
2551 real(real64) :: temp
2552
2553 ! Initialization
2554 if (present(err)) then
2555 errmgr => err
2556 else
2557 errmgr => deferr
2558 end if
2559 n = size(a, 2)
2560
2561 ! Input Checking
2562 if (kl < 0) go to 10
2563 if (ku < 0) go to 20
2564 if (left) then
2565 if (size(b) /= n) go to 30
2566 else
2567 if (size(b) < m) go to 30
2568 end if
2569
2570 ! Process
2571 if (left) then
2572 ! Compute A = A * B
2573 do j = 1, n
2574 k = ku + 1 - j
2575 i1 = max(1, j - ku) + k
2576 i2 = min(m, j + kl) + k
2577 if (alpha == one) then
2578 temp = b(j)
2579 else
2580 temp = alpha * b(j)
2581 end if
2582 do i = i1, i2
2583 a(i,j) = a(i,j) * temp
2584 end do
2585 end do
2586 else
2587 ! Compute A = B * A
2588 do j = 1, n
2589 k = ku + 1 - j
2590 i1 = max(1, j - ku)
2591 i2 = min(m, j + kl)
2592 if (alpha == 1.0d0) then
2593 do i = i1, i2
2594 a(i+k,j) = a(i+k,j) * b(i)
2595 end do
2596 else
2597 do i = i1, i2
2598 a(i+k,j) = alpha * a(i+k,j) * b(i)
2599 end do
2600 end if
2601 end do
2602 end if
2603
2604
2605 ! End
2606 return
2607
2608 ! KL < 0
260910 continue
2610 call errmgr%report_error("band_diag_mtx_mult_dbl", &
2611 "The number of subdiagonals must be at least 0.", &
2612 la_invalid_input_error)
2613 return
2614
2615 ! KU < 0
261620 continue
2617 call errmgr%report_error("band_diag_mtx_mult_dbl", &
2618 "The number of superdiagonals must be at least 0.", &
2619 la_invalid_input_error)
2620 return
2621
2622 ! B is not sized correctly
262330 continue
2624 call errmgr%report_error("band_diag_mtx_mult_dbl", &
2625 "Inner matrix dimensions do not agree.", &
2626 la_array_size_error)
2627 return
2628 end subroutine
2629
2630! ------------------------------------------------------------------------------
2631 module subroutine band_diag_mtx_mult_cmplx(left, m, kl, ku, alpha, a, b, err)
2632 ! Arguments
2633 logical, intent(in) :: left
2634 integer(int32), intent(in) :: m, kl, ku
2635 complex(real64), intent(in) :: alpha
2636 complex(real64), intent(inout), dimension(:,:) :: a
2637 complex(real64), intent(in), dimension(:) :: b
2638 class(errors), intent(inout), optional, target :: err
2639
2640 ! Parameters
2641 complex(real64), parameter :: one = (1.0d0, 0.0d0)
2642
2643 ! Local Variables
2644 class(errors), pointer :: errmgr
2645 type(errors), target :: deferr
2646 integer(int32) :: i, i1, i2, j, k, n
2647 complex(real64) :: temp
2648
2649 ! Initialization
2650 if (present(err)) then
2651 errmgr => err
2652 else
2653 errmgr => deferr
2654 end if
2655 n = size(a, 2)
2656
2657 ! Input Checking
2658 if (kl < 0) go to 10
2659 if (ku < 0) go to 20
2660 if (left) then
2661 if (size(b) /= n) go to 30
2662 else
2663 if (size(b) < m) go to 30
2664 end if
2665
2666 ! Process
2667 if (left) then
2668 ! Compute A = A * B
2669 do j = 1, n
2670 k = ku + 1 - j
2671 i1 = max(1, j - ku) + k
2672 i2 = min(m, j + kl) + k
2673 if (alpha == one) then
2674 temp = b(j)
2675 else
2676 temp = alpha * b(j)
2677 end if
2678 do i = i1, i2
2679 a(i,j) = a(i,j) * temp
2680 end do
2681 end do
2682 else
2683 ! Compute A = B * A
2684 do j = 1, n
2685 k = ku + 1 - j
2686 i1 = max(1, j - ku)
2687 i2 = min(m, j + kl)
2688 if (alpha == 1.0d0) then
2689 do i = i1, i2
2690 a(i+k,j) = a(i+k,j) * b(i)
2691 end do
2692 else
2693 do i = i1, i2
2694 a(i+k,j) = alpha * a(i+k,j) * b(i)
2695 end do
2696 end if
2697 end do
2698 end if
2699
2700
2701 ! End
2702 return
2703
2704 ! KL < 0
270510 continue
2706 call errmgr%report_error("band_diag_mtx_mult_cmplx", &
2707 "The number of subdiagonals must be at least 0.", &
2708 la_invalid_input_error)
2709 return
2710
2711 ! KU < 0
271220 continue
2713 call errmgr%report_error("band_diag_mtx_mult_cmplx", &
2714 "The number of superdiagonals must be at least 0.", &
2715 la_invalid_input_error)
2716 return
2717
2718 ! B is not sized correctly
271930 continue
2720 call errmgr%report_error("band_diag_mtx_mult_cmplx", &
2721 "Inner matrix dimensions do not agree.", &
2722 la_array_size_error)
2723 return
2724 end subroutine
2725
2726! ------------------------------------------------------------------------------
2727end submodule
A module providing explicit interfaces to BLAS routines.
Definition blas.f90:2
Provides a set of common linear algebra routines.
Definition linalg.f90:145