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 = 128) :: 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
67 "Matrix dimension mismatch. Input number ", flag, &
68 " was not sized correctly."
69 call errmgr%report_error(
"mtx_mult_mtx", errmsg, &
75 call dgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
82 module subroutine mtx_mult_vec(trans, alpha, a, b, beta, c, err)
84 logical,
intent(in) :: trans
85 real(real64),
intent(in) :: alpha, beta
86 real(real64),
intent(in),
dimension(:,:) :: a
87 real(real64),
intent(in),
dimension(:) :: b
88 real(real64),
intent(inout),
dimension(:) :: c
89 class(errors),
intent(inout),
optional,
target :: err
93 integer(int32) :: m, n, flag
94 class(errors),
pointer :: errmgr
95 type(errors),
target :: deferr
96 character(len = 128) :: errmsg
103 if (
present(err))
then
112 if (
size(b) /= m)
then
114 else if (
size(c) /= n)
then
118 if (
size(b) /= n)
then
120 else if (
size(c) /= m)
then
127 "Matrix dimension mismatch. Input number ", flag, &
128 " was not sized correctly."
129 call errmgr%report_error(
"mtx_mult_vec", errmsg, &
135 call dgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
144 module subroutine cmtx_mult_mtx(opa, opb, alpha, a, b, beta, c, err)
146 integer(int32),
intent(in) :: opa, opb
147 complex(real64),
intent(in) :: alpha, beta
148 complex(real64),
intent(in),
dimension(:,:) :: a, b
149 complex(real64),
intent(inout),
dimension(:,:) :: c
150 class(errors),
intent(inout),
optional,
target :: err
153 real(real64),
parameter :: zero = 0.0d0
154 real(real64),
parameter :: one = 1.0d0
158 integer(int32) :: m, n, k, lda, ldb, flag
159 class(errors),
pointer :: errmgr
160 type(errors),
target :: deferr
161 character(len = 128) :: errmsg
166 if (opa == la_transpose)
then
170 else if (opa == la_hermitian_transpose)
then
179 if (opb == la_transpose)
then
182 else if (opb == la_hermitian_transpose)
then
189 if (
present(err))
then
197 if (opa == la_transpose .or. opa == la_hermitian_transpose)
then
198 if (
size(a, 2) /= m) flag = 4
200 if (
size(a, 1) /= m) flag = 4
202 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
203 if (
size(b, 2) /= k .or.
size(b, 1) /= n) flag = 5
205 if (
size(b, 1) /= k .or.
size(b, 2) /= n) flag = 5
210 "Matrix dimension mismatch. Input number ", flag, &
211 " was not sized correctly."
212 call errmgr%report_error(
"cmtx_mult_mtx", errmsg, &
218 call zgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
225 module subroutine cmtx_mult_vec(opa, alpha, a, b, beta, c, err)
227 integer(int32),
intent(in) :: opa
228 complex(real64),
intent(in) :: alpha, beta
229 complex(real64),
intent(in),
dimension(:,:) :: a
230 complex(real64),
intent(in),
dimension(:) :: b
231 complex(real64),
intent(inout),
dimension(:) :: c
232 class(errors),
intent(inout),
optional,
target :: err
236 integer(int32) :: m, n, flag
237 class(errors),
pointer :: errmgr
238 type(errors),
target :: deferr
239 character(len = 128) :: errmsg
244 if (opa == la_transpose)
then
246 else if (opa == la_hermitian_transpose)
then
251 if (
present(err))
then
259 if (opa == la_transpose .or. opa == la_hermitian_transpose)
then
260 if (
size(b) /= m)
then
262 else if (
size(c) /= n)
then
266 if (
size(b) /= n)
then
268 else if (
size(c) /= m)
then
275 "Matrix dimension mismatch. Input number ", flag, &
276 " was not sized correctly."
277 call errmgr%report_error(
"cmtx_mult_vec", errmsg, &
283 call zgemv(t, m, n, alpha, a, m, b, 1, beta, c, 1)
292 module subroutine rank1_update_dbl(alpha, x, y, a, err)
294 real(real64),
intent(in) :: alpha
295 real(real64),
intent(in),
dimension(:) :: x, y
296 real(real64),
intent(inout),
dimension(:,:) :: a
297 class(errors),
intent(inout),
optional,
target :: err
300 real(real64),
parameter :: zero = 0.0d0
303 integer(int32) :: j, m, n
305 class(errors),
pointer :: errmgr
306 type(errors),
target :: deferr
311 if (
present(err))
then
318 if (
size(a, 1) /= m .or.
size(a, 2) /= n)
then
320 call errmgr%report_error(
"rank1_update_dbl", &
321 "Matrix dimension mismatch.", la_array_size_error)
327 if (y(j) /= zero)
then
329 a(:,j) = a(:,j) + temp * x
337 module subroutine rank1_update_cmplx(alpha, x, y, a, err)
339 complex(real64),
intent(in) :: alpha
340 complex(real64),
intent(in),
dimension(:) :: x, y
341 complex(real64),
intent(inout),
dimension(:,:) :: a
342 class(errors),
intent(inout),
optional,
target :: err
345 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
348 integer(int32) :: j, m, n
349 complex(real64) :: temp
350 class(errors),
pointer :: errmgr
351 type(errors),
target :: deferr
356 if (
present(err))
then
363 if (
size(a, 1) /= m .or.
size(a, 2) /= n)
then
365 call errmgr%report_error(
"rank1_update_cmplx", &
366 "Matrix dimension mismatch.", la_array_size_error)
372 if (y(j) /= zero)
then
373 temp = alpha * conjg(y(j))
374 a(:,j) = a(:,j) + temp * x
382 module subroutine diag_mtx_mult_mtx(lside, trans, alpha, a, b, beta, c, err)
384 logical,
intent(in) :: lside, trans
385 real(real64) :: alpha, beta
386 real(real64),
intent(in),
dimension(:) :: a
387 real(real64),
intent(in),
dimension(:,:) :: b
388 real(real64),
intent(inout),
dimension(:,:) :: c
389 class(errors),
intent(inout),
optional,
target :: err
392 real(real64),
parameter :: zero = 0.0d0
393 real(real64),
parameter :: one = 1.0d0
396 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
398 class(errors),
pointer :: errmgr
399 type(errors),
target :: deferr
400 character(len = 128) :: errmsg
408 if (
present(err))
then
422 if (nrowb /= n .or. ncolb < k) flag = 5
425 if (nrowb < k .or. ncolb /= n) flag = 5
434 if (ncolb /= m .or. nrowb < k) flag = 5
437 if (nrowb /= m .or. ncolb < k) flag = 5
443 write(errmsg, 100)
"Input number ", flag, &
444 " is not sized correctly."
445 call errmgr%report_error(
"diag_mtx_mult_mtx", trim(errmsg), &
452 if (beta == zero)
then
454 else if (beta /= one)
then
465 if (beta == zero)
then
467 else if (beta /= one)
then
468 c(i,:) = beta * c(i,:)
471 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
476 if (beta == zero)
then
478 else if (beta /= one)
then
479 c(i,:) = beta * c(i,:)
482 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
488 if (beta == zero)
then
491 c(k+1:m,:) = beta * c(k+1:m,:)
498 if (beta == zero)
then
500 else if (beta /= one)
then
501 c(:,i) = beta * c(:,i)
504 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
509 if (beta == zero)
then
511 else if (beta /= one)
then
512 c(:,i) = beta * c(:,i)
515 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
521 if (beta == zero)
then
523 else if (beta /= one)
then
524 c(:,k+1:m) = beta * c(:,k+1:m)
534 module subroutine diag_mtx_mult_mtx2(lside, alpha, a, b, err)
536 logical,
intent(in) :: lside
537 real(real64),
intent(in) :: alpha
538 real(real64),
intent(in),
dimension(:) :: a
539 real(real64),
intent(inout),
dimension(:,:) :: b
540 class(errors),
intent(inout),
optional,
target :: err
543 real(real64),
parameter :: zero = 0.0d0
544 real(real64),
parameter :: one = 1.0d0
547 integer(int32) :: i, m, n, k
549 class(errors),
pointer :: errmgr
550 type(errors),
target :: deferr
556 if (
present(err))
then
563 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
565 call errmgr%report_error(
"diag_mtx_mult_mtx2", &
566 "Input number 3 is not sized correctly.", &
576 if (temp /= one) b(i,:) = temp * b(i,:)
578 if (m > k) b(k+1:m,:) = zero
583 if (temp /= one) b(:,i) = temp * b(:,i)
585 if (n > k) b(:,k+1:n) = zero
590 module subroutine diag_mtx_mult_mtx3(lside, trans, alpha, a, b, beta, c, err)
592 logical,
intent(in) :: lside, trans
593 real(real64) :: alpha, beta
594 complex(real64),
intent(in),
dimension(:) :: a
595 real(real64),
intent(in),
dimension(:,:) :: b
596 complex(real64),
intent(inout),
dimension(:,:) :: c
597 class(errors),
intent(inout),
optional,
target :: err
600 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
601 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
604 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
605 complex(real64) :: temp
606 class(errors),
pointer :: errmgr
607 type(errors),
target :: deferr
608 character(len = 128) :: errmsg
616 if (
present(err))
then
630 if (nrowb /= n .or. ncolb < k) flag = 5
633 if (nrowb < k .or. ncolb /= n) flag = 5
642 if (ncolb /= m .or. nrowb < k) flag = 5
645 if (nrowb /= m .or. ncolb < k) flag = 5
651 write(errmsg, 100)
"Input number ", flag, &
652 " is not sized correctly."
653 call errmgr%report_error(
"diag_mtx_mult_mtx3", trim(errmsg), &
660 if (beta == zero)
then
662 else if (beta /= one)
then
673 if (beta == zero)
then
675 else if (beta /= one)
then
676 c(i,:) = beta * c(i,:)
679 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
684 if (beta == zero)
then
686 else if (beta /= one)
then
687 c(i,:) = beta * c(i,:)
690 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
696 if (beta == zero)
then
699 c(k+1:m,:) = beta * c(k+1:m,:)
706 if (beta == zero)
then
708 else if (beta /= one)
then
709 c(:,i) = beta * c(:,i)
712 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
717 if (beta == zero)
then
719 else if (beta /= one)
then
720 c(:,i) = beta * c(:,i)
723 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
729 if (beta == zero)
then
731 else if (beta /= one)
then
732 c(:,k+1:m) = beta * c(:,k+1:m)
742 module subroutine diag_mtx_mult_mtx4(lside, opb, alpha, a, b, beta, c, err)
744 logical,
intent(in) :: lside
745 integer(int32),
intent(in) :: opb
746 real(real64) :: alpha, beta
747 complex(real64),
intent(in),
dimension(:) :: a
748 complex(real64),
intent(in),
dimension(:,:) :: b
749 complex(real64),
intent(inout),
dimension(:,:) :: c
750 class(errors),
intent(inout),
optional,
target :: err
753 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
754 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
757 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
758 complex(real64) :: temp
759 class(errors),
pointer :: errmgr
760 type(errors),
target :: deferr
761 character(len = 128) :: errmsg
769 if (
present(err))
then
781 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
783 if (nrowb /= n .or. ncolb < k) flag = 5
786 if (nrowb < k .or. ncolb /= n) flag = 5
793 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
795 if (ncolb /= m .or. nrowb < k) flag = 5
798 if (nrowb /= m .or. ncolb < k) flag = 5
804 write(errmsg, 100)
"Input number ", flag, &
805 " is not sized correctly."
806 call errmgr%report_error(
"diag_mtx_mult_mtx4", trim(errmsg), &
813 if (beta == zero)
then
815 else if (beta /= one)
then
823 if (opb == la_transpose)
then
826 if (beta == zero)
then
828 else if (beta /= one)
then
829 c(i,:) = beta * c(i,:)
832 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
834 else if (opb == la_hermitian_transpose)
then
837 if (beta == zero)
then
839 else if (beta /= one)
then
840 c(i,:) = beta * c(i,:)
843 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
848 if (beta == zero)
then
850 else if (beta /= one)
then
851 c(i,:) = beta * c(i,:)
854 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
860 if (beta == zero)
then
863 c(k+1:m,:) = beta * c(k+1:m,:)
867 if (opb == la_transpose)
then
870 if (beta == zero)
then
872 else if (beta /= one)
then
873 c(:,i) = beta * c(:,i)
876 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
878 else if (opb == la_hermitian_transpose)
then
881 if (beta == zero)
then
883 else if (beta /= one)
then
884 c(:,i) = beta * c(:,i)
887 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
892 if (beta == zero)
then
894 else if (beta /= one)
then
895 c(:,i) = beta * c(:,i)
898 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
904 if (beta == zero)
then
906 else if (beta /= one)
then
907 c(:,k+1:m) = beta * c(:,k+1:m)
917 module subroutine diag_mtx_mult_mtx_cmplx(lside, opb, alpha, a, b, beta, c, err)
919 logical,
intent(in) :: lside
920 integer(int32),
intent(in) :: opb
921 complex(real64) :: alpha, beta
922 complex(real64),
intent(in),
dimension(:) :: a
923 complex(real64),
intent(in),
dimension(:,:) :: b
924 complex(real64),
intent(inout),
dimension(:,:) :: c
925 class(errors),
intent(inout),
optional,
target :: err
928 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
929 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
932 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
933 complex(real64) :: temp
934 class(errors),
pointer :: errmgr
935 type(errors),
target :: deferr
936 character(len = 128) :: errmsg
944 if (
present(err))
then
956 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
958 if (nrowb /= n .or. ncolb < k) flag = 5
961 if (nrowb < k .or. ncolb /= n) flag = 5
968 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
970 if (ncolb /= m .or. nrowb < k) flag = 5
973 if (nrowb /= m .or. ncolb < k) flag = 5
979 write(errmsg, 100)
"Input number ", flag, &
980 " is not sized correctly."
981 call errmgr%report_error(
"diag_mtx_mult_mtx_cmplx", trim(errmsg), &
988 if (beta == zero)
then
990 else if (beta /= one)
then
998 if (opb == la_transpose)
then
1001 if (beta == zero)
then
1003 else if (beta /= one)
then
1004 c(i,:) = beta * c(i,:)
1007 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
1009 else if (opb == la_hermitian_transpose)
then
1012 if (beta == zero)
then
1014 else if (beta /= one)
then
1015 c(i,:) = beta * c(i,:)
1018 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
1023 if (beta == zero)
then
1025 else if (beta /= one)
then
1026 c(i,:) = beta * c(i,:)
1029 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
1035 if (beta == zero)
then
1038 c(k+1:m,:) = beta * c(k+1:m,:)
1042 if (opb == la_transpose)
then
1045 if (beta == zero)
then
1047 else if (beta /= one)
then
1048 c(:,i) = beta * c(:,i)
1051 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
1053 else if (opb == la_hermitian_transpose)
then
1056 if (beta == zero)
then
1058 else if (beta /= one)
then
1059 c(:,i) = beta * c(:,i)
1062 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
1067 if (beta == zero)
then
1069 else if (beta /= one)
then
1070 c(:,i) = beta * c(:,i)
1073 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
1079 if (beta == zero)
then
1081 else if (beta /= one)
then
1082 c(:,k+1:m) = beta * c(:,k+1:m)
1092 module subroutine diag_mtx_mult_mtx2_cmplx(lside, alpha, a, b, err)
1094 logical,
intent(in) :: lside
1095 complex(real64),
intent(in) :: alpha
1096 complex(real64),
intent(in),
dimension(:) :: a
1097 complex(real64),
intent(inout),
dimension(:,:) :: b
1098 class(errors),
intent(inout),
optional,
target :: err
1101 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1102 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1105 integer(int32) :: i, m, n, k
1106 complex(real64) :: temp
1107 class(errors),
pointer :: errmgr
1108 type(errors),
target :: deferr
1114 if (
present(err))
then
1121 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
1123 call errmgr%report_error(
"diag_mtx_mult_mtx2_cmplx", &
1124 "Input number 3 is not sized correctly.", &
1125 la_array_size_error)
1134 if (temp /= one) b(i,:) = temp * b(i,:)
1136 if (m > k) b(k+1:m,:) = zero
1141 if (temp /= one) b(:,i) = temp * b(:,i)
1143 if (n > k) b(:,k+1:n) = zero
1148 module subroutine diag_mtx_mult_mtx_mix(lside, opb, alpha, a, b, beta, c, err)
1150 logical,
intent(in) :: lside
1151 integer(int32),
intent(in) :: opb
1152 complex(real64) :: alpha, beta
1153 real(real64),
intent(in),
dimension(:) :: a
1154 complex(real64),
intent(in),
dimension(:,:) :: b
1155 complex(real64),
intent(inout),
dimension(:,:) :: c
1156 class(errors),
intent(inout),
optional,
target :: err
1159 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1160 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1163 integer(int32) :: i, m, n, k, nrowb, ncolb, flag
1164 complex(real64) :: temp
1165 class(errors),
pointer :: errmgr
1166 type(errors),
target :: deferr
1167 character(len = 128) :: errmsg
1175 if (
present(err))
then
1187 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
1189 if (nrowb /= n .or. ncolb < k) flag = 5
1192 if (nrowb < k .or. ncolb /= n) flag = 5
1199 if (opb == la_transpose .or. opb == la_hermitian_transpose)
then
1201 if (ncolb /= m .or. nrowb < k) flag = 5
1204 if (nrowb /= m .or. ncolb < k) flag = 5
1210 write(errmsg, 100)
"Input number ", flag, &
1211 " is not sized correctly."
1212 call errmgr%report_error(
"diag_mtx_mult_mtx_mix", trim(errmsg), &
1213 la_array_size_error)
1218 if (alpha == 0)
then
1219 if (beta == zero)
then
1221 else if (beta /= one)
then
1229 if (opb == la_transpose)
then
1232 if (beta == zero)
then
1234 else if (beta /= one)
then
1235 c(i,:) = beta * c(i,:)
1238 if (temp /= one) c(i,:) = c(i,:) + temp * b(:,i)
1240 else if (opb == la_hermitian_transpose)
then
1243 if (beta == zero)
then
1245 else if (beta /= one)
then
1246 c(i,:) = beta * c(i,:)
1249 if (temp /= one) c(i,:) = c(i,:) + temp * conjg(b(:,i))
1254 if (beta == zero)
then
1256 else if (beta /= one)
then
1257 c(i,:) = beta * c(i,:)
1260 if (temp /= one) c(i,:) = c(i,:) + temp * b(i,:)
1266 if (beta == zero)
then
1269 c(k+1:m,:) = beta * c(k+1:m,:)
1273 if (opb == la_transpose)
then
1276 if (beta == zero)
then
1278 else if (beta /= one)
then
1279 c(:,i) = beta * c(:,i)
1282 if (temp /= one) c(:,i) = c(:,i) + temp * b(i,:)
1284 else if (opb == la_hermitian_transpose)
then
1287 if (beta == zero)
then
1289 else if (beta /= one)
then
1290 c(:,i) = beta * c(:,i)
1293 if (temp /= one) c(:,i) = c(:,i) + temp * conjg(b(i,:))
1298 if (beta == zero)
then
1300 else if (beta /= one)
then
1301 c(:,i) = beta * c(:,i)
1304 if (temp /= one) c(:,i) = c(:,i) + temp * b(:,i)
1310 if (beta == zero)
then
1312 else if (beta /= one)
then
1313 c(:,k+1:m) = beta * c(:,k+1:m)
1323 module subroutine diag_mtx_mult_mtx2_mix(lside, alpha, a, b, err)
1325 logical,
intent(in) :: lside
1326 complex(real64),
intent(in) :: alpha
1327 real(real64),
intent(in),
dimension(:) :: a
1328 complex(real64),
intent(inout),
dimension(:,:) :: b
1329 class(errors),
intent(inout),
optional,
target :: err
1332 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1333 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1336 integer(int32) :: i, m, n, k
1337 complex(real64) :: temp
1338 class(errors),
pointer :: errmgr
1339 type(errors),
target :: deferr
1345 if (
present(err))
then
1352 if ((lside .and. k > m) .or. (.not.lside .and. k > n))
then
1354 call errmgr%report_error(
"diag_mtx_mult_mtx2_cmplx", &
1355 "Input number 3 is not sized correctly.", &
1356 la_array_size_error)
1365 if (temp /= one) b(i,:) = temp * b(i,:)
1367 if (m > k) b(k+1:m,:) = zero
1372 if (temp /= one) b(:,i) = temp * b(:,i)
1374 if (n > k) b(:,k+1:n) = zero
1381 pure module function trace_dbl(x) result(y)
1383 real(real64),
intent(in),
dimension(:,:) :: x
1387 real(real64),
parameter :: zero = 0.0d0
1390 integer(int32) :: i, m, n, mn
1405 pure module function trace_cmplx(x) result(y)
1407 complex(real64),
intent(in),
dimension(:,:) :: x
1408 complex(real64) :: y
1411 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1414 integer(int32) :: i, m, n, mn
1429 module function mtx_rank_dbl(a, tol, work, olwork, err) result(rnk)
1431 real(real64),
intent(inout),
dimension(:,:) :: a
1432 real(real64),
intent(in),
optional :: tol
1433 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
1434 integer(int32),
intent(out),
optional :: olwork
1435 class(errors),
intent(inout),
optional,
target :: err
1436 integer(int32) :: rnk
1440 function dlamch(cmach)
result(x)
1441 use,
intrinsic :: iso_fortran_env, only : real64
1442 character,
intent(in) :: cmach
1448 integer(int32) :: i, m, n, mn, istat, lwork, flag
1449 real(real64),
pointer,
dimension(:) :: wptr, s, w
1450 real(real64),
allocatable,
target,
dimension(:) :: wrk
1451 real(real64) :: t, tref, smlnum
1452 real(real64),
dimension(1) :: dummy, temp
1453 class(errors),
pointer :: errmgr
1454 type(errors),
target :: deferr
1455 character(len = 128) :: errmsg
1461 smlnum = dlamch(
's')
1463 if (
present(err))
then
1471 call dgesvd(
'N',
'N', m, n, a, m, dummy, dummy, m, dummy, n, temp, &
1473 lwork = int(temp(1), int32) + mn
1474 if (
present(olwork))
then
1480 if (
present(work))
then
1481 if (
size(work) < lwork)
then
1483 call errmgr%report_error(
"mtx_rank", &
1484 "Incorrectly sized input array WORK, argument 5.", &
1485 la_array_size_error)
1488 wptr => work(1:lwork)
1490 allocate(wrk(lwork), stat = istat)
1491 if (istat /= 0)
then
1493 call errmgr%report_error(
"mtx_rank", &
1494 "Insufficient memory available.", &
1495 la_out_of_memory_error)
1501 w => wptr(mn+1:lwork)
1504 call dgesvd(
'N',
'N', m, n, a, m, s, dummy, m, dummy, n, w, &
1507 write(errmsg, 100) flag,
" superdiagonals could not " // &
1508 "converge to zero as part of the QR iteration process."
1509 call errmgr%report_warning(
"mtx_rank", errmsg, la_convergence_error)
1514 tref = max(m, n) * epsilon(t) * s(1)
1515 if (
present(tol))
then
1520 if (t < smlnum)
then
1541 module function mtx_rank_cmplx(a, tol, work, olwork, rwork, err) result(rnk)
1543 complex(real64),
intent(inout),
dimension(:,:) :: a
1544 real(real64),
intent(in),
optional :: tol
1545 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
1546 integer(int32),
intent(out),
optional :: olwork
1547 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
1548 class(errors),
intent(inout),
optional,
target :: err
1549 integer(int32) :: rnk
1553 function dlamch(cmach)
result(x)
1554 use,
intrinsic :: iso_fortran_env, only : real64
1555 character,
intent(in) :: cmach
1561 integer(int32) :: i, m, n, mn, istat, lwork, flag, lrwork
1562 real(real64),
pointer,
dimension(:) :: s, rwptr, rw
1563 real(real64),
allocatable,
target,
dimension(:) :: rwrk
1564 complex(real64),
allocatable,
target,
dimension(:) :: wrk
1565 complex(real64),
pointer,
dimension(:) :: wptr
1566 real(real64) :: t, tref, smlnum
1567 real(real64),
dimension(1) :: dummy
1568 complex(real64),
dimension(1) :: cdummy, temp
1569 class(errors),
pointer :: errmgr
1570 type(errors),
target :: deferr
1571 character(len = 128) :: errmsg
1578 smlnum = dlamch(
's')
1580 if (
present(err))
then
1587 call zgesvd(
'N',
'N', m, n, a, m, dummy, cdummy, m, cdummy, n, temp, &
1589 lwork = int(temp(1), int32)
1590 if (
present(olwork))
then
1596 if (
present(work))
then
1597 if (
size(work) < lwork)
then
1599 call errmgr%report_error(
"mtx_rank_cmplx", &
1600 "Incorrectly sized input array WORK, argument 5.", &
1601 la_array_size_error)
1604 wptr => work(1:lwork)
1606 allocate(wrk(lwork), stat = istat)
1607 if (istat /= 0)
then
1609 call errmgr%report_error(
"mtx_rank_cmplx", &
1610 "Insufficient memory available.", &
1611 la_out_of_memory_error)
1617 if (
present(rwork))
then
1618 if (
size(rwork) < lrwork)
then
1620 call errmgr%report_error(
"mtx_rank_cmplx", &
1621 "Incorrectly sized input array RWORK.", &
1622 la_array_size_error)
1625 rwptr => rwork(1:lrwork)
1627 allocate(rwrk(lrwork), stat = istat)
1628 if (istat /= 0)
then
1633 rw => rwptr(mn+1:lrwork)
1636 call zgesvd(
'N',
'N', m, n, a, m, s, cdummy, m, cdummy, n, wptr, &
1637 lwork - mn, rw, flag)
1639 write(errmsg, 100) flag,
" superdiagonals could not " // &
1640 "converge to zero as part of the QR iteration process."
1641 call errmgr%report_warning(
"mtx_rank_cmplx", errmsg, la_convergence_error)
1646 tref = max(m, n) * epsilon(t) * s(1)
1647 if (
present(tol))
then
1652 if (t < smlnum)
then
1673 module function det_dbl(a, iwork, err) result(x)
1675 real(real64),
intent(inout),
dimension(:,:) :: a
1676 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1677 class(errors),
intent(inout),
optional,
target :: err
1681 real(real64),
parameter :: zero = 0.0d0
1682 real(real64),
parameter :: one = 1.0d0
1683 real(real64),
parameter :: ten = 1.0d1
1684 real(real64),
parameter :: p1 = 1.0d-1
1687 integer(int32) :: i, ep, n, istat, flag
1688 integer(int32),
pointer,
dimension(:) :: ipvt
1689 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1690 real(real64) :: temp
1691 class(errors),
pointer :: errmgr
1692 type(errors),
target :: deferr
1697 if (
present(err))
then
1704 if (
size(a, 2) /= n)
then
1705 call errmgr%report_error(
"det", &
1706 "The supplied matrix must be square.", la_array_size_error)
1711 if (
present(iwork))
then
1712 if (
size(iwork) < n)
then
1714 call errmgr%report_error(
"det", &
1715 "Incorrectly sized input array IWORK, argument 2.", &
1716 la_array_size_error)
1721 allocate(iwrk(n), stat = istat)
1722 if (istat /= 0)
then
1724 call errmgr%report_error(
"det", &
1725 "Insufficient memory available.", &
1726 la_out_of_memory_error)
1733 call dgetrf(n, n, a, n, ipvt, flag)
1744 if (ipvt(i) /= i) temp = -temp
1746 temp = a(i,i) * temp
1747 if (temp == zero)
then
1752 do while (abs(temp) < one)
1757 do while (abs(temp) > ten)
1766 module function det_cmplx(a, iwork, err) result(x)
1768 complex(real64),
intent(inout),
dimension(:,:) :: a
1769 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1770 class(errors),
intent(inout),
optional,
target :: err
1771 complex(real64) :: x
1774 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1775 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1776 complex(real64),
parameter :: ten = (1.0d1, 0.0d0)
1777 complex(real64),
parameter :: p1 = (1.0d-1, 0.0d0)
1778 real(real64),
parameter :: real_one = 1.0d0
1779 real(real64),
parameter :: real_ten = 1.0d1
1782 integer(int32) :: i, ep, n, istat, flag
1783 integer(int32),
pointer,
dimension(:) :: ipvt
1784 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1785 complex(real64) :: temp
1786 class(errors),
pointer :: errmgr
1787 type(errors),
target :: deferr
1792 if (
present(err))
then
1799 if (
size(a, 2) /= n)
then
1800 call errmgr%report_error(
"det_cmplx", &
1801 "The supplied matrix must be square.", la_array_size_error)
1806 if (
present(iwork))
then
1807 if (
size(iwork) < n)
then
1809 call errmgr%report_error(
"det_cmplx", &
1810 "Incorrectly sized input array IWORK, argument 2.", &
1811 la_array_size_error)
1816 allocate(iwrk(n), stat = istat)
1817 if (istat /= 0)
then
1819 call errmgr%report_error(
"det_cmplx", &
1820 "Insufficient memory available.", &
1821 la_out_of_memory_error)
1828 call zgetrf(n, n, a, n, ipvt, flag)
1839 if (ipvt(i) /= i) temp = -temp
1841 temp = a(i,i) * temp
1842 if (temp == zero)
then
1847 do while (abs(temp) < real_one)
1852 do while (abs(temp) > real_ten)
1863 module subroutine swap_dbl(x, y, err)
1865 real(real64),
intent(inout),
dimension(:) :: x, y
1866 class(errors),
intent(inout),
optional,
target :: err
1869 integer(int32) :: i, n
1870 real(real64) :: temp
1871 class(errors),
pointer :: errmgr
1872 type(errors),
target :: deferr
1876 if (
present(err))
then
1883 if (
size(y) /= n)
then
1884 call errmgr%report_error(
"swap", &
1885 "The input arrays are not the same size.", &
1886 la_array_size_error)
1899 module subroutine swap_cmplx(x, y, err)
1901 complex(real64),
intent(inout),
dimension(:) :: x, y
1902 class(errors),
intent(inout),
optional,
target :: err
1905 integer(int32) :: i, n
1906 complex(real64) :: temp
1907 class(errors),
pointer :: errmgr
1908 type(errors),
target :: deferr
1912 if (
present(err))
then
1919 if (
size(y) /= n)
then
1920 call errmgr%report_error(
"swap_cmplx", &
1921 "The input arrays are not the same size.", &
1922 la_array_size_error)
1937 module subroutine recip_mult_array_dbl(a, x)
1939 real(real64),
intent(in) :: a
1940 real(real64),
intent(inout),
dimension(:) :: x
1944 function dlamch(cmach)
result(x)
1945 use,
intrinsic :: iso_fortran_env, only : real64
1946 character,
intent(in) :: cmach
1952 real(real64),
parameter :: zero = 0.0d0
1953 real(real64),
parameter :: one = 1.0d0
1954 real(real64),
parameter :: twotho = 2.0d3
1958 real(real64) :: bignum, cden, cden1, cnum, cnum1, mul, smlnum
1961 smlnum = dlamch(
's')
1962 bignum = one / smlnum
1963 if (log10(bignum) > twotho)
then
1964 smlnum = sqrt(smlnum)
1965 bignum = sqrt(bignum)
1974 cden1 = cden * smlnum
1975 cnum1 = cnum / bignum
1976 if (abs(cden1) > abs(cnum) .and. cnum /= zero)
then
1980 else if (abs(cnum1) > abs(cden))
then
2000 module subroutine tri_mtx_mult_dbl(upper, alpha, a, beta, b, err)
2002 logical,
intent(in) :: upper
2003 real(real64),
intent(in) :: alpha, beta
2004 real(real64),
intent(in),
dimension(:,:) :: a
2005 real(real64),
intent(inout),
dimension(:,:) :: b
2006 class(errors),
intent(inout),
optional,
target :: err
2009 real(real64),
parameter :: zero = 0.0d0
2012 integer(int32) :: i, j, k, n, d1, d2, flag
2013 real(real64) :: temp
2014 class(errors),
pointer :: errmgr
2015 type(errors),
target :: deferr
2016 character(len = 128) :: errmsg
2022 if (
present(err))
then
2030 if (
size(a, 2) /= n)
then
2033 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
2040 write(errmsg, 100)
"The matrix at input ", flag, &
2041 " was not sized appropriately. A matrix of ", n,
"-by-", n, &
2042 "was expected, but a matrix of ", d1,
"-by-", d2,
" was found."
2043 call errmgr%report_error(
"tri_mtx_mult_dbl", trim(errmsg), &
2044 la_array_size_error)
2051 if (beta == zero)
then
2056 temp = temp + a(k,i) * a(k,j)
2060 if (i /= j) b(j,i) = temp
2068 temp = temp + a(k,i) * a(k,j)
2071 b(i,j) = temp + beta * b(i,j)
2072 if (i /= j) b(j,i) = temp + beta * b(j,i)
2078 if (beta == zero)
then
2083 temp = temp + a(i,k) * a(j,k)
2087 if (i /= j) b(j,i) = temp
2095 temp = temp + a(i,k) * a(j,k)
2098 b(i,j) = temp + beta * b(i,j)
2099 if (i /= j) b(j,i) = temp + beta * b(j,i)
2106100
format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
2110 module subroutine tri_mtx_mult_cmplx(upper, alpha, a, beta, b, err)
2112 logical,
intent(in) :: upper
2113 complex(real64),
intent(in) :: alpha, beta
2114 complex(real64),
intent(in),
dimension(:,:) :: a
2115 complex(real64),
intent(inout),
dimension(:,:) :: b
2116 class(errors),
intent(inout),
optional,
target :: err
2119 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
2122 integer(int32) :: i, j, k, n, d1, d2, flag
2123 complex(real64) :: temp
2124 class(errors),
pointer :: errmgr
2125 type(errors),
target :: deferr
2126 character(len = 128) :: errmsg
2132 if (
present(err))
then
2140 if (
size(a, 2) /= n)
then
2143 else if (
size(b, 1) /= n .or.
size(b, 2) /= n)
then
2150 write(errmsg, 100)
"The matrix at input ", flag, &
2151 " was not sized appropriately. A matrix of ", n,
"-by-", n, &
2152 "was expected, but a matrix of ", d1,
"-by-", d2,
" was found."
2153 call errmgr%report_error(
"tri_mtx_mult_cmplx", trim(errmsg), &
2154 la_array_size_error)
2161 if (beta == zero)
then
2166 temp = temp + a(k,i) * a(k,j)
2170 if (i /= j) b(j,i) = temp
2178 temp = temp + a(k,i) * a(k,j)
2181 b(i,j) = temp + beta * b(i,j)
2182 if (i /= j) b(j,i) = temp + beta * b(j,i)
2188 if (beta == zero)
then
2193 temp = temp + a(i,k) * a(j,k)
2197 if (i /= j) b(j,i) = temp
2205 temp = temp + a(i,k) * a(j,k)
2208 b(i,j) = temp + beta * b(i,j)
2209 if (i /= j) b(j,i) = temp + beta * b(j,i)
2216100
format(a, i0, a, i0, a, i0, a, i0, a, i0, a)
Provides a set of common linear algebra routines.