3submodule(
linalg) linalg_basic
8 module subroutine mtx_mult_mtx(transa, transb, alpha, a, b, beta, c, err)
10 logical,
intent(in) :: transa, transb
11 real(real64),
intent(in) :: alpha, beta
12 real(real64),
intent(in),
dimension(:,:) :: a, b
13 real(real64),
intent(inout),
dimension(:,:) :: c
14 class(errors),
intent(inout),
optional,
target :: err
17 real(real64),
parameter :: zero = 0.0d0
18 real(real64),
parameter :: one = 1.0d0
22 integer(int32) :: m, n, k, lda, ldb, flag
23 class(errors),
pointer :: errmgr
24 type(errors),
target :: deferr
25 character(len = :),
allocatable :: errmsg
46 if (
present(err))
then
55 if (
size(a, 2) /= m) flag = 4
57 if (
size(a, 1) /= m) flag = 4
60 if (
size(b, 2) /= k .or.
size(b, 1) /= n) flag = 5
62 if (
size(b, 1) /= k .or.
size(b, 2) /= n) flag = 5
66 allocate(
character(len = 256) :: errmsg)
68 "Matrix dimension mismatch. Input number ", flag, &
69 " was not sized correctly."
70 call errmgr%report_error(
"mtx_mult_mtx", errmsg, &
76 call dgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
83 module subroutine mtx_mult_vec(trans, alpha, a, b, beta, c, err)
85 logical,
intent(in) :: trans
86 real(real64),
intent(in) :: alpha, beta
87 real(real64),
intent(in),
dimension(:,:) :: a
88 real(real64),
intent(in),
dimension(:) :: b
89 real(real64),
intent(inout),
dimension(:) :: c
90 class(errors),
intent(inout),
optional,
target :: err
94 integer(int32) :: m, n, flag
95 class(errors),
pointer :: errmgr
96 type(errors),
target :: deferr
97 character(len = :),
allocatable :: errmsg
104 if (
present(err))
then
113 if (
size(b) /= m)
then
115 else if (
size(c) /= n)
then
119 if (
size(b) /= n)
then
121 else if (
size(c) /= m)
then
127 allocate(
character(len = 256) :: errmsg)
129 "Matrix dimension mismatch. Input number ", flag, &
130 " was not sized correctly."
131 call errmgr%report_error(
"mtx_mult_vec", errmsg, &
137 call dgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
146 module subroutine cmtx_mult_mtx(opa, opb, alpha, a, b, beta, c, err)
148 integer(int32),
intent(in) :: opa, opb
149 complex(real64),
intent(in) :: alpha, beta
150 complex(real64),
intent(in),
dimension(:,:) :: a, b
151 complex(real64),
intent(inout),
dimension(:,:) :: c
152 class(errors),
intent(inout),
optional,
target :: err
155 real(real64),
parameter :: zero = 0.0d0
156 real(real64),
parameter :: one = 1.0d0
160 integer(int32) :: m, n, k, lda, ldb, flag
161 class(errors),
pointer :: errmgr
162 type(errors),
target :: deferr
163 character(len = :),
allocatable :: errmsg
168 if (opa == la_transpose)
then
172 else if (opa == la_hermitian_transpose)
then
181 if (opb == la_transpose)
then
184 else if (opb == la_hermitian_transpose)
then
191 if (
present(err))
then
199 if (opa == la_transpose .or. opa == la_hermitian_transpose)
then
200 if (
size(a, 2) /= m) flag = 4
202 if (
size(a, 1) /= m) flag = 4
204 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
205 if (
size(b, 2) /= k .or.
size(b, 1) /= n) flag = 5
207 if (
size(b, 1) /= k .or.
size(b, 2) /= n) flag = 5
211 allocate(
character(len = 256) :: errmsg)
213 "Matrix dimension mismatch. Input number ", flag, &
214 " was not sized correctly."
215 call errmgr%report_error(
"cmtx_mult_mtx", errmsg, &
221 call zgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
228 module subroutine cmtx_mult_vec(opa, alpha, a, b, beta, c, err)
230 integer(int32),
intent(in) :: opa
231 complex(real64),
intent(in) :: alpha, beta
232 complex(real64),
intent(in),
dimension(:,:) :: a
233 complex(real64),
intent(in),
dimension(:) :: b
234 complex(real64),
intent(inout),
dimension(:) :: c
235 class(errors),
intent(inout),
optional,
target :: err
239 integer(int32) :: m, n, flag
240 class(errors),
pointer :: errmgr
241 type(errors),
target :: deferr
242 character(len = :),
allocatable :: errmsg
247 if (opa == la_transpose)
then
249 else if (opa == la_hermitian_transpose)
then
254 if (
present(err))
then
262 if (opa == la_transpose .or. opa == la_hermitian_transpose)
then
263 if (
size(b) /= m)
then
265 else if (
size(c) /= n)
then
269 if (
size(b) /= n)
then
271 else if (
size(c) /= m)
then
277 allocate(
character(len = 256) :: errmsg)
279 "Matrix dimension mismatch. Input number ", flag, &
280 " was not sized correctly."
281 call errmgr%report_error(
"cmtx_mult_vec", errmsg, &
287 call zgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
296 module subroutine rank1_update_dbl(alpha, x, y, a, err)
298 real(real64),
intent(in) :: alpha
299 real(real64),
intent(in),
dimension(:) :: x, y
300 real(real64),
intent(inout),
dimension(:,:) :: a
301 class(errors),
intent(inout),
optional,
target :: err
304 real(real64),
parameter :: zero = 0.0d0
307 integer(int32) :: j, m, n
309 class(errors),
pointer :: errmgr
310 type(errors),
target :: deferr
315 if (
present(err))
then
322 if (
size(a, 1) /= m .or.
size(a, 2) /= n)
then
324 call errmgr%report_error(
"rank1_update_dbl", &
325 "Matrix dimension mismatch.", la_array_size_error)
331 if (y(j) /= zero)
then
333 a(:,j) = a(:,j) + temp * x
341 module subroutine rank1_update_cmplx(alpha, x, y, a, err)
343 complex(real64),
intent(in) :: alpha
344 complex(real64),
intent(in),
dimension(:) :: x, y
345 complex(real64),
intent(inout),
dimension(:,:) :: a
346 class(errors),
intent(inout),
optional,
target :: err
349 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
352 integer(int32) :: j, m, n
353 complex(real64) :: temp
354 class(errors),
pointer :: errmgr
355 type(errors),
target :: deferr
360 if (
present(err))
then
367 if (
size(a, 1) /= m .or.
size(a, 2) /= n)
then
369 call errmgr%report_error(
"rank1_update_cmplx", &
370 "Matrix dimension mismatch.", la_array_size_error)
376 if (y(j) /= zero)
then
377 temp = alpha * conjg(y(j))
378 a(:,j) = a(:,j) + temp * x
386 module subroutine diag_mtx_mult_mtx(lside, trans, alpha, a, b, beta, c, err)
388 logical,
intent(in) :: lside, trans
389 real(real64) :: alpha, beta
390 real(real64),
intent(in),
dimension(:) :: a
391 real(real64),
intent(in),
dimension(:,:) :: b
392 real(real64),
intent(inout),
dimension(:,:) :: c
393 class(errors),
intent(inout),
optional,
target :: err
396 real(real64),
parameter :: zero = 0.0d0
397 real(real64),
parameter :: one = 1.0d0
400 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
402 class(errors),
pointer :: errmgr
403 type(errors),
target :: deferr
404 character(len = :),
allocatable :: errmsg
412 if (
present(err))
then
426 if (nrowb /= n .or. ncolb < k) flag = 5
429 if (nrowb < k .or. ncolb /= n) flag = 5
438 if (ncolb /= m .or. nrowb < k) flag = 5
441 if (nrowb /= m .or. ncolb < k) flag = 5
447 allocate(
character(len = 256) :: errmsg)
448 write(errmsg, 100)
"Input number ", flag, &
449 " is not sized correctly."
450 call errmgr%report_error(
"diag_mtx_mult_mtx", trim(errmsg), &
457 if (beta == zero)
then
459 else if (beta /= one)
then
470 if (beta == zero)
then
472 else if (beta /= one)
then
473 c(i,:) = beta * c(i,:)
476 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
481 if (beta == zero)
then
483 else if (beta /= one)
then
484 c(i,:) = beta * c(i,:)
487 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
493 if (beta == zero)
then
496 c(k+1:m,:) = beta * c(k+1:m,:)
503 if (beta == zero)
then
505 else if (beta /= one)
then
506 c(:,i) = beta * c(:,i)
509 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
514 if (beta == zero)
then
516 else if (beta /= one)
then
517 c(:,i) = beta * c(:,i)
520 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
526 if (beta == zero)
then
528 else if (beta /= one)
then
529 c(:,k+1:m) = beta * c(:,k+1:m)
539 module subroutine diag_mtx_mult_mtx2(lside, alpha, a, b, err)
541 logical,
intent(in) :: lside
542 real(real64),
intent(in) :: alpha
543 real(real64),
intent(in),
dimension(:) :: a
544 real(real64),
intent(inout),
dimension(:,:) :: b
545 class(errors),
intent(inout),
optional,
target :: err
548 real(real64),
parameter :: zero = 0.0d0
549 real(real64),
parameter :: one = 1.0d0
552 integer(int32) :: i, m, n, k
554 class(errors),
pointer :: errmgr
555 type(errors),
target :: deferr
561 if (
present(err))
then
568 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
570 call errmgr%report_error(
"diag_mtx_mult_mtx2", &
571 "Input number 3 is not sized correctly.", &
581 if (temp /= one) b(i,:) = temp * b(i,:)
583 if (m > k) b(k+1:m,:) = zero
588 if (temp /= one) b(:,i) = temp * b(:,i)
590 if (n > k) b(:,k+1:n) = zero
595 module subroutine diag_mtx_mult_mtx3(lside, trans, alpha, a, b, beta, c, err)
597 logical,
intent(in) :: lside, trans
598 real(real64) :: alpha, beta
599 complex(real64),
intent(in),
dimension(:) :: a
600 real(real64),
intent(in),
dimension(:,:) :: b
601 complex(real64),
intent(inout),
dimension(:,:) :: c
602 class(errors),
intent(inout),
optional,
target :: err
605 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
606 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
609 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
610 complex(real64) :: temp
611 class(errors),
pointer :: errmgr
612 type(errors),
target :: deferr
613 character(len = :),
allocatable :: errmsg
621 if (
present(err))
then
635 if (nrowb /= n .or. ncolb < k) flag = 5
638 if (nrowb < k .or. ncolb /= n) flag = 5
647 if (ncolb /= m .or. nrowb < k) flag = 5
650 if (nrowb /= m .or. ncolb < k) flag = 5
656 allocate(
character(len = 256) :: errmsg)
657 write(errmsg, 100)
"Input number ", flag, &
658 " is not sized correctly."
659 call errmgr%report_error(
"diag_mtx_mult_mtx3", trim(errmsg), &
666 if (beta == zero)
then
668 else if (beta /= one)
then
679 if (beta == zero)
then
681 else if (beta /= one)
then
682 c(i,:) = beta * c(i,:)
685 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
690 if (beta == zero)
then
692 else if (beta /= one)
then
693 c(i,:) = beta * c(i,:)
696 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
702 if (beta == zero)
then
705 c(k+1:m,:) = beta * c(k+1:m,:)
712 if (beta == zero)
then
714 else if (beta /= one)
then
715 c(:,i) = beta * c(:,i)
718 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
723 if (beta == zero)
then
725 else if (beta /= one)
then
726 c(:,i) = beta * c(:,i)
729 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
735 if (beta == zero)
then
737 else if (beta /= one)
then
738 c(:,k+1:m) = beta * c(:,k+1:m)
748 module subroutine diag_mtx_mult_mtx4(lside, opb, alpha, a, b, beta, c, err)
750 logical,
intent(in) :: lside
751 integer(int32),
intent(in) :: opb
752 real(real64) :: alpha, beta
753 complex(real64),
intent(in),
dimension(:) :: a
754 complex(real64),
intent(in),
dimension(:,:) :: b
755 complex(real64),
intent(inout),
dimension(:,:) :: c
756 class(errors),
intent(inout),
optional,
target :: err
759 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
760 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
763 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
764 complex(real64) :: temp
765 class(errors),
pointer :: errmgr
766 type(errors),
target :: deferr
767 character(len = :),
allocatable :: errmsg
775 if (
present(err))
then
787 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
789 if (nrowb /= n .or. ncolb < k) flag = 5
792 if (nrowb < k .or. ncolb /= n) flag = 5
799 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
801 if (ncolb /= m .or. nrowb < k) flag = 5
804 if (nrowb /= m .or. ncolb < k) flag = 5
810 allocate(
character(len = 256) :: errmsg)
811 write(errmsg, 100)
"Input number ", flag, &
812 " is not sized correctly."
813 call errmgr%report_error(
"diag_mtx_mult_mtx4", trim(errmsg), &
820 if (beta == zero)
then
822 else if (beta /= one)
then
830 if (opb == la_transpose)
then
833 if (beta == zero)
then
835 else if (beta /= one)
then
836 c(i,:) = beta * c(i,:)
839 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
841 else if (opb == la_hermitian_transpose)
then
844 if (beta == zero)
then
846 else if (beta /= one)
then
847 c(i,:) = beta * c(i,:)
850 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
855 if (beta == zero)
then
857 else if (beta /= one)
then
858 c(i,:) = beta * c(i,:)
861 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
867 if (beta == zero)
then
870 c(k+1:m,:) = beta * c(k+1:m,:)
874 if (opb == la_transpose)
then
877 if (beta == zero)
then
879 else if (beta /= one)
then
880 c(:,i) = beta * c(:,i)
883 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
885 else if (opb == la_hermitian_transpose)
then
888 if (beta == zero)
then
890 else if (beta /= one)
then
891 c(:,i) = beta * c(:,i)
894 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
899 if (beta == zero)
then
901 else if (beta /= one)
then
902 c(:,i) = beta * c(:,i)
905 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
911 if (beta == zero)
then
913 else if (beta /= one)
then
914 c(:,k+1:m) = beta * c(:,k+1:m)
924 module subroutine diag_mtx_mult_mtx_cmplx(lside, opb, alpha, a, b, beta, c, err)
926 logical,
intent(in) :: lside
927 integer(int32),
intent(in) :: opb
928 complex(real64) :: alpha, beta
929 complex(real64),
intent(in),
dimension(:) :: a
930 complex(real64),
intent(in),
dimension(:,:) :: b
931 complex(real64),
intent(inout),
dimension(:,:) :: c
932 class(errors),
intent(inout),
optional,
target :: err
935 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
936 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
939 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
940 complex(real64) :: temp
941 class(errors),
pointer :: errmgr
942 type(errors),
target :: deferr
943 character(len = :),
allocatable :: errmsg
951 if (
present(err))
then
963 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
965 if (nrowb /= n .or. ncolb < k) flag = 5
968 if (nrowb < k .or. ncolb /= n) flag = 5
975 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
977 if (ncolb /= m .or. nrowb < k) flag = 5
980 if (nrowb /= m .or. ncolb < k) flag = 5
986 allocate(
character(len = 256) :: errmsg)
987 write(errmsg, 100)
"Input number ", flag, &
988 " is not sized correctly."
989 call errmgr%report_error(
"diag_mtx_mult_mtx_cmplx", trim(errmsg), &
996 if (beta == zero)
then
998 else if (beta /= one)
then
1006 if (opb == la_transpose)
then
1009 if (beta == zero)
then
1011 else if (beta /= one)
then
1012 c(i,:) = beta * c(i,:)
1015 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
1017 else if (opb == la_hermitian_transpose)
then
1020 if (beta == zero)
then
1022 else if (beta /= one)
then
1023 c(i,:) = beta * c(i,:)
1026 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
1031 if (beta == zero)
then
1033 else if (beta /= one)
then
1034 c(i,:) = beta * c(i,:)
1037 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
1043 if (beta == zero)
then
1046 c(k+1:m,:) = beta * c(k+1:m,:)
1050 if (opb == la_transpose)
then
1053 if (beta == zero)
then
1055 else if (beta /= one)
then
1056 c(:,i) = beta * c(:,i)
1059 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
1061 else if (opb == la_hermitian_transpose)
then
1064 if (beta == zero)
then
1066 else if (beta /= one)
then
1067 c(:,i) = beta * c(:,i)
1070 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
1075 if (beta == zero)
then
1077 else if (beta /= one)
then
1078 c(:,i) = beta * c(:,i)
1081 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
1087 if (beta == zero)
then
1089 else if (beta /= one)
then
1090 c(:,k+1:m) = beta * c(:,k+1:m)
1100 module subroutine diag_mtx_mult_mtx2_cmplx(lside, alpha, a, b, err)
1102 logical,
intent(in) :: lside
1103 complex(real64),
intent(in) :: alpha
1104 complex(real64),
intent(in),
dimension(:) :: a
1105 complex(real64),
intent(inout),
dimension(:,:) :: b
1106 class(errors),
intent(inout),
optional,
target :: err
1109 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1110 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1113 integer(int32) :: i, m, n, k
1114 complex(real64) :: temp
1115 class(errors),
pointer :: errmgr
1116 type(errors),
target :: deferr
1122 if (
present(err))
then
1129 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
1131 call errmgr%report_error(
"diag_mtx_mult_mtx2_cmplx", &
1132 "Input number 3 is not sized correctly.", &
1133 la_array_size_error)
1142 if (temp /= one) b(i,:) = temp * b(i,:)
1144 if (m > k) b(k+1:m,:) = zero
1149 if (temp /= one) b(:,i) = temp * b(:,i)
1151 if (n > k) b(:,k+1:n) = zero
1156 module subroutine diag_mtx_mult_mtx_mix(lside, opb, alpha, a, b, beta, c, err)
1158 logical,
intent(in) :: lside
1159 integer(int32),
intent(in) :: opb
1160 complex(real64) :: alpha, beta
1161 real(real64),
intent(in),
dimension(:) :: a
1162 complex(real64),
intent(in),
dimension(:,:) :: b
1163 complex(real64),
intent(inout),
dimension(:,:) :: c
1164 class(errors),
intent(inout),
optional,
target :: err
1167 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1168 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1171 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
1172 complex(real64) :: temp
1173 class(errors),
pointer :: errmgr
1174 type(errors),
target :: deferr
1175 character(len = :),
allocatable :: errmsg
1183 if (
present(err))
then
1195 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
1197 if (nrowb /= n .or. ncolb < k) flag = 5
1200 if (nrowb < k .or. ncolb /= n) flag = 5
1207 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
1209 if (ncolb /= m .or. nrowb < k) flag = 5
1212 if (nrowb /= m .or. ncolb < k) flag = 5
1218 allocate(
character(len = 256) :: errmsg)
1219 write(errmsg, 100)
"Input number ", flag, &
1220 " is not sized correctly."
1221 call errmgr%report_error(
"diag_mtx_mult_mtx_mix", trim(errmsg), &
1222 la_array_size_error)
1227 if (alpha == 0)
then
1228 if (beta == zero)
then
1230 else if (beta /= one)
then
1238 if (opb == la_transpose)
then
1241 if (beta == zero)
then
1243 else if (beta /= one)
then
1244 c(i,:) = beta * c(i,:)
1247 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
1249 else if (opb == la_hermitian_transpose)
then
1252 if (beta == zero)
then
1254 else if (beta /= one)
then
1255 c(i,:) = beta * c(i,:)
1258 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
1263 if (beta == zero)
then
1265 else if (beta /= one)
then
1266 c(i,:) = beta * c(i,:)
1269 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
1275 if (beta == zero)
then
1278 c(k+1:m,:) = beta * c(k+1:m,:)
1282 if (opb == la_transpose)
then
1285 if (beta == zero)
then
1287 else if (beta /= one)
then
1288 c(:,i) = beta * c(:,i)
1291 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
1293 else if (opb == la_hermitian_transpose)
then
1296 if (beta == zero)
then
1298 else if (beta /= one)
then
1299 c(:,i) = beta * c(:,i)
1302 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
1307 if (beta == zero)
then
1309 else if (beta /= one)
then
1310 c(:,i) = beta * c(:,i)
1313 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
1319 if (beta == zero)
then
1321 else if (beta /= one)
then
1322 c(:,k+1:m) = beta * c(:,k+1:m)
1332 module subroutine diag_mtx_mult_mtx2_mix(lside, alpha, a, b, err)
1334 logical,
intent(in) :: lside
1335 complex(real64),
intent(in) :: alpha
1336 real(real64),
intent(in),
dimension(:) :: a
1337 complex(real64),
intent(inout),
dimension(:,:) :: b
1338 class(errors),
intent(inout),
optional,
target :: err
1341 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1342 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1345 integer(int32) :: i, m, n, k
1346 complex(real64) :: temp
1347 class(errors),
pointer :: errmgr
1348 type(errors),
target :: deferr
1354 if (
present(err))
then
1361 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
1363 call errmgr%report_error(
"diag_mtx_mult_mtx2_cmplx", &
1364 "Input number 3 is not sized correctly.", &
1365 la_array_size_error)
1374 if (temp /= one) b(i,:) = temp * b(i,:)
1376 if (m > k) b(k+1:m,:) = zero
1381 if (temp /= one) b(:,i) = temp * b(:,i)
1383 if (n > k) b(:,k+1:n) = zero
1390 pure module function trace_dbl(x) result(y)
1392 real(real64),
intent(in),
dimension(:,:) :: x
1396 real(real64),
parameter :: zero = 0.0d0
1399 integer(int32) :: i, m, n, mn
1414 pure module function trace_cmplx(x) result(y)
1416 complex(real64),
intent(in),
dimension(:,:) :: x
1417 complex(real64) :: y
1420 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1423 integer(int32) :: i, m, n, mn
1438 module function mtx_rank_dbl(a, tol, work, olwork, err) result(rnk)
1440 real(real64),
intent(inout),
dimension(:,:) :: a
1441 real(real64),
intent(in),
optional :: tol
1442 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
1443 integer(int32),
intent(out),
optional :: olwork
1444 class(errors),
intent(inout),
optional,
target :: err
1445 integer(int32) :: rnk
1449 function dlamch(cmach)
result(x)
1450 use,
intrinsic :: iso_fortran_env, only : real64
1451 character,
intent(in) :: cmach
1457 integer(int32) :: i, m, n, mn, istat, lwork, flag
1458 real(real64),
pointer,
dimension(:) :: wptr, s, w
1459 real(real64),
allocatable,
target,
dimension(:) :: wrk
1460 real(real64) :: t, tref, smlnum
1461 real(real64),
dimension(1) :: dummy, temp
1462 class(errors),
pointer :: errmgr
1463 type(errors),
target :: deferr
1464 character(len = :),
allocatable :: errmsg
1470 smlnum = dlamch(
's')
1472 if (
present(err))
then
1480 call dgesvd(
'N',
'N', m, n, a, m, dummy, dummy, m, dummy, n, temp, &
1482 lwork = int(temp(1), int32) + mn
1483 if (
present(olwork))
then
1489 if (
present(work))
then
1490 if (
size(work) < lwork)
then
1492 call errmgr%report_error(
"mtx_rank", &
1493 "Incorrectly sized input array WORK, argument 5.", &
1494 la_array_size_error)
1497 wptr => work(1:lwork)
1499 allocate(wrk(lwork), stat = istat)
1500 if (istat /= 0)
then
1502 call errmgr%report_error(
"mtx_rank", &
1503 "Insufficient memory available.", &
1504 la_out_of_memory_error)
1510 w => wptr(mn+1:lwork)
1513 call dgesvd(
'N',
'N', m, n, a, m, s, dummy, m, dummy, n, w, &
1516 allocate(
character(len = 256) :: errmsg)
1517 write(errmsg, 100) flag,
" superdiagonals could not " // &
1518 "converge to zero as part of the QR iteration process."
1519 call errmgr%report_warning(
"mtx_rank", errmsg, la_convergence_error)
1524 tref = max(m, n) * epsilon(t) * s(1)
1525 if (
present(tol))
then
1530 if (t < smlnum)
then
1551 module function mtx_rank_cmplx(a, tol, work, olwork, rwork, err) result(rnk)
1553 complex(real64),
intent(inout),
dimension(:,:) :: a
1554 real(real64),
intent(in),
optional :: tol
1555 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
1556 integer(int32),
intent(out),
optional :: olwork
1557 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
1558 class(errors),
intent(inout),
optional,
target :: err
1559 integer(int32) :: rnk
1563 function dlamch(cmach)
result(x)
1564 use,
intrinsic :: iso_fortran_env, only : real64
1565 character,
intent(in) :: cmach
1571 integer(int32) :: i, m, n, mn, istat, lwork, flag, lrwork
1572 real(real64),
pointer,
dimension(:) :: s, rwptr, rw
1573 real(real64),
allocatable,
target,
dimension(:) :: rwrk
1574 complex(real64),
allocatable,
target,
dimension(:) :: wrk
1575 complex(real64),
pointer,
dimension(:) :: wptr
1576 real(real64) :: t, tref, smlnum
1577 real(real64),
dimension(1) :: dummy
1578 complex(real64),
dimension(1) :: cdummy, temp
1579 class(errors),
pointer :: errmgr
1580 type(errors),
target :: deferr
1581 character(len = :),
allocatable :: errmsg
1588 smlnum = dlamch(
's')
1590 if (
present(err))
then
1597 call zgesvd(
'N',
'N', m, n, a, m, dummy, cdummy, m, cdummy, n, temp, &
1599 lwork = int(temp(1), int32)
1600 if (
present(olwork))
then
1606 if (
present(work))
then
1607 if (
size(work) < lwork)
then
1609 call errmgr%report_error(
"mtx_rank_cmplx", &
1610 "Incorrectly sized input array WORK, argument 5.", &
1611 la_array_size_error)
1614 wptr => work(1:lwork)
1616 allocate(wrk(lwork), stat = istat)
1617 if (istat /= 0)
then
1619 call errmgr%report_error(
"mtx_rank_cmplx", &
1620 "Insufficient memory available.", &
1621 la_out_of_memory_error)
1627 if (
present(rwork))
then
1628 if (
size(rwork) < lrwork)
then
1630 call errmgr%report_error(
"mtx_rank_cmplx", &
1631 "Incorrectly sized input array RWORK.", &
1632 la_array_size_error)
1635 rwptr => rwork(1:lrwork)
1637 allocate(rwrk(lrwork), stat = istat)
1638 if (istat /= 0)
then
1643 rw => rwptr(mn+1:lrwork)
1646 call zgesvd(
'N',
'N', m, n, a, m, s, cdummy, m, cdummy, n, wptr, &
1647 lwork - mn, rw, flag)
1649 allocate(
character(len = 256) :: errmsg)
1650 write(errmsg, 100) flag,
" superdiagonals could not " // &
1651 "converge to zero as part of the QR iteration process."
1652 call errmgr%report_warning(
"mtx_rank_cmplx", errmsg, la_convergence_error)
1657 tref = max(m, n) * epsilon(t) * s(1)
1658 if (
present(tol))
then
1663 if (t < smlnum)
then
1684 module function det_dbl(a, iwork, err) result(x)
1686 real(real64),
intent(inout),
dimension(:,:) :: a
1687 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1688 class(errors),
intent(inout),
optional,
target :: err
1692 real(real64),
parameter :: zero = 0.0d0
1693 real(real64),
parameter :: one = 1.0d0
1694 real(real64),
parameter :: ten = 1.0d1
1695 real(real64),
parameter :: p1 = 1.0d-1
1698 integer(int32) :: i, ep, n, istat, flag
1699 integer(int32),
pointer,
dimension(:) :: ipvt
1700 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1701 real(real64) :: temp
1702 class(errors),
pointer :: errmgr
1703 type(errors),
target :: deferr
1708 if (
present(err))
then
1715 if (
size(a, 2) /= n)
then
1716 call errmgr%report_error(
"det", &
1717 "The supplied matrix must be square.", la_array_size_error)
1722 if (
present(iwork))
then
1723 if (
size(iwork) < n)
then
1725 call errmgr%report_error(
"det", &
1726 "Incorrectly sized input array IWORK, argument 2.", &
1727 la_array_size_error)
1732 allocate(iwrk(n), stat = istat)
1733 if (istat /= 0)
then
1735 call errmgr%report_error(
"det", &
1736 "Insufficient memory available.", &
1737 la_out_of_memory_error)
1744 call dgetrf(n, n, a, n, ipvt, flag)
1755 if (ipvt(i) /= i) temp = -temp
1757 temp = a(i,i) * temp
1758 if (temp == zero)
then
1763 do while (abs(temp) < one)
1768 do while (abs(temp) > ten)
1777 module function det_cmplx(a, iwork, err) result(x)
1779 complex(real64),
intent(inout),
dimension(:,:) :: a
1780 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1781 class(errors),
intent(inout),
optional,
target :: err
1782 complex(real64) :: x
1785 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1786 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1787 complex(real64),
parameter :: ten = (1.0d1, 0.0d0)
1788 complex(real64),
parameter :: p1 = (1.0d-1, 0.0d0)
1789 real(real64),
parameter :: real_one = 1.0d0
1790 real(real64),
parameter :: real_ten = 1.0d1
1793 integer(int32) :: i, ep, n, istat, flag
1794 integer(int32),
pointer,
dimension(:) :: ipvt
1795 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1796 complex(real64) :: temp
1797 class(errors),
pointer :: errmgr
1798 type(errors),
target :: deferr
1803 if (
present(err))
then
1810 if (
size(a, 2) /= n)
then
1811 call errmgr%report_error(
"det_cmplx", &
1812 "The supplied matrix must be square.", la_array_size_error)
1817 if (
present(iwork))
then
1818 if (
size(iwork) < n)
then
1820 call errmgr%report_error(
"det_cmplx", &
1821 "Incorrectly sized input array IWORK, argument 2.", &
1822 la_array_size_error)
1827 allocate(iwrk(n), stat = istat)
1828 if (istat /= 0)
then
1830 call errmgr%report_error(
"det_cmplx", &
1831 "Insufficient memory available.", &
1832 la_out_of_memory_error)
1839 call zgetrf(n, n, a, n, ipvt, flag)
1850 if (ipvt(i) /= i) temp = -temp
1852 temp = a(i,i) * temp
1853 if (temp == zero)
then
1858 do while (abs(temp) < real_one)
1863 do while (abs(temp) > real_ten)
1874 module subroutine swap_dbl(x, y, err)
1876 real(real64),
intent(inout),
dimension(:) :: x, y
1877 class(errors),
intent(inout),
optional,
target :: err
1880 integer(int32) :: i, n
1881 real(real64) :: temp
1882 class(errors),
pointer :: errmgr
1883 type(errors),
target :: deferr
1887 if (
present(err))
then
1894 if (
size(y) /= n)
then
1895 call errmgr%report_error(
"swap", &
1896 "The input arrays are not the same size.", &
1897 la_array_size_error)
1910 module subroutine swap_cmplx(x, y, err)
1912 complex(real64),
intent(inout),
dimension(:) :: x, y
1913 class(errors),
intent(inout),
optional,
target :: err
1916 integer(int32) :: i, n
1917 complex(real64) :: temp
1918 class(errors),
pointer :: errmgr
1919 type(errors),
target :: deferr
1923 if (
present(err))
then
1930 if (
size(y) /= n)
then
1931 call errmgr%report_error(
"swap_cmplx", &
1932 "The input arrays are not the same size.", &
1933 la_array_size_error)
1948 module subroutine recip_mult_array_dbl(a, x)
1950 real(real64),
intent(in) :: a
1951 real(real64),
intent(inout),
dimension(:) :: x
1955 function dlamch(cmach)
result(x)
1956 use,
intrinsic :: iso_fortran_env, only : real64
1957 character,
intent(in) :: cmach
1963 real(real64),
parameter :: zero = 0.0d0
1964 real(real64),
parameter :: one = 1.0d0
1965 real(real64),
parameter :: twotho = 2.0d3
1969 real(real64) :: bignum, cden, cden1, cnum, cnum1, mul, smlnum
1972 smlnum = dlamch(
's')
1973 bignum = one / smlnum
1974 if (log10(bignum) > twotho)
then
1975 smlnum = sqrt(smlnum)
1976 bignum = sqrt(bignum)
1985 cden1 = cden * smlnum
1986 cnum1 = cnum / bignum
1987 if (abs(cden1) > abs(cnum) .and. cnum /= zero)
then
1991 else if (abs(cnum1) > abs(cden))
then
2011 module subroutine tri_mtx_mult_dbl(upper, alpha, a, beta, b, err)
2013 logical,
intent(in) :: upper
2014 real(real64),
intent(in) :: alpha, beta
2015 real(real64),
intent(in),
dimension(:,:) :: a
2016 real(real64),
intent(inout),
dimension(:,:) :: b
2017 class(errors),
intent(inout),
optional,
target :: err
2020 real(real64),
parameter :: zero = 0.0d0
2023 integer(int32) :: i, j, k, n, d1, d2, flag
2024 real(real64) :: temp
2025 class(errors),
pointer :: errmgr
2026 type(errors),
target :: deferr
2027 character(len = :),
allocatable :: errmsg
2033 if (
present(err))
then
2041 if (
size(a, 2) /= n)
then
2044 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
2051 allocate(
character(len = 256) :: errmsg)
2052 write(errmsg, 100)
"The matrix at input ", flag, &
2053 " was not sized appropriately. A matrix of ", n,
"-by-", n, &
2054 "was expected, but a matrix of ", d1,
"-by-", d2,
" was found."
2055 call errmgr%report_error(
"tri_mtx_mult_dbl", trim(errmsg), &
2056 la_array_size_error)
2063 if (beta == zero)
then
2068 temp = temp + a(k,i) * a(k,j)
2072 if (i /= j) b(j,i) = temp
2080 temp = temp + a(k,i) * a(k,j)
2083 b(i,j) = temp + beta * b(i,j)
2084 if (i /= j) b(j,i) = temp + beta * b(j,i)
2090 if (beta == zero)
then
2095 temp = temp + a(i,k) * a(j,k)
2099 if (i /= j) b(j,i) = temp
2107 temp = temp + a(i,k) * a(j,k)
2110 b(i,j) = temp + beta * b(i,j)
2111 if (i /= j) b(j,i) = temp + beta * b(j,i)
2118100
format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2122 module subroutine tri_mtx_mult_cmplx(upper, alpha, a, beta, b, err)
2124 logical,
intent(in) :: upper
2125 complex(real64),
intent(in) :: alpha, beta
2126 complex(real64),
intent(in),
dimension(:,:) :: a
2127 complex(real64),
intent(inout),
dimension(:,:) :: b
2128 class(errors),
intent(inout),
optional,
target :: err
2131 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
2134 integer(int32) :: i, j, k, n, d1, d2, flag
2135 complex(real64) :: temp
2136 class(errors),
pointer :: errmgr
2137 type(errors),
target :: deferr
2138 character(len = :),
allocatable :: errmsg
2144 if (
present(err))
then
2152 if (
size(a, 2) /= n)
then
2155 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
2162 allocate(
character(len = 256) :: errmsg)
2163 write(errmsg, 100)
"The matrix at input ", flag, &
2164 " was not sized appropriately. A matrix of ", n,
"-by-", n, &
2165 "was expected, but a matrix of ", d1,
"-by-", d2,
" was found."
2166 call errmgr%report_error(
"tri_mtx_mult_cmplx", trim(errmsg), &
2167 la_array_size_error)
2174 if (beta == zero)
then
2179 temp = temp + a(k,i) * a(k,j)
2183 if (i /= j) b(j,i) = temp
2191 temp = temp + a(k,i) * a(k,j)
2194 b(i,j) = temp + beta * b(i,j)
2195 if (i /= j) b(j,i) = temp + beta * b(j,i)
2201 if (beta == zero)
then
2206 temp = temp + a(i,k) * a(j,k)
2210 if (i /= j) b(j,i) = temp
2218 temp = temp + a(i,k) * a(j,k)
2221 b(i,j) = temp + beta * b(i,j)
2222 if (i /= j) b(j,i) = temp + beta * b(j,i)
2229100
format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
Provides a set of common linear algebra routines.