3submodule(
linalg) linalg_basic
11 module subroutine mtx_mult_mtx(transa, transb, alpha, a, b, beta, c, err)
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
20 real(real64),
parameter :: zero = 0.0d0
21 real(real64),
parameter :: one = 1.0d0
25 integer(int32) :: m, n, k, lda, ldb, flag
26 class(errors),
pointer :: errmgr
27 type(errors),
target :: deferr
28 character(len = :),
allocatable :: errmsg
49 if (
present(err))
then
58 if (
size(a, 2) /= m) flag = 4
60 if (
size(a, 1) /= m) flag = 4
63 if (
size(b, 2) /= k .or.
size(b, 1) /= n) flag = 5
65 if (
size(b, 1) /= k .or.
size(b, 2) /= n) flag = 5
69 allocate(
character(len = 256) :: errmsg)
71 "Matrix dimension mismatch. Input number ", flag, &
72 " was not sized correctly."
73 call errmgr%report_error(
"mtx_mult_mtx", errmsg, &
79 call dgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
86 module subroutine mtx_mult_vec(trans, alpha, a, b, beta, c, err)
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
97 integer(int32) :: m, n, flag
98 class(errors),
pointer :: errmgr
99 type(errors),
target :: deferr
100 character(len = :),
allocatable :: errmsg
107 if (
present(err))
then
116 if (
size(b) /= m)
then
118 else if (
size(c) /= n)
then
122 if (
size(b) /= n)
then
124 else if (
size(c) /= m)
then
130 allocate(
character(len = 256) :: errmsg)
132 "Matrix dimension mismatch. Input number ", flag, &
133 " was not sized correctly."
134 call errmgr%report_error(
"mtx_mult_vec", errmsg, &
140 call dgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
149 module subroutine cmtx_mult_mtx(opa, opb, alpha, a, b, beta, c, err)
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
158 real(real64),
parameter :: zero = 0.0d0
159 real(real64),
parameter :: one = 1.0d0
163 integer(int32) :: m, n, k, lda, ldb, flag
164 class(errors),
pointer :: errmgr
165 type(errors),
target :: deferr
166 character(len = :),
allocatable :: errmsg
171 if (opa == la_transpose)
then
175 else if (opa == la_hermitian_transpose)
then
184 if (opb == la_transpose)
then
187 else if (opb == la_hermitian_transpose)
then
194 if (
present(err))
then
202 if (opa == la_transpose .or. opa == la_hermitian_transpose)
then
203 if (
size(a, 2) /= m) flag = 4
205 if (
size(a, 1) /= m) flag = 4
207 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
208 if (
size(b, 2) /= k .or.
size(b, 1) /= n) flag = 5
210 if (
size(b, 1) /= k .or.
size(b, 2) /= n) flag = 5
214 allocate(
character(len = 256) :: errmsg)
216 "Matrix dimension mismatch. Input number ", flag, &
217 " was not sized correctly."
218 call errmgr%report_error(
"cmtx_mult_mtx", errmsg, &
224 call zgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
231 module subroutine cmtx_mult_vec(opa, alpha, a, b, beta, c, err)
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
242 integer(int32) :: m, n, flag
243 class(errors),
pointer :: errmgr
244 type(errors),
target :: deferr
245 character(len = :),
allocatable :: errmsg
250 if (opa == la_transpose)
then
252 else if (opa == la_hermitian_transpose)
then
257 if (
present(err))
then
265 if (opa == la_transpose .or. opa == la_hermitian_transpose)
then
266 if (
size(b) /= m)
then
268 else if (
size(c) /= n)
then
272 if (
size(b) /= n)
then
274 else if (
size(c) /= m)
then
280 allocate(
character(len = 256) :: errmsg)
282 "Matrix dimension mismatch. Input number ", flag, &
283 " was not sized correctly."
284 call errmgr%report_error(
"cmtx_mult_vec", errmsg, &
290 call zgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
299 module subroutine rank1_update_dbl(alpha, x, y, a, err)
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
307 real(real64),
parameter :: zero = 0.0d0
310 integer(int32) :: j, m, n
312 class(errors),
pointer :: errmgr
313 type(errors),
target :: deferr
318 if (
present(err))
then
325 if (
size(a, 1) /= m .or.
size(a, 2) /= n)
then
327 call errmgr%report_error(
"rank1_update_dbl", &
328 "Matrix dimension mismatch.", la_array_size_error)
334 if (y(j) /= zero)
then
336 a(:,j) = a(:,j) + temp * x
344 module subroutine rank1_update_cmplx(alpha, x, y, a, err)
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
352 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
355 integer(int32) :: j, m, n
356 complex(real64) :: temp
357 class(errors),
pointer :: errmgr
358 type(errors),
target :: deferr
363 if (
present(err))
then
370 if (
size(a, 1) /= m .or.
size(a, 2) /= n)
then
372 call errmgr%report_error(
"rank1_update_cmplx", &
373 "Matrix dimension mismatch.", la_array_size_error)
379 if (y(j) /= zero)
then
380 temp = alpha * conjg(y(j))
381 a(:,j) = a(:,j) + temp * x
389 module subroutine diag_mtx_mult_mtx(lside, trans, alpha, a, b, beta, c, err)
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
399 real(real64),
parameter :: zero = 0.0d0
400 real(real64),
parameter :: one = 1.0d0
403 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
405 class(errors),
pointer :: errmgr
406 type(errors),
target :: deferr
407 character(len = :),
allocatable :: errmsg
415 if (
present(err))
then
429 if (nrowb /= n .or. ncolb < k) flag = 5
432 if (nrowb < k .or. ncolb /= n) flag = 5
441 if (ncolb /= m .or. nrowb < k) flag = 5
444 if (nrowb /= m .or. ncolb < k) flag = 5
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), &
460 if (beta == zero)
then
462 else if (beta /= one)
then
473 if (beta == zero)
then
475 else if (beta /= one)
then
476 c(i,:) = beta * c(i,:)
479 c(i,:) = c(i,:) + temp * b(:,i)
484 if (beta == zero)
then
486 else if (beta /= one)
then
487 c(i,:) = beta * c(i,:)
490 c(i,:) = c(i,:) + temp * b(i,:)
496 if (beta == zero)
then
499 c(k+1:m,:) = beta * c(k+1:m,:)
506 if (beta == zero)
then
508 else if (beta /= one)
then
509 c(:,i) = beta * c(:,i)
512 c(:,i) = c(:,i) + temp * b(i,:)
517 if (beta == zero)
then
519 else if (beta /= one)
then
520 c(:,i) = beta * c(:,i)
523 c(:,i) = c(:,i) + temp * b(:,i)
529 if (beta == zero)
then
531 else if (beta /= one)
then
532 c(:,k+1:m) = beta * c(:,k+1:m)
542 module subroutine diag_mtx_mult_mtx2(lside, alpha, a, b, err)
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
551 real(real64),
parameter :: zero = 0.0d0
552 real(real64),
parameter :: one = 1.0d0
555 integer(int32) :: i, m, n, k
557 class(errors),
pointer :: errmgr
558 type(errors),
target :: deferr
564 if (
present(err))
then
571 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
573 call errmgr%report_error(
"diag_mtx_mult_mtx2", &
574 "Input number 3 is not sized correctly.", &
584 b(i,:) = temp * b(i,:)
586 if (m > k) b(k+1:m,:) = zero
591 b(:,i) = temp * b(:,i)
593 if (n > k) b(:,k+1:n) = zero
598 module subroutine diag_mtx_mult_mtx3(lside, trans, alpha, a, b, beta, c, err)
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
608 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
609 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
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
624 if (
present(err))
then
638 if (nrowb /= n .or. ncolb < k) flag = 5
641 if (nrowb < k .or. ncolb /= n) flag = 5
650 if (ncolb /= m .or. nrowb < k) flag = 5
653 if (nrowb /= m .or. ncolb < k) flag = 5
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), &
669 if (beta == zero)
then
671 else if (beta /= one)
then
682 if (beta == zero)
then
684 else if (beta /= one)
then
685 c(i,:) = beta * c(i,:)
688 c(i,:) = c(i,:) + temp * b(:,i)
693 if (beta == zero)
then
695 else if (beta /= one)
then
696 c(i,:) = beta * c(i,:)
699 c(i,:) = c(i,:) + temp * b(i,:)
705 if (beta == zero)
then
708 c(k+1:m,:) = beta * c(k+1:m,:)
715 if (beta == zero)
then
717 else if (beta /= one)
then
718 c(:,i) = beta * c(:,i)
721 c(:,i) = c(:,i) + temp * b(i,:)
726 if (beta == zero)
then
728 else if (beta /= one)
then
729 c(:,i) = beta * c(:,i)
732 c(:,i) = c(:,i) + temp * b(:,i)
738 if (beta == zero)
then
740 else if (beta /= one)
then
741 c(:,k+1:m) = beta * c(:,k+1:m)
751 module subroutine diag_mtx_mult_mtx4(lside, opb, alpha, a, b, beta, c, err)
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
762 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
763 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
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
778 if (
present(err))
then
790 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
792 if (nrowb /= n .or. ncolb < k) flag = 5
795 if (nrowb < k .or. ncolb /= n) flag = 5
802 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
804 if (ncolb /= m .or. nrowb < k) flag = 5
807 if (nrowb /= m .or. ncolb < k) flag = 5
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), &
823 if (beta == zero)
then
825 else if (beta /= one)
then
833 if (opb == la_transpose)
then
836 if (beta == zero)
then
838 else if (beta /= one)
then
839 c(i,:) = beta * c(i,:)
842 c(i,:) = c(i,:) + temp * b(:,i)
844 else if (opb == la_hermitian_transpose)
then
847 if (beta == zero)
then
849 else if (beta /= one)
then
850 c(i,:) = beta * c(i,:)
853 c(i,:) = c(i,:) + temp * conjg(b(:,i))
858 if (beta == zero)
then
860 else if (beta /= one)
then
861 c(i,:) = beta * c(i,:)
864 c(i,:) = c(i,:) + temp * b(i,:)
870 if (beta == zero)
then
873 c(k+1:m,:) = beta * c(k+1:m,:)
877 if (opb == la_transpose)
then
880 if (beta == zero)
then
882 else if (beta /= one)
then
883 c(:,i) = beta * c(:,i)
886 c(:,i) = c(:,i) + temp * b(i,:)
888 else if (opb == la_hermitian_transpose)
then
891 if (beta == zero)
then
893 else if (beta /= one)
then
894 c(:,i) = beta * c(:,i)
897 c(:,i) = c(:,i) + temp * conjg(b(i,:))
902 if (beta == zero)
then
904 else if (beta /= one)
then
905 c(:,i) = beta * c(:,i)
908 c(:,i) = c(:,i) + temp * b(:,i)
914 if (beta == zero)
then
916 else if (beta /= one)
then
917 c(:,k+1:m) = beta * c(:,k+1:m)
927 module subroutine diag_mtx_mult_mtx_cmplx(lside, opb, alpha, a, b, beta, c, err)
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
938 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
939 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
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
954 if (
present(err))
then
966 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
968 if (nrowb /= n .or. ncolb < k) flag = 5
971 if (nrowb < k .or. ncolb /= n) flag = 5
978 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
980 if (ncolb /= m .or. nrowb < k) flag = 5
983 if (nrowb /= m .or. ncolb < k) flag = 5
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), &
999 if (beta == zero)
then
1001 else if (beta /= one)
then
1009 if (opb == la_transpose)
then
1012 if (beta == zero)
then
1014 else if (beta /= one)
then
1015 c(i,:) = beta * c(i,:)
1018 c(i,:) = c(i,:) + temp * b(:,i)
1020 else if (opb == la_hermitian_transpose)
then
1023 if (beta == zero)
then
1025 else if (beta /= one)
then
1026 c(i,:) = beta * c(i,:)
1029 c(i,:) = c(i,:) + temp * conjg(b(:,i))
1034 if (beta == zero)
then
1036 else if (beta /= one)
then
1037 c(i,:) = beta * c(i,:)
1040 c(i,:) = c(i,:) + temp * b(i,:)
1046 if (beta == zero)
then
1049 c(k+1:m,:) = beta * c(k+1:m,:)
1053 if (opb == la_transpose)
then
1056 if (beta == zero)
then
1058 else if (beta /= one)
then
1059 c(:,i) = beta * c(:,i)
1062 c(:,i) = c(:,i) + temp * b(i,:)
1064 else if (opb == la_hermitian_transpose)
then
1067 if (beta == zero)
then
1069 else if (beta /= one)
then
1070 c(:,i) = beta * c(:,i)
1073 c(:,i) = c(:,i) + temp * conjg(b(i,:))
1078 if (beta == zero)
then
1080 else if (beta /= one)
then
1081 c(:,i) = beta * c(:,i)
1084 c(:,i) = c(:,i) + temp * b(:,i)
1090 if (beta == zero)
then
1092 else if (beta /= one)
then
1093 c(:,k+1:m) = beta * c(:,k+1:m)
1103 module subroutine diag_mtx_mult_mtx2_cmplx(lside, alpha, a, b, err)
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
1112 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1113 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1116 integer(int32) :: i, m, n, k
1117 complex(real64) :: temp
1118 class(errors),
pointer :: errmgr
1119 type(errors),
target :: deferr
1125 if (
present(err))
then
1132 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
1134 call errmgr%report_error(
"diag_mtx_mult_mtx2_cmplx", &
1135 "Input number 3 is not sized correctly.", &
1136 la_array_size_error)
1145 b(i,:) = temp * b(i,:)
1147 if (m > k) b(k+1:m,:) = zero
1152 b(:,i) = temp * b(:,i)
1154 if (n > k) b(:,k+1:n) = zero
1159 module subroutine diag_mtx_mult_mtx_mix(lside, opb, alpha, a, b, beta, c, err)
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
1170 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1171 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
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
1186 if (
present(err))
then
1198 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
1200 if (nrowb /= n .or. ncolb < k) flag = 5
1203 if (nrowb < k .or. ncolb /= n) flag = 5
1210 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
1212 if (ncolb /= m .or. nrowb < k) flag = 5
1215 if (nrowb /= m .or. ncolb < k) flag = 5
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)
1230 if (alpha == 0)
then
1231 if (beta == zero)
then
1233 else if (beta /= one)
then
1241 if (opb == la_transpose)
then
1244 if (beta == zero)
then
1246 else if (beta /= one)
then
1247 c(i,:) = beta * c(i,:)
1250 c(i,:) = c(i,:) + temp * b(:,i)
1252 else if (opb == la_hermitian_transpose)
then
1255 if (beta == zero)
then
1257 else if (beta /= one)
then
1258 c(i,:) = beta * c(i,:)
1261 c(i,:) = c(i,:) + temp * conjg(b(:,i))
1266 if (beta == zero)
then
1268 else if (beta /= one)
then
1269 c(i,:) = beta * c(i,:)
1272 c(i,:) = c(i,:) + temp * b(i,:)
1278 if (beta == zero)
then
1281 c(k+1:m,:) = beta * c(k+1:m,:)
1285 if (opb == la_transpose)
then
1288 if (beta == zero)
then
1290 else if (beta /= one)
then
1291 c(:,i) = beta * c(:,i)
1294 c(:,i) = c(:,i) + temp * b(i,:)
1296 else if (opb == la_hermitian_transpose)
then
1299 if (beta == zero)
then
1301 else if (beta /= one)
then
1302 c(:,i) = beta * c(:,i)
1305 c(:,i) = c(:,i) + temp * conjg(b(i,:))
1310 if (beta == zero)
then
1312 else if (beta /= one)
then
1313 c(:,i) = beta * c(:,i)
1316 c(:,i) = c(:,i) + temp * b(:,i)
1322 if (beta == zero)
then
1324 else if (beta /= one)
then
1325 c(:,k+1:m) = beta * c(:,k+1:m)
1335 module subroutine diag_mtx_mult_mtx2_mix(lside, alpha, a, b, err)
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
1344 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1345 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1348 integer(int32) :: i, m, n, k
1349 complex(real64) :: temp
1350 class(errors),
pointer :: errmgr
1351 type(errors),
target :: deferr
1357 if (
present(err))
then
1364 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
1366 call errmgr%report_error(
"diag_mtx_mult_mtx2_cmplx", &
1367 "Input number 3 is not sized correctly.", &
1368 la_array_size_error)
1377 b(i,:) = temp * b(i,:)
1379 if (m > k) b(k+1:m,:) = zero
1384 b(:,i) = temp * b(:,i)
1386 if (n > k) b(:,k+1:n) = zero
1393 pure module function trace_dbl(x) result(y)
1395 real(real64),
intent(in),
dimension(:,:) :: x
1399 real(real64),
parameter :: zero = 0.0d0
1402 integer(int32) :: i, m, n, mn
1417 pure module function trace_cmplx(x) result(y)
1419 complex(real64),
intent(in),
dimension(:,:) :: x
1420 complex(real64) :: y
1423 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1426 integer(int32) :: i, m, n, mn
1441 module function mtx_rank_dbl(a, tol, work, olwork, err) result(rnk)
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
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
1466 if (
present(err))
then
1474 call dgesvd(
'N',
'N', m, n, a, m, dummy, dummy, m, dummy, n, temp, &
1476 lwork = int(temp(1), int32) + mn
1477 if (
present(olwork))
then
1483 if (
present(work))
then
1484 if (
size(work) < lwork)
then
1486 call errmgr%report_error(
"mtx_rank", &
1487 "Incorrectly sized input array WORK, argument 5.", &
1488 la_array_size_error)
1491 wptr => work(1:lwork)
1493 allocate(wrk(lwork), stat = istat)
1494 if (istat /= 0)
then
1496 call errmgr%report_error(
"mtx_rank", &
1497 "Insufficient memory available.", &
1498 la_out_of_memory_error)
1504 w => wptr(mn+1:lwork)
1507 call dgesvd(
'N',
'N', m, n, a, m, s, dummy, m, dummy, n, w, &
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)
1518 tref = max(m, n) * epsilon(t) * s(1)
1519 if (
present(tol))
then
1524 if (t < smlnum)
then
1545 module function mtx_rank_cmplx(a, tol, work, olwork, rwork, err) result(rnk)
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
1557 function dlamch(cmach)
result(x)
1558 use,
intrinsic :: iso_fortran_env, only : real64
1559 character,
intent(in) :: cmach
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
1584 if (
present(err))
then
1591 call zgesvd(
'N',
'N', m, n, a, m, dummy, cdummy, m, cdummy, n, temp, &
1593 lwork = int(temp(1), int32)
1594 if (
present(olwork))
then
1600 if (
present(work))
then
1601 if (
size(work) < lwork)
then
1603 call errmgr%report_error(
"mtx_rank_cmplx", &
1604 "Incorrectly sized input array WORK, argument 5.", &
1605 la_array_size_error)
1608 wptr => work(1:lwork)
1610 allocate(wrk(lwork), stat = istat)
1611 if (istat /= 0)
then
1613 call errmgr%report_error(
"mtx_rank_cmplx", &
1614 "Insufficient memory available.", &
1615 la_out_of_memory_error)
1621 if (
present(rwork))
then
1622 if (
size(rwork) < lrwork)
then
1624 call errmgr%report_error(
"mtx_rank_cmplx", &
1625 "Incorrectly sized input array RWORK.", &
1626 la_array_size_error)
1629 rwptr => rwork(1:lrwork)
1631 allocate(rwrk(lrwork), stat = istat)
1632 if (istat /= 0)
then
1637 rw => rwptr(mn+1:lrwork)
1640 call zgesvd(
'N',
'N', m, n, a, m, s, cdummy, m, cdummy, n, wptr, &
1641 lwork - mn, rw, flag)
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)
1651 tref = max(m, n) * epsilon(t) * s(1)
1652 if (
present(tol))
then
1657 if (t < smlnum)
then
1678 module function det_dbl(a, iwork, err) result(x)
1680 real(real64),
intent(inout),
dimension(:,:) :: a
1681 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1682 class(errors),
intent(inout),
optional,
target :: err
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
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
1702 if (
present(err))
then
1709 if (
size(a, 2) /= n)
then
1710 call errmgr%report_error(
"det", &
1711 "The supplied matrix must be square.", la_array_size_error)
1716 if (
present(iwork))
then
1717 if (
size(iwork) < n)
then
1719 call errmgr%report_error(
"det", &
1720 "Incorrectly sized input array IWORK, argument 2.", &
1721 la_array_size_error)
1726 allocate(iwrk(n), stat = istat)
1727 if (istat /= 0)
then
1729 call errmgr%report_error(
"det", &
1730 "Insufficient memory available.", &
1731 la_out_of_memory_error)
1738 call dgetrf(n, n, a, n, ipvt, flag)
1749 if (ipvt(i) /= i) temp = -temp
1751 temp = a(i,i) * temp
1752 if (temp == zero)
then
1757 do while (abs(temp) < one)
1762 do while (abs(temp) > ten)
1771 module function det_cmplx(a, iwork, err) result(x)
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
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
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
1797 if (
present(err))
then
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)
1811 if (
present(iwork))
then
1812 if (
size(iwork) < n)
then
1814 call errmgr%report_error(
"det_cmplx", &
1815 "Incorrectly sized input array IWORK, argument 2.", &
1816 la_array_size_error)
1821 allocate(iwrk(n), stat = istat)
1822 if (istat /= 0)
then
1824 call errmgr%report_error(
"det_cmplx", &
1825 "Insufficient memory available.", &
1826 la_out_of_memory_error)
1833 call zgetrf(n, n, a, n, ipvt, flag)
1844 if (ipvt(i) /= i) temp = -temp
1846 temp = a(i,i) * temp
1847 if (temp == zero)
then
1852 do while (abs(temp) < real_one)
1857 do while (abs(temp) > real_ten)
1868 module subroutine swap_dbl(x, y, err)
1870 real(real64),
intent(inout),
dimension(:) :: x, y
1871 class(errors),
intent(inout),
optional,
target :: err
1874 integer(int32) :: i, n
1875 real(real64) :: temp
1876 class(errors),
pointer :: errmgr
1877 type(errors),
target :: deferr
1881 if (
present(err))
then
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)
1904 module subroutine swap_cmplx(x, y, err)
1906 complex(real64),
intent(inout),
dimension(:) :: x, y
1907 class(errors),
intent(inout),
optional,
target :: err
1910 integer(int32) :: i, n
1911 complex(real64) :: temp
1912 class(errors),
pointer :: errmgr
1913 type(errors),
target :: deferr
1917 if (
present(err))
then
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)
1942 module subroutine recip_mult_array_dbl(a, x)
1944 real(real64),
intent(in) :: a
1945 real(real64),
intent(inout),
dimension(:) :: x
1949 function dlamch(cmach)
result(x)
1950 use,
intrinsic :: iso_fortran_env, only : real64
1951 character,
intent(in) :: cmach
1957 real(real64),
parameter :: zero = 0.0d0
1958 real(real64),
parameter :: one = 1.0d0
1959 real(real64),
parameter :: twotho = 2.0d3
1963 real(real64) :: bignum, cden, cden1, cnum, cnum1, mul, smlnum
1967 bignum = one / smlnum
1968 if (log10(bignum) > twotho)
then
1969 smlnum = sqrt(smlnum)
1970 bignum = sqrt(bignum)
1979 cden1 = cden * smlnum
1980 cnum1 = cnum / bignum
1981 if (abs(cden1) > abs(cnum) .and. cnum /= zero)
then
1985 else if (abs(cnum1) > abs(cden))
then
2005 module subroutine tri_mtx_mult_dbl(upper, alpha, a, beta, b, err)
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
2014 real(real64),
parameter :: zero = 0.0d0
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
2027 if (
present(err))
then
2035 if (
size(a, 2) /= n)
then
2038 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
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)
2057 if (beta == zero)
then
2062 temp = temp + a(k,i) * a(k,j)
2066 if (i /= j) b(j,i) = temp
2074 temp = temp + a(k,i) * a(k,j)
2077 b(i,j) = temp + beta * b(i,j)
2078 if (i /= j) b(j,i) = temp + beta * b(j,i)
2084 if (beta == zero)
then
2089 temp = temp + a(i,k) * a(j,k)
2093 if (i /= j) b(j,i) = temp
2101 temp = temp + a(i,k) * a(j,k)
2104 b(i,j) = temp + beta * b(i,j)
2105 if (i /= j) b(j,i) = temp + beta * b(j,i)
2112100
format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2116 module subroutine tri_mtx_mult_cmplx(upper, alpha, a, beta, b, err)
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
2125 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
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
2138 if (
present(err))
then
2146 if (
size(a, 2) /= n)
then
2149 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
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)
2168 if (beta == zero)
then
2173 temp = temp + a(k,i) * a(k,j)
2177 if (i /= j) b(j,i) = temp
2185 temp = temp + a(k,i) * a(k,j)
2188 b(i,j) = temp + beta * b(i,j)
2189 if (i /= j) b(j,i) = temp + beta * b(j,i)
2195 if (beta == zero)
then
2200 temp = temp + a(i,k) * a(j,k)
2204 if (i /= j) b(j,i) = temp
2212 temp = temp + a(i,k) * a(j,k)
2215 b(i,j) = temp + beta * b(i,j)
2216 if (i /= j) b(j,i) = temp + beta * b(j,i)
2223100
format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2229 module subroutine band_mtx_vec_mult_dbl(trans, kl, ku, alpha, a, x, beta, &
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
2241 integer(int32) :: m, n
2242 class(errors),
pointer :: errmgr
2243 type(errors),
target :: deferr
2246 if (
present(err))
then
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
2267 call dgbmv(
"T", m, n, kl, ku, alpha, a,
size(a, 1), x, 1, beta, y, 1)
2269 call dgbmv(
"N", m, n, kl, ku, alpha, a,
size(a, 1), x, 1, beta, y, 1)
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)
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)
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)
2298 module subroutine band_mtx_vec_mult_cmplx(trans, kl, ku, alpha, a, x, &
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
2312 integer(int32) :: m, n
2313 class(errors),
pointer :: errmgr
2314 type(errors),
target :: deferr
2317 if (
present(err))
then
2322 if (trans == la_transpose)
then
2325 else if (trans == la_hermitian_transpose)
then
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
2347 call zgbmv(op, m, n, kl, ku, alpha, a,
size(a, 1), x, 1, beta, y, 1)
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)
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)
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)
2375 module subroutine band_to_full_mtx_dbl(kl, ku, b, f, err)
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
2383 real(real64),
parameter :: zero = 0.0d0
2386 class(errors),
pointer :: errmgr
2387 type(errors),
target :: deferr
2388 integer(int32) :: i, j, k, m, n, i1, i2
2391 if (
present(err))
then
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
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)
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)
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)
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)
2455 module subroutine band_to_full_mtx_cmplx(kl, ku, b, f, err)
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
2463 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
2466 class(errors),
pointer :: errmgr
2467 type(errors),
target :: deferr
2468 integer(int32) :: i, j, k, m, n, i1, i2
2471 if (
present(err))
then
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
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)
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)
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)
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)
2535 module subroutine band_diag_mtx_mult_dbl(left, m, kl, ku, alpha, a, b, err)
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
2545 real(real64),
parameter :: one = 1.0d0
2548 class(errors),
pointer :: errmgr
2549 type(errors),
target :: deferr
2550 integer(int32) :: i, i1, i2, j, k, n
2551 real(real64) :: temp
2554 if (
present(err))
then
2562 if (kl < 0)
go to 10
2563 if (ku < 0)
go to 20
2565 if (
size(b) /= n)
go to 30
2567 if (
size(b) < m)
go to 30
2575 i1 = max(1, j - ku) + k
2576 i2 = min(m, j + kl) + k
2577 if (alpha == one)
then
2583 a(i,j) = a(i,j) * temp
2592 if (alpha == 1.0d0)
then
2594 a(i+k,j) = a(i+k,j) * b(i)
2598 a(i+k,j) = alpha * a(i+k,j) * b(i)
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)
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)
2624 call errmgr%report_error(
"band_diag_mtx_mult_dbl", &
2625 "Inner matrix dimensions do not agree.", &
2626 la_array_size_error)
2631 module subroutine band_diag_mtx_mult_cmplx(left, m, kl, ku, alpha, a, b, err)
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
2641 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
2644 class(errors),
pointer :: errmgr
2645 type(errors),
target :: deferr
2646 integer(int32) :: i, i1, i2, j, k, n
2647 complex(real64) :: temp
2650 if (
present(err))
then
2658 if (kl < 0)
go to 10
2659 if (ku < 0)
go to 20
2661 if (
size(b) /= n)
go to 30
2663 if (
size(b) < m)
go to 30
2671 i1 = max(1, j - ku) + k
2672 i2 = min(m, j + kl) + k
2673 if (alpha == one)
then
2679 a(i,j) = a(i,j) * temp
2688 if (alpha == 1.0d0)
then
2690 a(i+k,j) = a(i+k,j) * b(i)
2694 a(i+k,j) = alpha * a(i+k,j) * b(i)
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)
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)
2720 call errmgr%report_error(
"band_diag_mtx_mult_cmplx", &
2721 "Inner matrix dimensions do not agree.", &
2722 la_array_size_error)
A module providing explicit interfaces to BLAS routines.
Provides a set of common linear algebra routines.