7submodule(
linalg) linalg_solve
12 module subroutine solve_tri_mtx(lside, upper, trans, nounit, alpha, a, b, err)
14 logical,
intent(in) :: lside, upper, trans, nounit
15 real(real64),
intent(in) :: alpha
16 real(real64),
intent(in),
dimension(:,:) :: a
17 real(real64),
intent(inout),
dimension(:,:) :: b
18 class(errors),
intent(inout),
optional,
target :: err
21 character :: side, uplo, transa, diag
24 integer(int32) :: m, n, nrowa
25 class(errors),
pointer :: errmgr
26 type(errors),
target :: deferr
53 if (
present(err))
then
60 if (
size(a, 1) /= nrowa .or.
size(a, 2) /= nrowa)
then
62 call errmgr%report_error(
"solve_tri_mtx", &
63 "The input matrix must be square.", la_array_size_error)
68 call dtrsm(side, uplo, transa, diag, m, n, alpha, a, nrowa, b, m)
72 module subroutine solve_tri_mtx_cmplx(lside, upper, trans, nounit, alpha, a, b, err)
74 logical,
intent(in) :: lside, upper, trans, nounit
75 complex(real64),
intent(in) :: alpha
76 complex(real64),
intent(in),
dimension(:,:) :: a
77 complex(real64),
intent(inout),
dimension(:,:) :: b
78 class(errors),
intent(inout),
optional,
target :: err
81 character :: side, uplo, transa, diag
84 integer(int32) :: m, n, nrowa
85 class(errors),
pointer :: errmgr
86 type(errors),
target :: deferr
113 if (
present(err))
then
120 if (
size(a, 1) /= nrowa .or.
size(a, 2) /= nrowa)
then
122 call errmgr%report_error(
"solve_tri_mtx_cmplx", &
123 "The input matrix must be square.", la_array_size_error)
128 call ztrsm(side, uplo, transa, diag, m, n, alpha, a, nrowa, b, m)
132 module subroutine solve_tri_vec(upper, trans, nounit, a, x, err)
134 logical,
intent(in) :: upper, trans, nounit
135 real(real64),
intent(in),
dimension(:,:) :: a
136 real(real64),
intent(inout),
dimension(:) :: x
137 class(errors),
intent(inout),
optional,
target :: err
140 real(real64),
parameter :: zero = 0.0d0
143 character :: uplo, t, diag
145 class(errors),
pointer :: errmgr
146 type(errors),
target :: deferr
165 if (
present(err))
then
172 if (
size(a, 2) /= n)
then
174 call errmgr%report_error(
"solve_tri_vec", &
175 "The input matrix must be square.", la_array_size_error)
177 else if (
size(x) /= n)
then
179 call errmgr%report_error(
"solve_tri_vec", &
180 "The inner matrix dimensions must be equal.", &
186 call dtrsv(uplo, t, diag, n, a, n, x, 1)
190 module subroutine solve_tri_vec_cmplx(upper, trans, nounit, a, x, err)
192 logical,
intent(in) :: upper, trans, nounit
193 complex(real64),
intent(in),
dimension(:,:) :: a
194 complex(real64),
intent(inout),
dimension(:) :: x
195 class(errors),
intent(inout),
optional,
target :: err
198 real(real64),
parameter :: zero = 0.0d0
201 character :: uplo, t, diag
203 class(errors),
pointer :: errmgr
204 type(errors),
target :: deferr
223 if (
present(err))
then
230 if (
size(a, 2) /= n)
then
232 call errmgr%report_error(
"solve_tri_vec_cmplx", &
233 "The input matrix must be square.", la_array_size_error)
235 else if (
size(x) /= n)
then
237 call errmgr%report_error(
"solve_tri_vec_cmplx", &
238 "The inner matrix dimensions must be equal.", &
244 call ztrsv(uplo, t, diag, n, a, n, x, 1)
250 module subroutine solve_lu_mtx(a, ipvt, b, err)
252 real(real64),
intent(in),
dimension(:,:) :: a
253 integer(int32),
intent(in),
dimension(:) :: ipvt
254 real(real64),
intent(inout),
dimension(:,:) :: b
255 class(errors),
intent(inout),
optional,
target :: err
258 integer(int32) :: n, nrhs, flag
259 class(errors),
pointer :: errmgr
260 type(errors),
target :: deferr
261 character(len = 128) :: errmsg
266 if (
present(err))
then
274 if (
size(a, 2) /= n)
then
276 else if (
size(ipvt) /= n)
then
278 else if (
size(b, 1) /= n)
then
283 write(errmsg, 100)
"Input number ", flag, &
284 " is not sized correctly."
285 call errmgr%report_error(
"solve_lu_mtx", trim(errmsg), &
291 call dgetrs(
"N", n, nrhs, a, n, ipvt, b, n, flag)
298 module subroutine solve_lu_mtx_cmplx(a, ipvt, b, err)
300 complex(real64),
intent(in),
dimension(:,:) :: a
301 integer(int32),
intent(in),
dimension(:) :: ipvt
302 complex(real64),
intent(inout),
dimension(:,:) :: b
303 class(errors),
intent(inout),
optional,
target :: err
306 integer(int32) :: n, nrhs, flag
307 class(errors),
pointer :: errmgr
308 type(errors),
target :: deferr
309 character(len = 128) :: errmsg
314 if (
present(err))
then
322 if (
size(a, 2) /= n)
then
324 else if (
size(ipvt) /= n)
then
326 else if (
size(b, 1) /= n)
then
331 write(errmsg, 100)
"Input number ", flag, &
332 " is not sized correctly."
333 call errmgr%report_error(
"solve_lu_mtx_cmplx", trim(errmsg), &
339 call zgetrs(
"N", n, nrhs, a, n, ipvt, b, n, flag)
346 module subroutine solve_lu_vec(a, ipvt, b, err)
348 real(real64),
intent(in),
dimension(:,:) :: a
349 integer(int32),
intent(in),
dimension(:) :: ipvt
350 real(real64),
intent(inout),
dimension(:) :: b
351 class(errors),
intent(inout),
optional,
target :: err
354 integer(int32) :: n, flag
355 class(errors),
pointer :: errmgr
356 type(errors),
target :: deferr
357 character(len = 128) :: errmsg
361 if (
present(err))
then
369 if (
size(a, 2) /= n)
then
371 else if (
size(ipvt) /= n)
then
373 else if (
size(b) /= n)
then
378 write(errmsg, 100)
"Input number ", flag, &
379 " is not sized correctly."
380 call errmgr%report_error(
"solve_lu_vec", trim(errmsg), &
386 call dgetrs(
"N", n, 1, a, n, ipvt, b, n, flag)
393 module subroutine solve_lu_vec_cmplx(a, ipvt, b, err)
395 complex(real64),
intent(in),
dimension(:,:) :: a
396 integer(int32),
intent(in),
dimension(:) :: ipvt
397 complex(real64),
intent(inout),
dimension(:) :: b
398 class(errors),
intent(inout),
optional,
target :: err
401 integer(int32) :: n, flag
402 class(errors),
pointer :: errmgr
403 type(errors),
target :: deferr
404 character(len = 128) :: errmsg
408 if (
present(err))
then
416 if (
size(a, 2) /= n)
then
418 else if (
size(ipvt) /= n)
then
420 else if (
size(b) /= n)
then
425 write(errmsg, 100)
"Input number ", flag, &
426 " is not sized correctly."
427 call errmgr%report_error(
"solve_lu_vec_cmplx", trim(errmsg), &
433 call zgetrs(
"N", n, 1, a, n, ipvt, b, n, flag)
442 module subroutine solve_qr_no_pivot_mtx(a, tau, b, work, olwork, err)
444 real(real64),
intent(inout),
dimension(:,:) :: a, b
445 real(real64),
intent(in),
dimension(:) :: tau
446 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
447 integer(int32),
intent(out),
optional :: olwork
448 class(errors),
intent(inout),
optional,
target :: err
451 real(real64),
parameter :: one = 1.0d0
454 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
455 real(real64),
pointer,
dimension(:) :: wptr
456 real(real64),
allocatable,
target,
dimension(:) :: wrk
457 class(errors),
pointer :: errmgr
458 type(errors),
target :: deferr
459 character(len = 128) :: errmsg
466 if (
present(err))
then
476 else if (
size(tau) /= k)
then
478 else if (
size(b, 1) /= m)
then
483 write(errmsg, 100)
"Input number ", flag, &
484 " is not sized correctly."
485 call errmgr%report_error(
"solve_qr_no_pivot_mtx", trim(errmsg), &
491 call mult_qr(.true., .true., a, tau, b, olwork = lwork)
492 if (
present(olwork))
then
498 if (
present(work))
then
499 if (
size(work) < lwork)
then
501 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
502 "Incorrectly sized input array WORK, argument 4.", &
506 wptr => work(1:lwork)
508 allocate(wrk(lwork), stat = istat)
511 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
512 "Insufficient memory available.", &
513 la_out_of_memory_error)
520 call mult_qr(.true., .true., a, tau, b, wptr)
523 call solve_triangular_system(.true., .true., .false., .true., one, &
524 a(1:n,1:n), b(1:n,:))
531 module subroutine solve_qr_no_pivot_mtx_cmplx(a, tau, b, work, olwork, err)
533 complex(real64),
intent(inout),
dimension(:,:) :: a, b
534 complex(real64),
intent(in),
dimension(:) :: tau
535 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
536 integer(int32),
intent(out),
optional :: olwork
537 class(errors),
intent(inout),
optional,
target :: err
540 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
543 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
544 complex(real64),
pointer,
dimension(:) :: wptr
545 complex(real64),
allocatable,
target,
dimension(:) :: wrk
546 class(errors),
pointer :: errmgr
547 type(errors),
target :: deferr
548 character(len = 128) :: errmsg
555 if (
present(err))
then
565 else if (
size(tau) /= k)
then
567 else if (
size(b, 1) /= m)
then
572 write(errmsg, 100)
"Input number ", flag, &
573 " is not sized correctly."
574 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
575 trim(errmsg), la_array_size_error)
580 call mult_qr(.true., .true., a, tau, b, olwork = lwork)
581 if (
present(olwork))
then
587 if (
present(work))
then
588 if (
size(work) < lwork)
then
590 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
591 "Incorrectly sized input array WORK, argument 4.", &
595 wptr => work(1:lwork)
597 allocate(wrk(lwork), stat = istat)
600 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
601 "Insufficient memory available.", &
602 la_out_of_memory_error)
609 call mult_qr(.true., .true., a, tau, b, wptr)
612 call solve_triangular_system(.true., .true., .false., .true., one, &
613 a(1:n,1:n), b(1:n,:))
620 module subroutine solve_qr_no_pivot_vec(a, tau, b, work, olwork, err)
622 real(real64),
intent(inout),
dimension(:,:) :: a
623 real(real64),
intent(in),
dimension(:) :: tau
624 real(real64),
intent(inout),
dimension(:) :: b
625 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
626 integer(int32),
intent(out),
optional :: olwork
627 class(errors),
intent(inout),
optional,
target :: err
630 integer(int32) :: m, n, k, flag, lwork, istat
631 real(real64),
pointer,
dimension(:) :: wptr
632 real(real64),
allocatable,
target,
dimension(:) :: wrk
633 class(errors),
pointer :: errmgr
634 type(errors),
target :: deferr
635 character(len = 128) :: errmsg
641 if (
present(err))
then
651 else if (
size(tau) /= k)
then
653 else if (
size(b) /= m)
then
658 write(errmsg, 100)
"Input number ", flag, &
659 " is not sized correctly."
660 call errmgr%report_error(
"solve_qr_no_pivot_vec", trim(errmsg), &
666 call mult_qr(.true., a, tau, b, olwork = lwork)
667 if (
present(olwork))
then
673 if (
present(work))
then
674 if (
size(work) < lwork)
then
676 call errmgr%report_error(
"solve_qr_no_pivot_vec", &
677 "Incorrectly sized input array WORK, argument 4.", &
681 wptr => work(1:lwork)
683 allocate(wrk(lwork), stat = istat)
686 call errmgr%report_error(
"solve_qr_no_pivot_vec", &
687 "Insufficient memory available.", &
688 la_out_of_memory_error)
695 call mult_qr(.true., a, tau, b, work = wptr)
698 call solve_triangular_system(.true., .false., .true., a(1:n,1:n), b)
705 module subroutine solve_qr_no_pivot_vec_cmplx(a, tau, b, work, olwork, err)
707 complex(real64),
intent(inout),
dimension(:,:) :: a
708 complex(real64),
intent(in),
dimension(:) :: tau
709 complex(real64),
intent(inout),
dimension(:) :: b
710 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
711 integer(int32),
intent(out),
optional :: olwork
712 class(errors),
intent(inout),
optional,
target :: err
715 integer(int32) :: m, n, k, flag, lwork, istat
716 complex(real64),
pointer,
dimension(:) :: wptr
717 complex(real64),
allocatable,
target,
dimension(:) :: wrk
718 class(errors),
pointer :: errmgr
719 type(errors),
target :: deferr
720 character(len = 128) :: errmsg
726 if (
present(err))
then
736 else if (
size(tau) /= k)
then
738 else if (
size(b) /= m)
then
743 write(errmsg, 100)
"Input number ", flag, &
744 " is not sized correctly."
745 call errmgr%report_error(
"solve_qr_no_pivot_vec_cmplx", &
746 trim(errmsg), la_array_size_error)
751 call mult_qr(.true., a, tau, b, olwork = lwork)
752 if (
present(olwork))
then
758 if (
present(work))
then
759 if (
size(work) < lwork)
then
761 call errmgr%report_error(
"solve_qr_no_pivot_vec_cmplx", &
762 "Incorrectly sized input array WORK, argument 4.", &
766 wptr => work(1:lwork)
768 allocate(wrk(lwork), stat = istat)
771 call errmgr%report_error(
"solve_qr_no_pivot_vec_cmplx", &
772 "Insufficient memory available.", &
773 la_out_of_memory_error)
780 call mult_qr(.true., a, tau, b, work = wptr)
783 call solve_triangular_system(.true., .false., .true., a(1:n,1:n), b)
790 module subroutine solve_qr_pivot_mtx(a, tau, jpvt, b, work, olwork, err)
792 real(real64),
intent(inout),
dimension(:,:) :: a
793 real(real64),
intent(in),
dimension(:) :: tau
794 integer(int32),
intent(in),
dimension(:) :: jpvt
795 real(real64),
intent(inout),
dimension(:,:) :: b
796 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
797 integer(int32),
intent(out),
optional :: olwork
798 class(errors),
intent(inout),
optional,
target :: err
801 integer(int32),
parameter :: imin = 2
802 integer(int32),
parameter :: imax = 1
803 real(real64),
parameter :: zero = 0.0d0
804 real(real64),
parameter :: one = 1.0d0
807 integer(int32) :: i, j, m, n, mn, nrhs, lwork, ismin, ismax, &
808 rnk, maxmn, flag, istat, lwork1, lwork2, lwork3
809 real(real64) :: rcond, smax, smin, smaxpr, sminpr, s1, c1, s2, c2
810 real(real64),
pointer,
dimension(:) :: wptr, w, tau2
811 real(real64),
allocatable,
target,
dimension(:) :: wrk
812 class(errors),
pointer :: errmgr
813 type(errors),
target :: deferr
814 character(len = 128) :: errmsg
824 rcond = epsilon(rcond)
825 if (
present(err))
then
833 if (
size(tau) /= mn)
then
835 else if (
size(jpvt) /= n)
then
837 else if (
size(b, 1) /= maxmn)
then
842 write(errmsg, 100)
"Input number ", flag, &
843 " is not sized correctly."
844 call errmgr%report_error(
"solve_qr_pivot_mtx", trim(errmsg), &
850 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
851 call mult_qr(.true., .true., a, tau, b(1:m,:), olwork = lwork2)
852 call mult_rz(.true., .true., n, a(1:mn,:), a(1:mn,1), b(1:n,:), &
854 lwork = max(lwork1, lwork2, lwork3, 2 * mn + 1) + mn
855 if (
present(olwork))
then
861 if (
present(work))
then
862 if (
size(work) < lwork)
then
864 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
865 "Incorrectly sized input array WORK, argument 5.", &
869 wptr => work(1:lwork)
871 allocate(wrk(lwork), stat = istat)
874 call errmgr%report_error(
"solve_qr_pivot_mtx", &
875 "Insufficient memory available.", &
876 la_out_of_memory_error)
887 if (abs(a(1,1)) == zero)
then
897 call dlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
898 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
899 call dlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
900 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
901 if (smaxpr * rcond <= sminpr)
then
903 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
904 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
920 w => wptr(rnk+1:lwork)
921 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
924 call mult_qr(.true., .true., a, tau, b(1:m,:), w)
927 call solve_triangular_system(.true., .true., .false., .true., one, &
928 a(1:rnk,1:rnk), b(1:rnk,:))
929 if (n > rnk) b(rnk+1:n,:) = zero
933 call mult_rz(.true., .true., n - rnk, a(1:rnk,:), tau2, b(1:n,:), w)
939 wptr(jpvt(i)) = b(i,j)
949 module subroutine solve_qr_pivot_mtx_cmplx(a, tau, jpvt, b, work, olwork, err)
951 complex(real64),
intent(inout),
dimension(:,:) :: a
952 complex(real64),
intent(in),
dimension(:) :: tau
953 integer(int32),
intent(in),
dimension(:) :: jpvt
954 complex(real64),
intent(inout),
dimension(:,:) :: b
955 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
956 integer(int32),
intent(out),
optional :: olwork
957 class(errors),
intent(inout),
optional,
target :: err
960 integer(int32),
parameter :: imin = 2
961 integer(int32),
parameter :: imax = 1
962 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
963 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
966 integer(int32) :: i, j, m, n, mn, nrhs, lwork, ismin, ismax, &
967 rnk, maxmn, flag, istat, lwork1, lwork2, lwork3
968 real(real64) :: rcond, smax, smin, smaxpr, sminpr
969 complex(real64) :: s1, c1, s2, c2
970 complex(real64),
pointer,
dimension(:) :: wptr, w, tau2
971 complex(real64),
allocatable,
target,
dimension(:) :: wrk
972 class(errors),
pointer :: errmgr
973 type(errors),
target :: deferr
974 character(len = 128) :: errmsg
984 rcond = epsilon(rcond)
985 if (
present(err))
then
993 if (
size(tau) /= mn)
then
995 else if (
size(jpvt) /= n)
then
997 else if (
size(b, 1) /= maxmn)
then
1002 write(errmsg, 100)
"Input number ", flag, &
1003 " is not sized correctly."
1004 call errmgr%report_error(
"solve_qr_pivot_mtx_cmplx", &
1005 trim(errmsg), la_array_size_error)
1010 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1011 call mult_qr(.true., .true., a, tau, b(1:m,:), olwork = lwork2)
1012 call mult_rz(.true., .true., n, a(1:mn,:), a(1:mn,1), b(1:n,:), &
1014 lwork = max(lwork1, lwork2, lwork3, 2 * mn + 1) + mn
1015 if (
present(olwork))
then
1021 if (
present(work))
then
1022 if (
size(work) < lwork)
then
1024 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
1025 "Incorrectly sized input array WORK, argument 5.", &
1026 la_array_size_error)
1029 wptr => work(1:lwork)
1031 allocate(wrk(lwork), stat = istat)
1032 if (istat /= 0)
then
1034 call errmgr%report_error(
"solve_qr_pivot_mtx_cmplx", &
1035 "Insufficient memory available.", &
1036 la_out_of_memory_error)
1047 if (abs(a(1,1)) == zero)
then
1057 call zlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1058 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1059 call zlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1060 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1061 if (smaxpr * rcond <= sminpr)
then
1063 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1064 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1066 wptr(ismin+rnk) = c1
1067 wptr(ismax+rnk) = c2
1080 w => wptr(rnk+1:lwork)
1081 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
1084 call mult_qr(.true., .true., a, tau, b(1:m,:), w)
1087 call solve_triangular_system(.true., .true., .false., .true., one, &
1088 a(1:rnk,1:rnk), b(1:rnk,:))
1089 if (n > rnk) b(rnk+1:n,:) = zero
1093 call mult_rz(.true., .true., n - rnk, a(1:rnk,:), tau2, b(1:n,:), w)
1099 wptr(jpvt(i)) = b(i,j)
1109 module subroutine solve_qr_pivot_vec(a, tau, jpvt, b, work, olwork, err)
1111 real(real64),
intent(inout),
dimension(:,:) :: a
1112 real(real64),
intent(in),
dimension(:) :: tau
1113 integer(int32),
intent(in),
dimension(:) :: jpvt
1114 real(real64),
intent(inout),
dimension(:) :: b
1115 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
1116 integer(int32),
intent(out),
optional :: olwork
1117 class(errors),
intent(inout),
optional,
target :: err
1120 integer(int32),
parameter :: imin = 2
1121 integer(int32),
parameter :: imax = 1
1122 real(real64),
parameter :: zero = 0.0d0
1123 real(real64),
parameter :: one = 1.0d0
1126 integer(int32) :: i, m, n, mn, lwork, ismin, ismax, rnk, maxmn, flag, &
1127 istat, lwork1, lwork2
1128 real(real64) :: rcond, smax, smin, smaxpr, sminpr, s1, c1, s2, c2
1129 real(real64),
pointer,
dimension(:) :: wptr, w, tau2
1130 real(real64),
allocatable,
target,
dimension(:) :: wrk
1131 class(errors),
pointer :: errmgr
1132 type(errors),
target :: deferr
1133 character(len = 128) :: errmsg
1142 rcond = epsilon(rcond)
1143 if (
present(err))
then
1151 if (
size(tau) /= mn)
then
1153 else if (
size(jpvt) /= n)
then
1155 else if (
size(b) /= maxmn)
then
1160 write(errmsg, 100)
"Input number ", flag, &
1161 " is not sized correctly."
1162 call errmgr%report_error(
"solve_qr_pivot_vec", trim(errmsg), &
1163 la_array_size_error)
1168 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1169 call mult_rz(.true., n, a(1:mn,:), a(1:mn,1), b(1:n), olwork = lwork2)
1170 lwork = max(lwork1, lwork2, 2 * mn + 1) + mn
1171 if (
present(olwork))
then
1177 if (
present(work))
then
1178 if (
size(work) < lwork)
then
1180 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
1181 "Incorrectly sized input array WORK, argument 5.", &
1182 la_array_size_error)
1185 wptr => work(1:lwork)
1187 allocate(wrk(lwork), stat = istat)
1188 if (istat /= 0)
then
1190 call errmgr%report_error(
"solve_qr_pivot_vec", &
1191 "Insufficient memory available.", &
1192 la_out_of_memory_error)
1203 if (abs(a(1,1)) == zero)
then
1213 call dlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1214 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1215 call dlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1216 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1217 if (smaxpr * rcond <= sminpr)
then
1219 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1220 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1222 wptr(ismin+rnk) = c1
1223 wptr(ismax+rnk) = c2
1236 w => wptr(rnk+1:lwork)
1237 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
1240 call mult_qr(.true., a, tau, b(1:m))
1243 call solve_triangular_system(.true., .false., .true., a(1:rnk,1:rnk), &
1245 if (n > rnk) b(rnk+1:n) = zero
1249 call mult_rz(.true., n - rnk, a(1:rnk,:), tau2, b(1:n), w)
1254 wptr(jpvt(i)) = b(i)
1263 module subroutine solve_qr_pivot_vec_cmplx(a, tau, jpvt, b, work, olwork, err)
1265 complex(real64),
intent(inout),
dimension(:,:) :: a
1266 complex(real64),
intent(in),
dimension(:) :: tau
1267 integer(int32),
intent(in),
dimension(:) :: jpvt
1268 complex(real64),
intent(inout),
dimension(:) :: b
1269 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
1270 integer(int32),
intent(out),
optional :: olwork
1271 class(errors),
intent(inout),
optional,
target :: err
1274 integer(int32),
parameter :: imin = 2
1275 integer(int32),
parameter :: imax = 1
1276 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1277 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1280 integer(int32) :: i, m, n, mn, lwork, ismin, ismax, rnk, maxmn, flag, &
1281 istat, lwork1, lwork2
1282 real(real64) :: rcond, smax, smin, smaxpr, sminpr
1283 complex(real64) :: s1, c1, s2, c2
1284 complex(real64),
pointer,
dimension(:) :: wptr, w, tau2
1285 complex(real64),
allocatable,
target,
dimension(:) :: wrk
1286 class(errors),
pointer :: errmgr
1287 type(errors),
target :: deferr
1288 character(len = 128) :: errmsg
1297 rcond = epsilon(rcond)
1298 if (
present(err))
then
1306 if (
size(tau) /= mn)
then
1308 else if (
size(jpvt) /= n)
then
1310 else if (
size(b) /= maxmn)
then
1315 write(errmsg, 100)
"Input number ", flag, &
1316 " is not sized correctly."
1317 call errmgr%report_error(
"solve_qr_pivot_vec_cmplx", trim(errmsg), &
1318 la_array_size_error)
1323 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1324 call mult_rz(.true., n, a(1:mn,:), a(1:mn,1), b(1:n), olwork = lwork2)
1325 lwork = max(lwork1, lwork2, 2 * mn + 1) + mn
1326 if (
present(olwork))
then
1332 if (
present(work))
then
1333 if (
size(work) < lwork)
then
1335 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
1336 "Incorrectly sized input array WORK, argument 5.", &
1337 la_array_size_error)
1340 wptr => work(1:lwork)
1342 allocate(wrk(lwork), stat = istat)
1343 if (istat /= 0)
then
1345 call errmgr%report_error(
"solve_qr_pivot_vec_cmplx", &
1346 "Insufficient memory available.", &
1347 la_out_of_memory_error)
1358 if (abs(a(1,1)) == zero)
then
1368 call zlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1369 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1370 call zlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1371 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1372 if (smaxpr * rcond <= sminpr)
then
1374 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1375 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1377 wptr(ismin+rnk) = c1
1378 wptr(ismax+rnk) = c2
1391 w => wptr(rnk+1:lwork)
1392 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
1395 call mult_qr(.true., a, tau, b(1:m))
1398 call solve_triangular_system(.true., .false., .true., a(1:rnk,1:rnk), &
1400 if (n > rnk) b(rnk+1:n) = zero
1404 call mult_rz(.true., n - rnk, a(1:rnk,:), tau2, b(1:n), w)
1409 wptr(jpvt(i)) = b(i)
1420 module subroutine solve_cholesky_mtx(upper, a, b, err)
1422 logical,
intent(in) :: upper
1423 real(real64),
intent(in),
dimension(:,:) :: a
1424 real(real64),
intent(inout),
dimension(:,:) :: b
1425 class(errors),
intent(inout),
optional,
target :: err
1429 integer(int32) :: n, nrhs, flag
1430 class(errors),
pointer :: errmgr
1431 type(errors),
target :: deferr
1432 character(len = 128) :: errmsg
1442 if (
present(err))
then
1450 if (
size(a, 2) /= n)
then
1452 else if (
size(b, 1) /= n)
then
1456 write(errmsg, 100)
"Input number ", flag, &
1457 " is not sized correctly."
1458 call errmgr%report_error(
"solve_cholesky_mtx", trim(errmsg), &
1459 la_array_size_error)
1464 call dpotrs(uplo, n, nrhs, a, n, b, n, flag)
1471 module subroutine solve_cholesky_mtx_cmplx(upper, a, b, err)
1473 logical,
intent(in) :: upper
1474 complex(real64),
intent(in),
dimension(:,:) :: a
1475 complex(real64),
intent(inout),
dimension(:,:) :: b
1476 class(errors),
intent(inout),
optional,
target :: err
1480 integer(int32) :: n, nrhs, flag
1481 class(errors),
pointer :: errmgr
1482 type(errors),
target :: deferr
1483 character(len = 128) :: errmsg
1493 if (
present(err))
then
1501 if (
size(a, 2) /= n)
then
1503 else if (
size(b, 1) /= n)
then
1507 write(errmsg, 100)
"Input number ", flag, &
1508 " is not sized correctly."
1509 call errmgr%report_error(
"solve_cholesky_mtx_cmplx", trim(errmsg), &
1510 la_array_size_error)
1515 call zpotrs(uplo, n, nrhs, a, n, b, n, flag)
1522 module subroutine solve_cholesky_vec(upper, a, b, err)
1524 logical,
intent(in) :: upper
1525 real(real64),
intent(in),
dimension(:,:) :: a
1526 real(real64),
intent(inout),
dimension(:) :: b
1527 class(errors),
intent(inout),
optional,
target :: err
1531 integer(int32) :: n, flag
1532 class(errors),
pointer :: errmgr
1533 type(errors),
target :: deferr
1534 character(len = 128) :: errmsg
1543 if (
present(err))
then
1551 if (
size(a, 2) /= n)
then
1553 else if (
size(b) /= n)
then
1557 write(errmsg, 100)
"Input number ", flag, &
1558 " is not sized correctly."
1559 call errmgr%report_error(
"solve_cholesky_vec", trim(errmsg), &
1560 la_array_size_error)
1565 call dpotrs(uplo, n, 1, a, n, b, n, flag)
1572 module subroutine solve_cholesky_vec_cmplx(upper, a, b, err)
1574 logical,
intent(in) :: upper
1575 complex(real64),
intent(in),
dimension(:,:) :: a
1576 complex(real64),
intent(inout),
dimension(:) :: b
1577 class(errors),
intent(inout),
optional,
target :: err
1581 integer(int32) :: n, flag
1582 class(errors),
pointer :: errmgr
1583 type(errors),
target :: deferr
1584 character(len = 128) :: errmsg
1593 if (
present(err))
then
1601 if (
size(a, 2) /= n)
then
1603 else if (
size(b) /= n)
then
1607 write(errmsg, 100)
"Input number ", flag, &
1608 " is not sized correctly."
1609 call errmgr%report_error(
"solve_cholesky_vec_cmplx", trim(errmsg), &
1610 la_array_size_error)
1615 call zpotrs(uplo, n, 1, a, n, b, n, flag)
1624 module subroutine mtx_inverse_dbl(a, iwork, work, olwork, err)
1626 real(real64),
intent(inout),
dimension(:,:) :: a
1627 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1628 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
1629 integer(int32),
intent(out),
optional :: olwork
1630 class(errors),
intent(inout),
optional,
target :: err
1633 integer(int32) :: n, liwork, lwork, istat, flag, itemp(1)
1634 integer(int32),
pointer,
dimension(:) :: iptr
1635 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1636 real(real64),
pointer,
dimension(:) :: wptr
1637 real(real64),
allocatable,
target,
dimension(:) :: wrk
1638 real(real64),
dimension(1) :: temp
1639 class(errors),
pointer :: errmgr
1640 type(errors),
target :: deferr
1645 if (
present(err))
then
1652 if (
size(a, 2) /= n)
then
1653 call errmgr%report_error(
"mtx_inverse", &
1654 "The matrix must be squre to invert.", la_array_size_error)
1659 call dgetri(n, a, n, itemp, temp, -1, flag)
1660 lwork = int(temp(1), int32)
1661 if (
present(olwork))
then
1667 if (
present(work))
then
1668 if (
size(work) < lwork)
then
1670 call errmgr%report_error(
"mtx_inverse_dbl", &
1671 "Incorrectly sized input array WORK, argument 3.", &
1672 la_array_size_error)
1675 wptr => work(1:lwork)
1677 allocate(wrk(lwork), stat = istat)
1678 if (istat /= 0)
then
1680 call errmgr%report_error(
"mtx_inverse_dbl", &
1681 "Insufficient memory available.", &
1682 la_out_of_memory_error)
1689 if (
present(iwork))
then
1690 if (
size(iwork) < liwork)
then
1692 call errmgr%report_error(
"mtx_inverse_dbl", &
1693 "Incorrectly sized input array IWORK, argument 2.", &
1694 la_array_size_error)
1697 iptr => iwork(1:liwork)
1699 allocate(iwrk(liwork), stat = istat)
1700 if (istat /= 0)
then
1702 call errmgr%report_error(
"mtx_inverse_dbl", &
1703 "Insufficient memory available.", &
1704 la_out_of_memory_error)
1711 call dgetrf(n, n, a, n, iptr, flag)
1714 call dgetri(n, a, n, iptr, wptr, lwork, flag)
1718 call errmgr%report_error(
"mtx_inverse_dbl", &
1719 "The matrix is singular; therefore, the inverse could " // &
1720 "not be computed.", la_singular_matrix_error)
1725 module subroutine mtx_inverse_cmplx(a, iwork, work, olwork, err)
1727 complex(real64),
intent(inout),
dimension(:,:) :: a
1728 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1729 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
1730 integer(int32),
intent(out),
optional :: olwork
1731 class(errors),
intent(inout),
optional,
target :: err
1734 integer(int32) :: n, liwork, lwork, istat, flag, itemp(1)
1735 integer(int32),
pointer,
dimension(:) :: iptr
1736 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1737 complex(real64),
pointer,
dimension(:) :: wptr
1738 complex(real64),
allocatable,
target,
dimension(:) :: wrk
1739 complex(real64),
dimension(1) :: temp
1740 class(errors),
pointer :: errmgr
1741 type(errors),
target :: deferr
1746 if (
present(err))
then
1753 if (
size(a, 2) /= n)
then
1754 call errmgr%report_error(
"mtx_inverse_cmplx", &
1755 "The matrix must be squre to invert.", la_array_size_error)
1760 call zgetri(n, a, n, itemp, temp, -1, flag)
1761 lwork = int(temp(1), int32)
1762 if (
present(olwork))
then
1768 if (
present(work))
then
1769 if (
size(work) < lwork)
then
1771 call errmgr%report_error(
"mtx_inverse_cmplx", &
1772 "Incorrectly sized input array WORK, argument 3.", &
1773 la_array_size_error)
1776 wptr => work(1:lwork)
1778 allocate(wrk(lwork), stat = istat)
1779 if (istat /= 0)
then
1781 call errmgr%report_error(
"mtx_inverse_cmplx", &
1782 "Insufficient memory available.", &
1783 la_out_of_memory_error)
1790 if (
present(iwork))
then
1791 if (
size(iwork) < liwork)
then
1793 call errmgr%report_error(
"mtx_inverse_cmplx", &
1794 "Incorrectly sized input array IWORK, argument 2.", &
1795 la_array_size_error)
1798 iptr => iwork(1:liwork)
1800 allocate(iwrk(liwork), stat = istat)
1801 if (istat /= 0)
then
1803 call errmgr%report_error(
"mtx_inverse_cmplx", &
1804 "Insufficient memory available.", &
1805 la_out_of_memory_error)
1812 call zgetrf(n, n, a, n, iptr, flag)
1815 call zgetri(n, a, n, iptr, wptr, lwork, flag)
1819 call errmgr%report_error(
"mtx_inverse_cmplx", &
1820 "The matrix is singular; therefore, the inverse could " // &
1821 "not be computed.", la_singular_matrix_error)
1826 module subroutine mtx_pinverse_dbl(a, ainv, tol, work, olwork, err)
1828 real(real64),
intent(inout),
dimension(:,:) :: a
1829 real(real64),
intent(out),
dimension(:,:) :: ainv
1830 real(real64),
intent(in),
optional :: tol
1831 real(real64),
intent(out),
target,
dimension(:),
optional :: work
1832 integer(int32),
intent(out),
optional :: olwork
1833 class(errors),
intent(inout),
optional,
target :: err
1837 function dlamch(cmach)
result(x)
1838 use,
intrinsic :: iso_fortran_env, only : real64
1839 character,
intent(in) :: cmach
1845 real(real64),
parameter :: zero = 0.0d0
1846 real(real64),
parameter :: one = 1.0d0
1849 integer(int32) :: i, m, n, mn, lwork, istat, flag, i1, i2a, i2b, i3a, &
1851 real(real64),
pointer,
dimension(:) :: s, wptr, w
1852 real(real64),
pointer,
dimension(:,:) :: u, vt
1853 real(real64),
allocatable,
target,
dimension(:) :: wrk
1854 real(real64),
dimension(1) :: temp
1855 real(real64) :: t, tref, tolcheck
1856 class(errors),
pointer :: errmgr
1857 type(errors),
target :: deferr
1858 character(len = 128) :: errmsg
1866 i2b = i2a + n * n - 1
1870 tolcheck = dlamch(
's')
1871 if (
present(err))
then
1878 if (
size(ainv, 1) /= n .or.
size(ainv, 2) /= m)
then
1879 write(errmsg, 100) &
1880 "The output matrix AINV is not sized appropriately. " // &
1881 "It is expected to be ", n,
"-by-", m,
"."
1882 call errmgr%report_error(
"mtx_pinverse", errmsg, &
1883 la_array_size_error)
1888 call dgesvd(
'S',
'A', m, n, a, m, a(1:mn,:), a, m, a, n, temp, -1, flag)
1889 lwork = int(temp(1), int32)
1890 lwork = lwork + m * mn + n * n + mn
1891 if (
present(olwork))
then
1897 if (
present(work))
then
1898 if (
size(work) < lwork)
then
1900 call errmgr%report_error(
"mtx_pinverse", &
1901 "Incorrectly sized input array WORK, argument 4.", &
1902 la_array_size_error)
1905 wptr => work(1:lwork)
1907 allocate(wrk(lwork), stat = istat)
1908 if (istat /= 0)
then
1910 call errmgr%report_error(
"mtx_pinverse", &
1911 "Insufficient memory available.", &
1912 la_out_of_memory_error)
1917 u(1:m,1:mn) => wptr(1:i1)
1918 vt(1:n,1:n) => wptr(i2a:i2b)
1923 call dgesvd(
'S',
'A', m, n, a, m, s, u, m, vt, n, w,
size(w), flag)
1927 write(errmsg, 101) flag,
" superdiagonals could not " // &
1928 "converge to zero as part of the QR iteration process."
1929 call errmgr%report_warning(
"mtx_pinverse", errmsg, &
1930 la_convergence_error)
1936 tref = max(m, n) * epsilon(t) * s(1)
1937 if (
present(tol))
then
1943 if (t < tolcheck)
then
1961 call recip_mult_array(s(i), vt(i,1:n))
1966 call mtx_mult(.true., .true., one, vt(1:mn,:), u, zero, ainv)
1969100
format(a, i0, a, i0, a)
1974 module subroutine mtx_pinverse_cmplx(a, ainv, tol, work, olwork, rwork, err)
1976 complex(real64),
intent(inout),
dimension(:,:) :: a
1977 complex(real64),
intent(out),
dimension(:,:) :: ainv
1978 real(real64),
intent(in),
optional :: tol
1979 complex(real64),
intent(out),
target,
dimension(:),
optional :: work
1980 integer(int32),
intent(out),
optional :: olwork
1981 real(real64),
intent(out),
target,
dimension(:),
optional :: rwork
1982 class(errors),
intent(inout),
optional,
target :: err
1986 function dlamch(cmach)
result(x)
1987 use,
intrinsic :: iso_fortran_env, only : real64
1988 character,
intent(in) :: cmach
1994 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1995 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1998 integer(int32) :: i, m, n, mn, lwork, istat, flag, i1, i2a, i2b, i3, &
2000 real(real64),
pointer,
dimension(:) :: s, rwptr, rw
2001 real(real64),
allocatable,
target,
dimension(:) :: rwrk
2002 complex(real64),
pointer,
dimension(:) :: wptr, w
2003 complex(real64),
pointer,
dimension(:,:) :: u, vt
2004 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2005 complex(real64) :: temp(1), val
2006 real(real64) :: t, tref, tolcheck, rtemp(1)
2007 class(errors),
pointer :: errmgr
2008 type(errors),
target :: deferr
2009 character(len = 128) :: errmsg
2018 i2b = i2a + n * n - 1
2020 tolcheck = dlamch(
's')
2021 if (
present(err))
then
2028 if (
size(ainv, 1) /= n .or.
size(ainv, 2) /= m)
then
2029 write(errmsg, 100) &
2030 "The output matrix AINV is not sized appropriately. " // &
2031 "It is expected to be ", n,
"-by-", m,
"."
2032 call errmgr%report_error(
"mtx_pinverse_cmplx", errmsg, &
2033 la_array_size_error)
2038 call zgesvd(
'S',
'A', m, n, a, m, rtemp, a, m, a, n, temp, -1, &
2040 lwork = int(temp(1), int32)
2041 lwork = lwork + m * mn + n * n
2042 if (
present(olwork))
then
2048 if (
present(work))
then
2049 if (
size(work) < lwork)
then
2051 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2052 "Incorrectly sized input array WORK, argument 4.", &
2053 la_array_size_error)
2056 wptr => work(1:lwork)
2058 allocate(wrk(lwork), stat = istat)
2059 if (istat /= 0)
then
2061 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2062 "Insufficient memory available.", &
2063 la_out_of_memory_error)
2069 if (
present(rwork))
then
2070 if (
size(rwork) < lrwork)
then
2072 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2073 "Incorrectly sized input array RWORK, argument 6.", &
2074 la_array_size_error)
2077 rwptr => rwork(1:lrwork)
2079 allocate(rwrk(lrwork), stat = istat)
2080 if (istat /= 0)
then
2082 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2083 "Insufficient memory available.", &
2084 la_out_of_memory_error)
2089 u(1:m,1:mn) => wptr(1:i1)
2090 vt(1:n,1:n) => wptr(i2a:i2b)
2093 rw => rwptr(mn+1:lrwork)
2096 call zgesvd(
'S',
'A', m, n, a, m, s, u, m, vt, n, w,
size(w), rw, flag)
2100 write(errmsg, 101) flag,
" superdiagonals could not " // &
2101 "converge to zero as part of the QR iteration process."
2102 call errmgr%report_warning(
"mtx_pinverse_cmplx", errmsg, &
2103 la_convergence_error)
2109 tref = max(m, n) * epsilon(t) * s(1)
2110 if (
present(tol))
then
2116 if (t < tolcheck)
then
2135 vt(i,1:n) = conjg(vt(i,1:n)) / s(i)
2148 val = val + vt(k,i) * conjg(u(j,k))
2155100
format(a, i0, a, i0, a)
2162 module subroutine solve_least_squares_mtx(a, b, work, olwork, err)
2164 real(real64),
intent(inout),
dimension(:,:) :: a, b
2165 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2166 integer(int32),
intent(out),
optional :: olwork
2167 class(errors),
intent(inout),
optional,
target :: err
2170 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag
2171 real(real64),
pointer,
dimension(:) :: wptr
2172 real(real64),
allocatable,
target,
dimension(:) :: wrk
2173 real(real64),
dimension(1) :: temp
2174 class(errors),
pointer :: errmgr
2175 type(errors),
target :: deferr
2182 if (
present(err))
then
2189 if (
size(b, 1) /= maxmn)
then
2190 call errmgr%report_error(
"solve_least_squares_mtx", &
2191 "Input 2 is not sized correctly.", la_array_size_error)
2196 call dgels(
'N', m, n, nrhs, a, m, b, maxmn, temp, -1, flag)
2197 lwork = int(temp(1), int32)
2198 if (
present(olwork))
then
2204 if (
present(work))
then
2205 if (
size(work) < lwork)
then
2207 call errmgr%report_error(
"solve_least_squares_mtx", &
2208 "Incorrectly sized input array WORK, argument 3.", &
2209 la_array_size_error)
2212 wptr => work(1:lwork)
2214 allocate(wrk(lwork), stat = istat)
2215 if (istat /= 0)
then
2217 call errmgr%report_error(
"solve_least_squares_mtx", &
2218 "Insufficient memory available.", &
2219 la_out_of_memory_error)
2226 call dgels(
'N', m, n, nrhs, a, m, b, maxmn, wptr, lwork, flag)
2228 call errmgr%report_error(
"solve_least_squares_mtx", &
2229 "The supplied matrix is not of full rank; therefore, " // &
2230 "the solution could not be computed via this routine. " // &
2231 "Try a routine that utilizes column pivoting.", &
2232 la_invalid_operation_error)
2237 module subroutine solve_least_squares_mtx_cmplx(a, b, work, olwork, err)
2239 complex(real64),
intent(inout),
dimension(:,:) :: a, b
2240 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2241 integer(int32),
intent(out),
optional :: olwork
2242 class(errors),
intent(inout),
optional,
target :: err
2245 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag
2246 complex(real64),
pointer,
dimension(:) :: wptr
2247 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2248 complex(real64),
dimension(1) :: temp
2249 class(errors),
pointer :: errmgr
2250 type(errors),
target :: deferr
2257 if (
present(err))
then
2264 if (
size(b, 1) /= maxmn)
then
2265 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2266 "Input 2 is not sized correctly.", la_array_size_error)
2271 call zgels(
'N', m, n, nrhs, a, m, b, maxmn, temp, -1, flag)
2272 lwork = int(temp(1), int32)
2273 if (
present(olwork))
then
2279 if (
present(work))
then
2280 if (
size(work) < lwork)
then
2282 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2283 "Incorrectly sized input array WORK, argument 3.", &
2284 la_array_size_error)
2287 wptr => work(1:lwork)
2289 allocate(wrk(lwork), stat = istat)
2290 if (istat /= 0)
then
2292 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2293 "Insufficient memory available.", &
2294 la_out_of_memory_error)
2301 call zgels(
'N', m, n, nrhs, a, m, b, maxmn, wptr, lwork, flag)
2303 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2304 "The supplied matrix is not of full rank; therefore, " // &
2305 "the solution could not be computed via this routine. " // &
2306 "Try a routine that utilizes column pivoting.", &
2307 la_invalid_operation_error)
2312 module subroutine solve_least_squares_vec(a, b, work, olwork, err)
2314 real(real64),
intent(inout),
dimension(:,:) :: a
2315 real(real64),
intent(inout),
dimension(:) :: b
2316 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2317 integer(int32),
intent(out),
optional :: olwork
2318 class(errors),
intent(inout),
optional,
target :: err
2321 integer(int32) :: m, n, maxmn, lwork, istat, flag
2322 real(real64),
pointer,
dimension(:) :: wptr
2323 real(real64),
allocatable,
target,
dimension(:) :: wrk
2324 real(real64),
dimension(1) :: temp
2325 class(errors),
pointer :: errmgr
2326 type(errors),
target :: deferr
2332 if (
present(err))
then
2339 if (
size(b) /= maxmn)
then
2340 call errmgr%report_error(
"solve_least_squares_vec", &
2341 "Input 2 is not sized correctly.", la_array_size_error)
2346 call dgels(
'N', m, n, 1, a, m, b, maxmn, temp, -1, flag)
2347 lwork = int(temp(1), int32)
2348 if (
present(olwork))
then
2354 if (
present(work))
then
2355 if (
size(work) < lwork)
then
2357 call errmgr%report_error(
"solve_least_squares_vec", &
2358 "Incorrectly sized input array WORK, argument 3.", &
2359 la_array_size_error)
2362 wptr => work(1:lwork)
2364 allocate(wrk(lwork), stat = istat)
2365 if (istat /= 0)
then
2367 call errmgr%report_error(
"solve_least_squares_vec", &
2368 "Insufficient memory available.", &
2369 la_out_of_memory_error)
2376 call dgels(
'N', m, n, 1, a, m, b, maxmn, wptr, lwork, flag)
2378 call errmgr%report_error(
"solve_least_squares_mtx", &
2379 "The supplied matrix is not of full rank; therefore, " // &
2380 "the solution could not be computed via this routine. " // &
2381 "Try a routine that utilizes column pivoting.", &
2382 la_invalid_operation_error)
2387 module subroutine solve_least_squares_vec_cmplx(a, b, work, olwork, err)
2389 complex(real64),
intent(inout),
dimension(:,:) :: a
2390 complex(real64),
intent(inout),
dimension(:) :: b
2391 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2392 integer(int32),
intent(out),
optional :: olwork
2393 class(errors),
intent(inout),
optional,
target :: err
2396 integer(int32) :: m, n, maxmn, lwork, istat, flag
2397 complex(real64),
pointer,
dimension(:) :: wptr
2398 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2399 complex(real64),
dimension(1) :: temp
2400 class(errors),
pointer :: errmgr
2401 type(errors),
target :: deferr
2407 if (
present(err))
then
2414 if (
size(b) /= maxmn)
then
2415 call errmgr%report_error(
"solve_least_squares_vec_cmplx", &
2416 "Input 2 is not sized correctly.", la_array_size_error)
2421 call zgels(
'N', m, n, 1, a, m, b, maxmn, temp, -1, flag)
2422 lwork = int(temp(1), int32)
2423 if (
present(olwork))
then
2429 if (
present(work))
then
2430 if (
size(work) < lwork)
then
2432 call errmgr%report_error(
"solve_least_squares_vec_cmplx", &
2433 "Incorrectly sized input array WORK, argument 3.", &
2434 la_array_size_error)
2437 wptr => work(1:lwork)
2439 allocate(wrk(lwork), stat = istat)
2440 if (istat /= 0)
then
2442 call errmgr%report_error(
"solve_least_squares_vec_cmplx", &
2443 "Insufficient memory available.", &
2444 la_out_of_memory_error)
2451 call zgels(
'N', m, n, 1, a, m, b, maxmn, wptr, lwork, flag)
2453 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2454 "The supplied matrix is not of full rank; therefore, " // &
2455 "the solution could not be computed via this routine. " // &
2456 "Try a routine that utilizes column pivoting.", &
2457 la_invalid_operation_error)
2462 module subroutine solve_least_squares_mtx_pvt(a, b, ipvt, arnk, work, olwork, err)
2464 real(real64),
intent(inout),
dimension(:,:) :: a, b
2465 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2466 integer(int32),
intent(out),
optional :: arnk
2467 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2468 integer(int32),
intent(out),
optional :: olwork
2469 class(errors),
intent(inout),
optional,
target :: err
2472 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag, rnk
2473 real(real64),
pointer,
dimension(:) :: wptr
2474 real(real64),
allocatable,
target,
dimension(:) :: wrk
2475 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2476 integer(int32),
pointer,
dimension(:) :: iptr
2477 real(real64),
dimension(1) :: temp
2478 integer(int32),
dimension(1) :: itemp
2480 class(errors),
pointer :: errmgr
2481 type(errors),
target :: deferr
2482 character(len = 128) :: errmsg
2490 if (
present(arnk)) arnk = 0
2491 if (
present(err))
then
2499 if (
size(b, 1) /= maxmn)
then
2503 write(errmsg, 100)
"Input number ", flag, &
2504 " is not sized correctly."
2505 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2506 trim(errmsg), la_array_size_error)
2511 call dgelsy(m, n, nrhs, a, m, b, maxmn, itemp, rc, rnk, temp, -1, flag)
2512 lwork = int(temp(1), int32)
2513 if (
present(olwork))
then
2519 if (
present(ipvt))
then
2520 if (
size(ipvt) < n)
then
2522 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2523 "Incorrectly sized pivot array, argument 3.", &
2524 la_array_size_error)
2529 allocate(iwrk(n), stat = istat)
2530 if (istat /= 0)
then
2532 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2533 "Insufficient memory available.", &
2534 la_out_of_memory_error)
2541 if (
present(work))
then
2542 if (
size(work) < lwork)
then
2544 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2545 "Incorrectly sized input array WORK, argument 5.", &
2546 la_array_size_error)
2549 wptr => work(1:lwork)
2551 allocate(wrk(lwork), stat = istat)
2552 if (istat /= 0)
then
2554 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2555 "Insufficient memory available.", &
2556 la_out_of_memory_error)
2563 call dgelsy(m, n, nrhs, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2565 if (
present(arnk)) arnk = rnk
2572 module subroutine solve_least_squares_mtx_pvt_cmplx(a, b, ipvt, arnk, &
2573 work, olwork, rwork, err)
2575 complex(real64),
intent(inout),
dimension(:,:) :: a, b
2576 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2577 integer(int32),
intent(out),
optional :: arnk
2578 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2579 integer(int32),
intent(out),
optional :: olwork
2580 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
2581 class(errors),
intent(inout),
optional,
target :: err
2584 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag, rnk, lrwork
2585 complex(real64),
pointer,
dimension(:) :: wptr
2586 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2587 real(real64),
pointer,
dimension(:) :: rwptr
2588 real(real64),
allocatable,
target,
dimension(:) :: rwrk
2589 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2590 integer(int32),
pointer,
dimension(:) :: iptr
2591 complex(real64),
dimension(1) :: temp
2592 real(real64),
dimension(1) :: rtemp
2593 integer(int32),
dimension(1) :: itemp
2595 class(errors),
pointer :: errmgr
2596 type(errors),
target :: deferr
2597 character(len = 128) :: errmsg
2606 if (
present(arnk)) arnk = 0
2607 if (
present(err))
then
2615 if (
size(b, 1) /= maxmn)
then
2619 write(errmsg, 100)
"Input number ", flag, &
2620 " is not sized correctly."
2621 call errmgr%report_error(
"solve_least_squares_mtx_pvt_cmplx", &
2622 trim(errmsg), la_array_size_error)
2627 call zgelsy(m, n, nrhs, a, m, b, maxmn, itemp, rc, rnk, temp, -1, &
2629 lwork = int(temp(1), int32)
2630 if (
present(olwork))
then
2636 if (
present(ipvt))
then
2637 if (
size(ipvt) < n)
then
2639 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2640 "Incorrectly sized pivot array, argument 3.", &
2641 la_array_size_error)
2646 allocate(iwrk(n), stat = istat)
2647 if (istat /= 0)
then
2649 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2650 "Insufficient memory available.", &
2651 la_out_of_memory_error)
2658 if (
present(work))
then
2659 if (
size(work) < lwork)
then
2661 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2662 "Incorrectly sized input array WORK, argument 5.", &
2663 la_array_size_error)
2666 wptr => work(1:lwork)
2668 allocate(wrk(lwork), stat = istat)
2669 if (istat /= 0)
then
2671 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2672 "Insufficient memory available.", &
2673 la_out_of_memory_error)
2679 if (
present(rwork))
then
2680 if (
size(rwork) < lrwork)
then
2682 call errmgr%report_error(
"solve_least_squares_mtx_pvt_cmplx", &
2683 "Incorrectly sized input array RWORK, argument 7.", &
2684 la_array_size_error)
2687 rwptr => rwork(1:lrwork)
2689 allocate(rwrk(lrwork), stat = istat)
2690 if (istat /= 0)
then
2692 call errmgr%report_error(
"solve_least_squares_mtx_pvt_cmplx", &
2693 "Insufficient memory available.", &
2694 la_out_of_memory_error)
2701 call zgelsy(m, n, nrhs, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2703 if (
present(arnk)) arnk = rnk
2710 module subroutine solve_least_squares_vec_pvt(a, b, ipvt, arnk, work, olwork, err)
2712 real(real64),
intent(inout),
dimension(:,:) :: a
2713 real(real64),
intent(inout),
dimension(:) :: b
2714 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2715 integer(int32),
intent(out),
optional :: arnk
2716 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2717 integer(int32),
intent(out),
optional :: olwork
2718 class(errors),
intent(inout),
optional,
target :: err
2721 integer(int32) :: m, n, maxmn, lwork, istat, flag, rnk
2722 real(real64),
pointer,
dimension(:) :: wptr
2723 real(real64),
allocatable,
target,
dimension(:) :: wrk
2724 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2725 integer(int32),
pointer,
dimension(:) :: iptr
2726 real(real64),
dimension(1) :: temp
2727 integer(int32),
dimension(1) :: itemp
2729 class(errors),
pointer :: errmgr
2730 type(errors),
target :: deferr
2731 character(len = 128) :: errmsg
2738 if (
present(arnk)) arnk = 0
2739 if (
present(err))
then
2747 if (
size(b, 1) /= maxmn)
then
2751 write(errmsg, 100)
"Input number ", flag, &
2752 " is not sized correctly."
2753 call errmgr%report_error(
"solve_least_squares_vec_pvt", &
2754 trim(errmsg), la_array_size_error)
2759 call dgelsy(m, n, 1, a, m, b, maxmn, itemp, rc, rnk, temp, -1, flag)
2760 lwork = int(temp(1), int32)
2761 if (
present(olwork))
then
2767 if (
present(ipvt))
then
2768 if (
size(ipvt) < n)
then
2770 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2771 "Incorrectly sized pivot array, argument 3.", &
2772 la_array_size_error)
2777 allocate(iwrk(n), stat = istat)
2778 if (istat /= 0)
then
2780 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2781 "Insufficient memory available.", &
2782 la_out_of_memory_error)
2789 if (
present(work))
then
2790 if (
size(work) < lwork)
then
2792 call errmgr%report_error(
"solve_least_squares_vec_pvt", &
2793 "Incorrectly sized input array WORK, argument 5.", &
2794 la_array_size_error)
2797 wptr => work(1:lwork)
2799 allocate(wrk(lwork), stat = istat)
2800 if (istat /= 0)
then
2802 call errmgr%report_error(
"solve_least_squares_vec_pvt", &
2803 "Insufficient memory available.", &
2804 la_out_of_memory_error)
2811 call dgelsy(m, n, 1, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, flag)
2812 if (
present(arnk)) arnk = rnk
2819 module subroutine solve_least_squares_vec_pvt_cmplx(a, b, ipvt, arnk, &
2820 work, olwork, rwork, err)
2822 complex(real64),
intent(inout),
dimension(:,:) :: a
2823 complex(real64),
intent(inout),
dimension(:) :: b
2824 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2825 integer(int32),
intent(out),
optional :: arnk
2826 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2827 integer(int32),
intent(out),
optional :: olwork
2828 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
2829 class(errors),
intent(inout),
optional,
target :: err
2832 integer(int32) :: m, n, maxmn, lwork, istat, flag, rnk
2833 complex(real64),
pointer,
dimension(:) :: wptr
2834 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2835 real(real64),
pointer,
dimension(:) :: rwptr
2836 real(real64),
allocatable,
target,
dimension(:) :: rwrk
2837 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2838 integer(int32),
pointer,
dimension(:) :: iptr
2839 complex(real64),
dimension(1) :: temp
2840 real(real64),
dimension(1) :: rtemp
2841 integer(int32),
dimension(1) :: itemp
2843 class(errors),
pointer :: errmgr
2844 type(errors),
target :: deferr
2845 character(len = 128) :: errmsg
2853 if (
present(arnk)) arnk = 0
2854 if (
present(err))
then
2862 if (
size(b, 1) /= maxmn)
then
2866 write(errmsg, 100)
"Input number ", flag, &
2867 " is not sized correctly."
2868 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2869 trim(errmsg), la_array_size_error)
2874 call zgelsy(m, n, 1, a, m, b, maxmn, itemp, rc, rnk, temp, -1, rtemp, &
2876 lwork = int(temp(1), int32)
2877 if (
present(olwork))
then
2883 if (
present(ipvt))
then
2884 if (
size(ipvt) < n)
then
2886 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2887 "Incorrectly sized pivot array, argument 3.", &
2888 la_array_size_error)
2893 allocate(iwrk(n), stat = istat)
2894 if (istat /= 0)
then
2896 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2897 "Insufficient memory available.", &
2898 la_out_of_memory_error)
2905 if (
present(work))
then
2906 if (
size(work) < lwork)
then
2908 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2909 "Incorrectly sized input array WORK, argument 5.", &
2910 la_array_size_error)
2913 wptr => work(1:lwork)
2915 allocate(wrk(lwork), stat = istat)
2916 if (istat /= 0)
then
2918 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2919 "Insufficient memory available.", &
2920 la_out_of_memory_error)
2926 if (
present(rwork))
then
2927 if (
size(rwork) < lrwork)
then
2929 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2930 "Incorrectly sized input array RWORK, argument 7.", &
2931 la_array_size_error)
2934 rwptr => rwork(1:lrwork)
2936 allocate(rwrk(lrwork), stat = istat)
2937 if (istat /= 0)
then
2939 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2940 "Insufficient memory available.", &
2941 la_out_of_memory_error)
2948 call zgelsy(m, n, 1, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2950 if (
present(arnk)) arnk = rnk
2957 module subroutine solve_least_squares_mtx_svd(a, b, s, arnk, work, olwork, err)
2959 real(real64),
intent(inout),
dimension(:,:) :: a, b
2960 integer(int32),
intent(out),
optional :: arnk
2961 real(real64),
intent(out),
target,
optional,
dimension(:) :: work, s
2962 integer(int32),
intent(out),
optional :: olwork
2963 class(errors),
intent(inout),
optional,
target :: err
2966 integer(int32) :: m, n, nrhs, mn, maxmn, istat, flag, lwork, rnk
2967 real(real64),
pointer,
dimension(:) :: wptr, sptr
2968 real(real64),
allocatable,
target,
dimension(:) :: wrk, sing
2969 real(real64),
dimension(1) :: temp
2970 real(real64) :: rcond
2971 class(errors),
pointer :: errmgr
2972 type(errors),
target :: deferr
2973 character(len = 128) :: errmsg
2981 rcond = epsilon(rcond)
2982 if (
present(arnk)) arnk = 0
2983 if (
present(err))
then
2991 if (
size(b, 1) /= maxmn)
then
2996 write(errmsg, 100)
"Input number ", flag, &
2997 " is not sized correctly."
2998 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
2999 trim(errmsg), la_array_size_error)
3004 call dgelss(m, n, nrhs, a, m, b, maxmn, temp, rcond, rnk, temp, -1, &
3006 lwork = int(temp(1), int32)
3007 if (
present(olwork))
then
3013 if (
present(s))
then
3014 if (
size(s) < mn)
then
3016 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3017 "Incorrectly sized input array S, argument 3.", &
3018 la_array_size_error)
3023 allocate(sing(mn), stat = istat)
3024 if (istat /= 0)
then
3026 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3027 "Insufficient memory available.", &
3028 la_out_of_memory_error)
3034 if (
present(work))
then
3035 if (
size(work) < lwork)
then
3037 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3038 "Incorrectly sized input array WORK, argument 5.", &
3039 la_array_size_error)
3042 wptr => work(1:lwork)
3044 allocate(wrk(lwork), stat = istat)
3045 if (istat /= 0)
then
3047 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3048 "Insufficient memory available.", &
3049 la_out_of_memory_error)
3056 call dgelss(m, n, nrhs, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3058 if (
present(arnk)) arnk = rnk
3060 write(errmsg, 101) flag,
" superdiagonals could not " // &
3061 "converge to zero as part of the QR iteration process."
3062 call errmgr%report_warning(
"solve_least_squares_mtx_svd", errmsg, &
3063 la_convergence_error)
3072 module subroutine solve_least_squares_mtx_svd_cmplx(a, b, s, arnk, work, &
3075 complex(real64),
intent(inout),
dimension(:,:) :: a, b
3076 integer(int32),
intent(out),
optional :: arnk
3077 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3078 real(real64),
intent(out),
target,
optional,
dimension(:) :: s, rwork
3079 integer(int32),
intent(out),
optional :: olwork
3080 class(errors),
intent(inout),
optional,
target :: err
3083 integer(int32) :: m, n, nrhs, mn, maxmn, istat, flag, lwork, rnk, lrwork
3084 complex(real64),
pointer,
dimension(:) :: wptr
3085 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3086 real(real64),
pointer,
dimension(:) :: rwptr, sptr
3087 real(real64),
allocatable,
target,
dimension(:) :: rwrk, sing
3088 complex(real64),
dimension(1) :: temp
3089 real(real64),
dimension(1) :: rtemp
3090 real(real64) :: rcond
3091 class(errors),
pointer :: errmgr
3092 type(errors),
target :: deferr
3093 character(len = 128) :: errmsg
3102 rcond = epsilon(rcond)
3103 if (
present(arnk)) arnk = 0
3104 if (
present(err))
then
3112 if (
size(b, 1) /= maxmn)
then
3117 write(errmsg, 100)
"Input number ", flag, &
3118 " is not sized correctly."
3119 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3120 trim(errmsg), la_array_size_error)
3125 call zgelss(m, n, nrhs, a, m, b, maxmn, rtemp, rcond, rnk, temp, -1, &
3127 lwork = int(temp(1), int32)
3128 if (
present(olwork))
then
3134 if (
present(s))
then
3135 if (
size(s) < mn)
then
3137 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3138 "Incorrectly sized input array S, argument 3.", &
3139 la_array_size_error)
3144 allocate(sing(mn), stat = istat)
3145 if (istat /= 0)
then
3147 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3148 "Insufficient memory available.", &
3149 la_out_of_memory_error)
3155 if (
present(work))
then
3156 if (
size(work) < lwork)
then
3158 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3159 "Incorrectly sized input array WORK, argument 5.", &
3160 la_array_size_error)
3163 wptr => work(1:lwork)
3165 allocate(wrk(lwork), stat = istat)
3166 if (istat /= 0)
then
3168 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3169 "Insufficient memory available.", &
3170 la_out_of_memory_error)
3176 if (
present(rwork))
then
3177 if (
size(rwork) < lrwork)
then
3179 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3180 "Incorrectly sized input array RWORK, argument 7.", &
3181 la_array_size_error)
3184 rwptr => rwork(1:lrwork)
3186 allocate(rwrk(lrwork), stat = istat)
3187 if (istat /= 0)
then
3189 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3190 "Insufficient memory available.", &
3191 la_out_of_memory_error)
3198 call zgelss(m, n, nrhs, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3200 if (
present(arnk)) arnk = rnk
3202 write(errmsg, 101) flag,
" superdiagonals could not " // &
3203 "converge to zero as part of the QR iteration process."
3204 call errmgr%report_warning(
"solve_least_squares_mtx_svd_cmplx", &
3205 errmsg, la_convergence_error)
3214 module subroutine solve_least_squares_vec_svd(a, b, s, arnk, work, olwork, err)
3216 real(real64),
intent(inout),
dimension(:,:) :: a
3217 real(real64),
intent(inout),
dimension(:) :: b
3218 integer(int32),
intent(out),
optional :: arnk
3219 real(real64),
intent(out),
target,
optional,
dimension(:) :: work, s
3220 integer(int32),
intent(out),
optional :: olwork
3221 class(errors),
intent(inout),
optional,
target :: err
3224 integer(int32) :: m, n, mn, maxmn, istat, flag, lwork, rnk
3225 real(real64),
pointer,
dimension(:) :: wptr, sptr
3226 real(real64),
allocatable,
target,
dimension(:) :: wrk, sing
3227 real(real64),
dimension(1) :: temp
3228 real(real64) :: rcond
3229 class(errors),
pointer :: errmgr
3230 type(errors),
target :: deferr
3231 character(len = 128) :: errmsg
3238 rcond = epsilon(rcond)
3239 if (
present(arnk)) arnk = 0
3240 if (
present(err))
then
3248 if (
size(b) /= maxmn)
then
3253 write(errmsg, 100)
"Input number ", flag, &
3254 " is not sized correctly."
3255 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3256 trim(errmsg), la_array_size_error)
3261 call dgelss(m, n, 1, a, m, b, maxmn, temp, rcond, rnk, temp, -1, flag)
3262 lwork = int(temp(1), int32)
3263 if (
present(olwork))
then
3269 if (
present(s))
then
3270 if (
size(s) < mn)
then
3272 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3273 "Incorrectly sized input array S, argument 3.", &
3274 la_array_size_error)
3279 allocate(sing(mn), stat = istat)
3280 if (istat /= 0)
then
3282 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3283 "Insufficient memory available.", &
3284 la_out_of_memory_error)
3290 if (
present(work))
then
3291 if (
size(work) < lwork)
then
3293 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3294 "Incorrectly sized input array WORK, argument 5.", &
3295 la_array_size_error)
3298 wptr => work(1:lwork)
3300 allocate(wrk(lwork), stat = istat)
3301 if (istat /= 0)
then
3303 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3304 "Insufficient memory available.", &
3305 la_out_of_memory_error)
3312 call dgelss(m, n, 1, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3314 if (
present(arnk)) arnk = rnk
3316 write(errmsg, 101) flag,
" superdiagonals could not " // &
3317 "converge to zero as part of the QR iteration process."
3318 call errmgr%report_warning(
"solve_least_squares_vec_svd", errmsg, &
3319 la_convergence_error)
3328 module subroutine solve_least_squares_vec_svd_cmplx(a, b, s, arnk, work, &
3331 complex(real64),
intent(inout),
dimension(:,:) :: a
3332 complex(real64),
intent(inout),
dimension(:) :: b
3333 integer(int32),
intent(out),
optional :: arnk
3334 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3335 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork, s
3336 integer(int32),
intent(out),
optional :: olwork
3337 class(errors),
intent(inout),
optional,
target :: err
3340 integer(int32) :: m, n, mn, maxmn, istat, flag, lwork, rnk, lrwork
3341 real(real64),
pointer,
dimension(:) :: rwptr, sptr
3342 real(real64),
allocatable,
target,
dimension(:) :: rwrk, sing
3343 complex(real64),
pointer,
dimension(:) :: wptr
3344 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3345 complex(real64),
dimension(1) :: temp
3346 real(real64),
dimension(1) :: rtemp
3347 real(real64) :: rcond
3348 class(errors),
pointer :: errmgr
3349 type(errors),
target :: deferr
3350 character(len = 128) :: errmsg
3358 rcond = epsilon(rcond)
3359 if (
present(arnk)) arnk = 0
3360 if (
present(err))
then
3368 if (
size(b) /= maxmn)
then
3373 write(errmsg, 100)
"Input number ", flag, &
3374 " is not sized correctly."
3375 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3376 trim(errmsg), la_array_size_error)
3381 call zgelss(m, n, 1, a, m, b, maxmn, rtemp, rcond, rnk, temp, -1, &
3383 lwork = int(temp(1), int32)
3384 if (
present(olwork))
then
3390 if (
present(s))
then
3391 if (
size(s) < mn)
then
3393 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3394 "Incorrectly sized input array S, argument 3.", &
3395 la_array_size_error)
3400 allocate(sing(mn), stat = istat)
3401 if (istat /= 0)
then
3403 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3404 "Insufficient memory available.", &
3405 la_out_of_memory_error)
3411 if (
present(work))
then
3412 if (
size(work) < lwork)
then
3414 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3415 "Incorrectly sized input array WORK, argument 5.", &
3416 la_array_size_error)
3419 wptr => work(1:lwork)
3421 allocate(wrk(lwork), stat = istat)
3422 if (istat /= 0)
then
3424 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3425 "Insufficient memory available.", &
3426 la_out_of_memory_error)
3432 if (
present(rwork))
then
3433 if (
size(rwork) < lrwork)
then
3435 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3436 "Incorrectly sized input array RWORK, argument 7.", &
3437 la_array_size_error)
3440 rwptr => rwork(1:lrwork)
3442 allocate(rwrk(lrwork), stat = istat)
3443 if (istat /= 0)
then
3445 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3446 "Insufficient memory available.", &
3447 la_out_of_memory_error)
3454 call zgelss(m, n, 1, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3456 if (
present(arnk)) arnk = rnk
3458 write(errmsg, 101) flag,
" superdiagonals could not " // &
3459 "converge to zero as part of the QR iteration process."
3460 call errmgr%report_warning(
"solve_least_squares_vec_svd_cmplx", &
3461 errmsg, la_convergence_error)
3472 module subroutine solve_lq_mtx(a, tau, b, work, olwork, err)
3474 real(real64),
intent(in),
dimension(:,:) :: a
3475 real(real64),
intent(in),
dimension(:) :: tau
3476 real(real64),
intent(inout),
dimension(:,:) :: b
3477 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
3478 integer(int32),
intent(out),
optional :: olwork
3479 class(errors),
intent(inout),
optional,
target :: err
3482 real(real64),
parameter :: one = 1.0d0
3485 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
3486 real(real64),
pointer,
dimension(:) :: wptr
3487 real(real64),
allocatable,
target,
dimension(:) :: wrk
3488 class(errors),
pointer :: errmgr
3489 type(errors),
target :: deferr
3490 character(len = 128) :: errmsg
3497 if (
present(err))
then
3507 else if (
size(tau) /= k)
then
3509 else if (
size(b, 1) /= n)
then
3515 write(errmsg, 100)
"Input number ", flag, &
3516 " is not sized correctly."
3517 call errmgr%report_error(
"solve_lq_mtx", trim(errmsg), &
3518 la_array_size_error)
3523 call mult_lq(.true., .true., a, tau, b, olwork = lwork)
3525 if (
present(olwork))
then
3531 if (
present(work))
then
3532 if (
size(work) < lwork)
then
3534 call errmgr%report_error(
"solve_lq_mtx", &
3535 "Incorrectly sized input array WORK, argument 4.", &
3536 la_array_size_error)
3539 wptr => work(1:lwork)
3541 allocate(wrk(lwork), stat = istat)
3542 if (istat /= 0)
then
3544 call errmgr%report_error(
"solve_lq_mtx", &
3545 "Insufficient memory available.", &
3546 la_out_of_memory_error)
3554 call solve_triangular_system(.true., .false., .false., .true., one, &
3555 a(1:m,1:m), b(1:m,:), errmgr)
3556 if (errmgr%has_error_occurred())
return
3559 call mult_lq(.true., .true., a, tau, b, work = wptr, err = errmgr)
3560 if (errmgr%has_error_occurred())
return
3567 module subroutine solve_lq_mtx_cmplx(a, tau, b, work, olwork, err)
3569 complex(real64),
intent(in),
dimension(:,:) :: a
3570 complex(real64),
intent(in),
dimension(:) :: tau
3571 complex(real64),
intent(inout),
dimension(:,:) :: b
3572 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3573 integer(int32),
intent(out),
optional :: olwork
3574 class(errors),
intent(inout),
optional,
target :: err
3577 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
3580 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
3581 complex(real64),
pointer,
dimension(:) :: wptr
3582 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3583 class(errors),
pointer :: errmgr
3584 type(errors),
target :: deferr
3585 character(len = 128) :: errmsg
3592 if (
present(err))
then
3602 else if (
size(tau) /= k)
then
3604 else if (
size(b, 1) /= n)
then
3610 write(errmsg, 100)
"Input number ", flag, &
3611 " is not sized correctly."
3612 call errmgr%report_error(
"solve_lq_mtx_cmplx", trim(errmsg), &
3613 la_array_size_error)
3618 call mult_lq(.true., .true., a, tau, b, olwork = lwork)
3620 if (
present(olwork))
then
3626 if (
present(work))
then
3627 if (
size(work) < lwork)
then
3629 call errmgr%report_error(
"solve_lq_mtx_cmplx", &
3630 "Incorrectly sized input array WORK, argument 4.", &
3631 la_array_size_error)
3634 wptr => work(1:lwork)
3636 allocate(wrk(lwork), stat = istat)
3637 if (istat /= 0)
then
3639 call errmgr%report_error(
"solve_lq_mtx_cmplx", &
3640 "Insufficient memory available.", &
3641 la_out_of_memory_error)
3649 call solve_triangular_system(.true., .false., .false., .true., one, &
3650 a(1:m,1:m), b(1:m,:), errmgr)
3651 if (errmgr%has_error_occurred())
return
3654 call mult_lq(.true., .true., a, tau, b, work = wptr, err = errmgr)
3655 if (errmgr%has_error_occurred())
return
3662 module subroutine solve_lq_vec(a, tau, b, work, olwork, err)
3664 real(real64),
intent(in),
dimension(:,:) :: a
3665 real(real64),
intent(in),
dimension(:) :: tau
3666 real(real64),
intent(inout),
dimension(:) :: b
3667 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
3668 integer(int32),
intent(out),
optional :: olwork
3669 class(errors),
intent(inout),
optional,
target :: err
3672 integer(int32) :: m, n, k, lwork, flag, istat
3673 real(real64),
pointer,
dimension(:) :: wptr
3674 real(real64),
allocatable,
target,
dimension(:) :: wrk
3675 class(errors),
pointer :: errmgr
3676 type(errors),
target :: deferr
3677 character(len = 128) :: errmsg
3683 if (
present(err))
then
3693 else if (
size(tau) /= k)
then
3695 else if (
size(b) /= n)
then
3701 write(errmsg, 100)
"Input number ", flag, &
3702 " is not sized correctly."
3703 call errmgr%report_error(
"solve_lq_vec", trim(errmsg), &
3704 la_array_size_error)
3709 call mult_lq(.true., a, tau, b, olwork = lwork)
3711 if (
present(olwork))
then
3717 if (
present(work))
then
3718 if (
size(work) < lwork)
then
3720 call errmgr%report_error(
"solve_lq_vec", &
3721 "Incorrectly sized input array WORK, argument 4.", &
3722 la_array_size_error)
3725 wptr => work(1:lwork)
3727 allocate(wrk(lwork), stat = istat)
3728 if (istat /= 0)
then
3730 call errmgr%report_error(
"solve_lq_vec", &
3731 "Insufficient memory available.", &
3732 la_out_of_memory_error)
3740 call solve_triangular_system(.false., .false., .true., a(1:m,1:m), &
3742 if (errmgr%has_error_occurred())
return
3745 call mult_lq(.true., a, tau, b, work = wptr, err = errmgr)
3746 if (errmgr%has_error_occurred())
return
3753 module subroutine solve_lq_vec_cmplx(a, tau, b, work, olwork, err)
3755 complex(real64),
intent(in),
dimension(:,:) :: a
3756 complex(real64),
intent(in),
dimension(:) :: tau
3757 complex(real64),
intent(inout),
dimension(:) :: b
3758 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3759 integer(int32),
intent(out),
optional :: olwork
3760 class(errors),
intent(inout),
optional,
target :: err
3763 integer(int32) :: m, n, k, lwork, flag, istat
3764 complex(real64),
pointer,
dimension(:) :: wptr
3765 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3766 class(errors),
pointer :: errmgr
3767 type(errors),
target :: deferr
3768 character(len = 128) :: errmsg
3774 if (
present(err))
then
3784 else if (
size(tau) /= k)
then
3786 else if (
size(b) /= n)
then
3792 write(errmsg, 100)
"Input number ", flag, &
3793 " is not sized correctly."
3794 call errmgr%report_error(
"solve_lq_vec_cmplx", trim(errmsg), &
3795 la_array_size_error)
3800 call mult_lq(.true., a, tau, b, olwork = lwork)
3802 if (
present(olwork))
then
3808 if (
present(work))
then
3809 if (
size(work) < lwork)
then
3811 call errmgr%report_error(
"solve_lq_vec_cmplx", &
3812 "Incorrectly sized input array WORK, argument 4.", &
3813 la_array_size_error)
3816 wptr => work(1:lwork)
3818 allocate(wrk(lwork), stat = istat)
3819 if (istat /= 0)
then
3821 call errmgr%report_error(
"solve_lq_vec_cmplx", &
3822 "Insufficient memory available.", &
3823 la_out_of_memory_error)
3831 call solve_triangular_system(.false., .false., .true., a(1:m,1:m), &
3833 if (errmgr%has_error_occurred())
return
3836 call mult_lq(.true., a, tau, b, work = wptr, err = errmgr)
3837 if (errmgr%has_error_occurred())
return
Provides a set of common linear algebra routines.