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 = :),
allocatable :: 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 allocate(
character(len = 256) :: errmsg)
284 write(errmsg, 100)
"Input number ", flag, &
285 " is not sized correctly."
286 call errmgr%report_error(
"solve_lu_mtx", trim(errmsg), &
292 call dgetrs(
"N", n, nrhs, a, n, ipvt, b, n, flag)
299 module subroutine solve_lu_mtx_cmplx(a, ipvt, b, err)
301 complex(real64),
intent(in),
dimension(:,:) :: a
302 integer(int32),
intent(in),
dimension(:) :: ipvt
303 complex(real64),
intent(inout),
dimension(:,:) :: b
304 class(errors),
intent(inout),
optional,
target :: err
307 integer(int32) :: n, nrhs, flag
308 class(errors),
pointer :: errmgr
309 type(errors),
target :: deferr
310 character(len = :),
allocatable :: errmsg
315 if (
present(err))
then
323 if (
size(a, 2) /= n)
then
325 else if (
size(ipvt) /= n)
then
327 else if (
size(b, 1) /= n)
then
332 allocate(
character(len = 256) :: errmsg)
333 write(errmsg, 100)
"Input number ", flag, &
334 " is not sized correctly."
335 call errmgr%report_error(
"solve_lu_mtx_cmplx", trim(errmsg), &
341 call zgetrs(
"N", n, nrhs, a, n, ipvt, b, n, flag)
348 module subroutine solve_lu_vec(a, ipvt, b, err)
350 real(real64),
intent(in),
dimension(:,:) :: a
351 integer(int32),
intent(in),
dimension(:) :: ipvt
352 real(real64),
intent(inout),
dimension(:) :: b
353 class(errors),
intent(inout),
optional,
target :: err
356 integer(int32) :: n, flag
357 class(errors),
pointer :: errmgr
358 type(errors),
target :: deferr
359 character(len = :),
allocatable :: errmsg
363 if (
present(err))
then
371 if (
size(a, 2) /= n)
then
373 else if (
size(ipvt) /= n)
then
375 else if (
size(b) /= n)
then
380 allocate(
character(len = 256) :: errmsg)
381 write(errmsg, 100)
"Input number ", flag, &
382 " is not sized correctly."
383 call errmgr%report_error(
"solve_lu_vec", trim(errmsg), &
389 call dgetrs(
"N", n, 1, a, n, ipvt, b, n, flag)
396 module subroutine solve_lu_vec_cmplx(a, ipvt, b, err)
398 complex(real64),
intent(in),
dimension(:,:) :: a
399 integer(int32),
intent(in),
dimension(:) :: ipvt
400 complex(real64),
intent(inout),
dimension(:) :: b
401 class(errors),
intent(inout),
optional,
target :: err
404 integer(int32) :: n, flag
405 class(errors),
pointer :: errmgr
406 type(errors),
target :: deferr
407 character(len = :),
allocatable :: errmsg
411 if (
present(err))
then
419 if (
size(a, 2) /= n)
then
421 else if (
size(ipvt) /= n)
then
423 else if (
size(b) /= n)
then
428 allocate(
character(len = 256) :: errmsg)
429 write(errmsg, 100)
"Input number ", flag, &
430 " is not sized correctly."
431 call errmgr%report_error(
"solve_lu_vec_cmplx", trim(errmsg), &
437 call zgetrs(
"N", n, 1, a, n, ipvt, b, n, flag)
446 module subroutine solve_qr_no_pivot_mtx(a, tau, b, work, olwork, err)
448 real(real64),
intent(inout),
dimension(:,:) :: a, b
449 real(real64),
intent(in),
dimension(:) :: tau
450 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
451 integer(int32),
intent(out),
optional :: olwork
452 class(errors),
intent(inout),
optional,
target :: err
455 real(real64),
parameter :: one = 1.0d0
458 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
459 real(real64),
pointer,
dimension(:) :: wptr
460 real(real64),
allocatable,
target,
dimension(:) :: wrk
461 class(errors),
pointer :: errmgr
462 type(errors),
target :: deferr
463 character(len = :),
allocatable :: errmsg
470 if (
present(err))
then
480 else if (
size(tau) /= k)
then
482 else if (
size(b, 1) /= m)
then
487 allocate(
character(len = 256) :: errmsg)
488 write(errmsg, 100)
"Input number ", flag, &
489 " is not sized correctly."
490 call errmgr%report_error(
"solve_qr_no_pivot_mtx", trim(errmsg), &
496 call mult_qr(.true., .true., a, tau, b, olwork = lwork)
497 if (
present(olwork))
then
503 if (
present(work))
then
504 if (
size(work) < lwork)
then
506 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
507 "Incorrectly sized input array WORK, argument 4.", &
511 wptr => work(1:lwork)
513 allocate(wrk(lwork), stat = istat)
516 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
517 "Insufficient memory available.", &
518 la_out_of_memory_error)
525 call mult_qr(.true., .true., a, tau, b, wptr)
528 call solve_triangular_system(.true., .true., .false., .true., one, &
529 a(1:n,1:n), b(1:n,:))
536 module subroutine solve_qr_no_pivot_mtx_cmplx(a, tau, b, work, olwork, err)
538 complex(real64),
intent(inout),
dimension(:,:) :: a, b
539 complex(real64),
intent(in),
dimension(:) :: tau
540 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
541 integer(int32),
intent(out),
optional :: olwork
542 class(errors),
intent(inout),
optional,
target :: err
545 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
548 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
549 complex(real64),
pointer,
dimension(:) :: wptr
550 complex(real64),
allocatable,
target,
dimension(:) :: wrk
551 class(errors),
pointer :: errmgr
552 type(errors),
target :: deferr
553 character(len = :),
allocatable :: errmsg
560 if (
present(err))
then
570 else if (
size(tau) /= k)
then
572 else if (
size(b, 1) /= m)
then
577 allocate(
character(len = 256) :: errmsg)
578 write(errmsg, 100)
"Input number ", flag, &
579 " is not sized correctly."
580 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
581 trim(errmsg), la_array_size_error)
586 call mult_qr(.true., .true., a, tau, b, olwork = lwork)
587 if (
present(olwork))
then
593 if (
present(work))
then
594 if (
size(work) < lwork)
then
596 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
597 "Incorrectly sized input array WORK, argument 4.", &
601 wptr => work(1:lwork)
603 allocate(wrk(lwork), stat = istat)
606 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
607 "Insufficient memory available.", &
608 la_out_of_memory_error)
615 call mult_qr(.true., .true., a, tau, b, wptr)
618 call solve_triangular_system(.true., .true., .false., .true., one, &
619 a(1:n,1:n), b(1:n,:))
626 module subroutine solve_qr_no_pivot_vec(a, tau, b, work, olwork, err)
628 real(real64),
intent(inout),
dimension(:,:) :: a
629 real(real64),
intent(in),
dimension(:) :: tau
630 real(real64),
intent(inout),
dimension(:) :: b
631 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
632 integer(int32),
intent(out),
optional :: olwork
633 class(errors),
intent(inout),
optional,
target :: err
636 integer(int32) :: m, n, k, flag, lwork, istat
637 real(real64),
pointer,
dimension(:) :: wptr
638 real(real64),
allocatable,
target,
dimension(:) :: wrk
639 class(errors),
pointer :: errmgr
640 type(errors),
target :: deferr
641 character(len = :),
allocatable :: errmsg
647 if (
present(err))
then
657 else if (
size(tau) /= k)
then
659 else if (
size(b) /= m)
then
664 allocate(
character(len = 256) :: errmsg)
665 write(errmsg, 100)
"Input number ", flag, &
666 " is not sized correctly."
667 call errmgr%report_error(
"solve_qr_no_pivot_vec", trim(errmsg), &
673 call mult_qr(.true., a, tau, b, olwork = lwork)
674 if (
present(olwork))
then
680 if (
present(work))
then
681 if (
size(work) < lwork)
then
683 call errmgr%report_error(
"solve_qr_no_pivot_vec", &
684 "Incorrectly sized input array WORK, argument 4.", &
688 wptr => work(1:lwork)
690 allocate(wrk(lwork), stat = istat)
693 call errmgr%report_error(
"solve_qr_no_pivot_vec", &
694 "Insufficient memory available.", &
695 la_out_of_memory_error)
702 call mult_qr(.true., a, tau, b, work = wptr)
705 call solve_triangular_system(.true., .false., .true., a(1:n,1:n), b)
712 module subroutine solve_qr_no_pivot_vec_cmplx(a, tau, b, work, olwork, err)
714 complex(real64),
intent(inout),
dimension(:,:) :: a
715 complex(real64),
intent(in),
dimension(:) :: tau
716 complex(real64),
intent(inout),
dimension(:) :: b
717 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
718 integer(int32),
intent(out),
optional :: olwork
719 class(errors),
intent(inout),
optional,
target :: err
722 integer(int32) :: m, n, k, flag, lwork, istat
723 complex(real64),
pointer,
dimension(:) :: wptr
724 complex(real64),
allocatable,
target,
dimension(:) :: wrk
725 class(errors),
pointer :: errmgr
726 type(errors),
target :: deferr
727 character(len = :),
allocatable :: errmsg
733 if (
present(err))
then
743 else if (
size(tau) /= k)
then
745 else if (
size(b) /= m)
then
750 allocate(
character(len = 256) :: errmsg)
751 write(errmsg, 100)
"Input number ", flag, &
752 " is not sized correctly."
753 call errmgr%report_error(
"solve_qr_no_pivot_vec_cmplx", &
754 trim(errmsg), la_array_size_error)
759 call mult_qr(.true., a, tau, b, olwork = lwork)
760 if (
present(olwork))
then
766 if (
present(work))
then
767 if (
size(work) < lwork)
then
769 call errmgr%report_error(
"solve_qr_no_pivot_vec_cmplx", &
770 "Incorrectly sized input array WORK, argument 4.", &
774 wptr => work(1:lwork)
776 allocate(wrk(lwork), stat = istat)
779 call errmgr%report_error(
"solve_qr_no_pivot_vec_cmplx", &
780 "Insufficient memory available.", &
781 la_out_of_memory_error)
788 call mult_qr(.true., a, tau, b, work = wptr)
791 call solve_triangular_system(.true., .false., .true., a(1:n,1:n), b)
798 module subroutine solve_qr_pivot_mtx(a, tau, jpvt, b, work, olwork, err)
800 real(real64),
intent(inout),
dimension(:,:) :: a
801 real(real64),
intent(in),
dimension(:) :: tau
802 integer(int32),
intent(in),
dimension(:) :: jpvt
803 real(real64),
intent(inout),
dimension(:,:) :: b
804 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
805 integer(int32),
intent(out),
optional :: olwork
806 class(errors),
intent(inout),
optional,
target :: err
809 integer(int32),
parameter :: imin = 2
810 integer(int32),
parameter :: imax = 1
811 real(real64),
parameter :: zero = 0.0d0
812 real(real64),
parameter :: one = 1.0d0
815 integer(int32) :: i, j, m, n, mn, nrhs, lwork, ismin, ismax, &
816 rnk, maxmn, flag, istat, lwork1, lwork2, lwork3
817 real(real64) :: rcond, smax, smin, smaxpr, sminpr, s1, c1, s2, c2
818 real(real64),
pointer,
dimension(:) :: wptr, w, tau2
819 real(real64),
allocatable,
target,
dimension(:) :: wrk
820 class(errors),
pointer :: errmgr
821 type(errors),
target :: deferr
822 character(len = :),
allocatable :: errmsg
832 rcond = epsilon(rcond)
833 if (
present(err))
then
841 if (
size(tau) /= mn)
then
843 else if (
size(jpvt) /= n)
then
845 else if (
size(b, 1) /= maxmn)
then
850 allocate(
character(len = 256) :: errmsg)
851 write(errmsg, 100)
"Input number ", flag, &
852 " is not sized correctly."
853 call errmgr%report_error(
"solve_qr_pivot_mtx", trim(errmsg), &
859 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
860 call mult_qr(.true., .true., a, tau, b(1:m,:), olwork = lwork2)
861 call mult_rz(.true., .true., n, a(1:mn,:), a(1:mn,1), b(1:n,:), &
863 lwork = max(lwork1, lwork2, lwork3, 2 * mn + 1) + mn
864 if (
present(olwork))
then
870 if (
present(work))
then
871 if (
size(work) < lwork)
then
873 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
874 "Incorrectly sized input array WORK, argument 5.", &
878 wptr => work(1:lwork)
880 allocate(wrk(lwork), stat = istat)
883 call errmgr%report_error(
"solve_qr_pivot_mtx", &
884 "Insufficient memory available.", &
885 la_out_of_memory_error)
896 if (abs(a(1,1)) == zero)
then
906 call dlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
907 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
908 call dlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
909 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
910 if (smaxpr * rcond <= sminpr)
then
912 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
913 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
929 w => wptr(rnk+1:lwork)
930 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
933 call mult_qr(.true., .true., a, tau, b(1:m,:), w)
936 call solve_triangular_system(.true., .true., .false., .true., one, &
937 a(1:rnk,1:rnk), b(1:rnk,:))
938 if (n > rnk) b(rnk+1:n,:) = zero
942 call mult_rz(.true., .true., n - rnk, a(1:rnk,:), tau2, b(1:n,:), w)
948 wptr(jpvt(i)) = b(i,j)
958 module subroutine solve_qr_pivot_mtx_cmplx(a, tau, jpvt, b, work, olwork, err)
960 complex(real64),
intent(inout),
dimension(:,:) :: a
961 complex(real64),
intent(in),
dimension(:) :: tau
962 integer(int32),
intent(in),
dimension(:) :: jpvt
963 complex(real64),
intent(inout),
dimension(:,:) :: b
964 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
965 integer(int32),
intent(out),
optional :: olwork
966 class(errors),
intent(inout),
optional,
target :: err
969 integer(int32),
parameter :: imin = 2
970 integer(int32),
parameter :: imax = 1
971 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
972 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
975 integer(int32) :: i, j, m, n, mn, nrhs, lwork, ismin, ismax, &
976 rnk, maxmn, flag, istat, lwork1, lwork2, lwork3
977 real(real64) :: rcond, smax, smin, smaxpr, sminpr
978 complex(real64) :: s1, c1, s2, c2
979 complex(real64),
pointer,
dimension(:) :: wptr, w, tau2
980 complex(real64),
allocatable,
target,
dimension(:) :: wrk
981 class(errors),
pointer :: errmgr
982 type(errors),
target :: deferr
983 character(len = :),
allocatable :: errmsg
993 rcond = epsilon(rcond)
994 if (
present(err))
then
1002 if (
size(tau) /= mn)
then
1004 else if (
size(jpvt) /= n)
then
1006 else if (
size(b, 1) /= maxmn)
then
1011 allocate(
character(len = 256) :: errmsg)
1012 write(errmsg, 100)
"Input number ", flag, &
1013 " is not sized correctly."
1014 call errmgr%report_error(
"solve_qr_pivot_mtx_cmplx", &
1015 trim(errmsg), la_array_size_error)
1020 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1021 call mult_qr(.true., .true., a, tau, b(1:m,:), olwork = lwork2)
1022 call mult_rz(.true., .true., n, a(1:mn,:), a(1:mn,1), b(1:n,:), &
1024 lwork = max(lwork1, lwork2, lwork3, 2 * mn + 1) + mn
1025 if (
present(olwork))
then
1031 if (
present(work))
then
1032 if (
size(work) < lwork)
then
1034 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
1035 "Incorrectly sized input array WORK, argument 5.", &
1036 la_array_size_error)
1039 wptr => work(1:lwork)
1041 allocate(wrk(lwork), stat = istat)
1042 if (istat /= 0)
then
1044 call errmgr%report_error(
"solve_qr_pivot_mtx_cmplx", &
1045 "Insufficient memory available.", &
1046 la_out_of_memory_error)
1057 if (abs(a(1,1)) == zero)
then
1067 call zlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1068 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1069 call zlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1070 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1071 if (smaxpr * rcond <= sminpr)
then
1073 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1074 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1076 wptr(ismin+rnk) = c1
1077 wptr(ismax+rnk) = c2
1090 w => wptr(rnk+1:lwork)
1091 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
1094 call mult_qr(.true., .true., a, tau, b(1:m,:), w)
1097 call solve_triangular_system(.true., .true., .false., .true., one, &
1098 a(1:rnk,1:rnk), b(1:rnk,:))
1099 if (n > rnk) b(rnk+1:n,:) = zero
1103 call mult_rz(.true., .true., n - rnk, a(1:rnk,:), tau2, b(1:n,:), w)
1109 wptr(jpvt(i)) = b(i,j)
1119 module subroutine solve_qr_pivot_vec(a, tau, jpvt, b, work, olwork, err)
1121 real(real64),
intent(inout),
dimension(:,:) :: a
1122 real(real64),
intent(in),
dimension(:) :: tau
1123 integer(int32),
intent(in),
dimension(:) :: jpvt
1124 real(real64),
intent(inout),
dimension(:) :: b
1125 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
1126 integer(int32),
intent(out),
optional :: olwork
1127 class(errors),
intent(inout),
optional,
target :: err
1130 integer(int32),
parameter :: imin = 2
1131 integer(int32),
parameter :: imax = 1
1132 real(real64),
parameter :: zero = 0.0d0
1133 real(real64),
parameter :: one = 1.0d0
1136 integer(int32) :: i, m, n, mn, lwork, ismin, ismax, rnk, maxmn, flag, &
1137 istat, lwork1, lwork2
1138 real(real64) :: rcond, smax, smin, smaxpr, sminpr, s1, c1, s2, c2
1139 real(real64),
pointer,
dimension(:) :: wptr, w, tau2
1140 real(real64),
allocatable,
target,
dimension(:) :: wrk
1141 class(errors),
pointer :: errmgr
1142 type(errors),
target :: deferr
1143 character(len = :),
allocatable :: errmsg
1152 rcond = epsilon(rcond)
1153 if (
present(err))
then
1161 if (
size(tau) /= mn)
then
1163 else if (
size(jpvt) /= n)
then
1165 else if (
size(b) /= maxmn)
then
1170 allocate(
character(len = 256) :: errmsg)
1171 write(errmsg, 100)
"Input number ", flag, &
1172 " is not sized correctly."
1173 call errmgr%report_error(
"solve_qr_pivot_vec", trim(errmsg), &
1174 la_array_size_error)
1179 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1180 call mult_rz(.true., n, a(1:mn,:), a(1:mn,1), b(1:n), olwork = lwork2)
1181 lwork = max(lwork1, lwork2, 2 * mn + 1) + mn
1182 if (
present(olwork))
then
1188 if (
present(work))
then
1189 if (
size(work) < lwork)
then
1191 call errmgr%report_error(
"solve_qr_no_pivot_mtx", &
1192 "Incorrectly sized input array WORK, argument 5.", &
1193 la_array_size_error)
1196 wptr => work(1:lwork)
1198 allocate(wrk(lwork), stat = istat)
1199 if (istat /= 0)
then
1201 call errmgr%report_error(
"solve_qr_pivot_vec", &
1202 "Insufficient memory available.", &
1203 la_out_of_memory_error)
1214 if (abs(a(1,1)) == zero)
then
1224 call dlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1225 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1226 call dlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1227 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1228 if (smaxpr * rcond <= sminpr)
then
1230 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1231 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1233 wptr(ismin+rnk) = c1
1234 wptr(ismax+rnk) = c2
1247 w => wptr(rnk+1:lwork)
1248 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
1251 call mult_qr(.true., a, tau, b(1:m))
1254 call solve_triangular_system(.true., .false., .true., a(1:rnk,1:rnk), &
1256 if (n > rnk) b(rnk+1:n) = zero
1260 call mult_rz(.true., n - rnk, a(1:rnk,:), tau2, b(1:n), w)
1265 wptr(jpvt(i)) = b(i)
1274 module subroutine solve_qr_pivot_vec_cmplx(a, tau, jpvt, b, work, olwork, err)
1276 complex(real64),
intent(inout),
dimension(:,:) :: a
1277 complex(real64),
intent(in),
dimension(:) :: tau
1278 integer(int32),
intent(in),
dimension(:) :: jpvt
1279 complex(real64),
intent(inout),
dimension(:) :: b
1280 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
1281 integer(int32),
intent(out),
optional :: olwork
1282 class(errors),
intent(inout),
optional,
target :: err
1285 integer(int32),
parameter :: imin = 2
1286 integer(int32),
parameter :: imax = 1
1287 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
1288 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
1291 integer(int32) :: i, m, n, mn, lwork, ismin, ismax, rnk, maxmn, flag, &
1292 istat, lwork1, lwork2
1293 real(real64) :: rcond, smax, smin, smaxpr, sminpr
1294 complex(real64) :: s1, c1, s2, c2
1295 complex(real64),
pointer,
dimension(:) :: wptr, w, tau2
1296 complex(real64),
allocatable,
target,
dimension(:) :: wrk
1297 class(errors),
pointer :: errmgr
1298 type(errors),
target :: deferr
1299 character(len = :),
allocatable :: errmsg
1308 rcond = epsilon(rcond)
1309 if (
present(err))
then
1317 if (
size(tau) /= mn)
then
1319 else if (
size(jpvt) /= n)
then
1321 else if (
size(b) /= maxmn)
then
1326 allocate(
character(len = 256) :: errmsg)
1327 write(errmsg, 100)
"Input number ", flag, &
1328 " is not sized correctly."
1329 call errmgr%report_error(
"solve_qr_pivot_vec_cmplx", trim(errmsg), &
1330 la_array_size_error)
1335 call rz_factor(a(1:mn,:), a(1:mn,1), olwork = lwork1)
1336 call mult_rz(.true., n, a(1:mn,:), a(1:mn,1), b(1:n), olwork = lwork2)
1337 lwork = max(lwork1, lwork2, 2 * mn + 1) + mn
1338 if (
present(olwork))
then
1344 if (
present(work))
then
1345 if (
size(work) < lwork)
then
1347 call errmgr%report_error(
"solve_qr_no_pivot_mtx_cmplx", &
1348 "Incorrectly sized input array WORK, argument 5.", &
1349 la_array_size_error)
1352 wptr => work(1:lwork)
1354 allocate(wrk(lwork), stat = istat)
1355 if (istat /= 0)
then
1357 call errmgr%report_error(
"solve_qr_pivot_vec_cmplx", &
1358 "Insufficient memory available.", &
1359 la_out_of_memory_error)
1370 if (abs(a(1,1)) == zero)
then
1380 call zlaic1(imin, rnk, wptr(ismin:ismin+rnk-1), smin, &
1381 a(1:rnk-1,i), a(i,i), sminpr, s1, c1)
1382 call zlaic1(imax, rnk, wptr(ismax:ismax+rnk-1), smax, &
1383 a(1:rnk-1,i), a(i,i), smaxpr, s2, c2)
1384 if (smaxpr * rcond <= sminpr)
then
1386 wptr(ismin+i-1) = s1 * wptr(ismin+i-1)
1387 wptr(ismax+i-1) = s2 * wptr(ismax+i-1)
1389 wptr(ismin+rnk) = c1
1390 wptr(ismax+rnk) = c2
1403 w => wptr(rnk+1:lwork)
1404 if (rnk < n)
call rz_factor(a(1:rnk,:), tau2, w)
1407 call mult_qr(.true., a, tau, b(1:m))
1410 call solve_triangular_system(.true., .false., .true., a(1:rnk,1:rnk), &
1412 if (n > rnk) b(rnk+1:n) = zero
1416 call mult_rz(.true., n - rnk, a(1:rnk,:), tau2, b(1:n), w)
1421 wptr(jpvt(i)) = b(i)
1432 module subroutine solve_cholesky_mtx(upper, a, b, err)
1434 logical,
intent(in) :: upper
1435 real(real64),
intent(in),
dimension(:,:) :: a
1436 real(real64),
intent(inout),
dimension(:,:) :: b
1437 class(errors),
intent(inout),
optional,
target :: err
1441 integer(int32) :: n, nrhs, flag
1442 class(errors),
pointer :: errmgr
1443 type(errors),
target :: deferr
1444 character(len = :),
allocatable :: errmsg
1454 if (
present(err))
then
1462 if (
size(a, 2) /= n)
then
1464 else if (
size(b, 1) /= n)
then
1468 allocate(
character(len = 256) :: errmsg)
1469 write(errmsg, 100)
"Input number ", flag, &
1470 " is not sized correctly."
1471 call errmgr%report_error(
"solve_cholesky_mtx", trim(errmsg), &
1472 la_array_size_error)
1477 call dpotrs(uplo, n, nrhs, a, n, b, n, flag)
1484 module subroutine solve_cholesky_mtx_cmplx(upper, a, b, err)
1486 logical,
intent(in) :: upper
1487 complex(real64),
intent(in),
dimension(:,:) :: a
1488 complex(real64),
intent(inout),
dimension(:,:) :: b
1489 class(errors),
intent(inout),
optional,
target :: err
1493 integer(int32) :: n, nrhs, flag
1494 class(errors),
pointer :: errmgr
1495 type(errors),
target :: deferr
1496 character(len = :),
allocatable :: errmsg
1506 if (
present(err))
then
1514 if (
size(a, 2) /= n)
then
1516 else if (
size(b, 1) /= n)
then
1520 allocate(
character(len = 256) :: errmsg)
1521 write(errmsg, 100)
"Input number ", flag, &
1522 " is not sized correctly."
1523 call errmgr%report_error(
"solve_cholesky_mtx_cmplx", trim(errmsg), &
1524 la_array_size_error)
1529 call zpotrs(uplo, n, nrhs, a, n, b, n, flag)
1536 module subroutine solve_cholesky_vec(upper, a, b, err)
1538 logical,
intent(in) :: upper
1539 real(real64),
intent(in),
dimension(:,:) :: a
1540 real(real64),
intent(inout),
dimension(:) :: b
1541 class(errors),
intent(inout),
optional,
target :: err
1545 integer(int32) :: n, flag
1546 class(errors),
pointer :: errmgr
1547 type(errors),
target :: deferr
1548 character(len = :),
allocatable :: errmsg
1557 if (
present(err))
then
1565 if (
size(a, 2) /= n)
then
1567 else if (
size(b) /= n)
then
1571 allocate(
character(len = 256) :: errmsg)
1572 write(errmsg, 100)
"Input number ", flag, &
1573 " is not sized correctly."
1574 call errmgr%report_error(
"solve_cholesky_vec", trim(errmsg), &
1575 la_array_size_error)
1580 call dpotrs(uplo, n, 1, a, n, b, n, flag)
1587 module subroutine solve_cholesky_vec_cmplx(upper, a, b, err)
1589 logical,
intent(in) :: upper
1590 complex(real64),
intent(in),
dimension(:,:) :: a
1591 complex(real64),
intent(inout),
dimension(:) :: b
1592 class(errors),
intent(inout),
optional,
target :: err
1596 integer(int32) :: n, flag
1597 class(errors),
pointer :: errmgr
1598 type(errors),
target :: deferr
1599 character(len = :),
allocatable :: errmsg
1608 if (
present(err))
then
1616 if (
size(a, 2) /= n)
then
1618 else if (
size(b) /= n)
then
1622 allocate(
character(len = 256) :: errmsg)
1623 write(errmsg, 100)
"Input number ", flag, &
1624 " is not sized correctly."
1625 call errmgr%report_error(
"solve_cholesky_vec_cmplx", trim(errmsg), &
1626 la_array_size_error)
1631 call zpotrs(uplo, n, 1, a, n, b, n, flag)
1640 module subroutine mtx_inverse_dbl(a, iwork, work, olwork, err)
1642 real(real64),
intent(inout),
dimension(:,:) :: a
1643 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1644 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
1645 integer(int32),
intent(out),
optional :: olwork
1646 class(errors),
intent(inout),
optional,
target :: err
1649 integer(int32) :: n, liwork, lwork, istat, flag, itemp(1)
1650 integer(int32),
pointer,
dimension(:) :: iptr
1651 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1652 real(real64),
pointer,
dimension(:) :: wptr
1653 real(real64),
allocatable,
target,
dimension(:) :: wrk
1654 real(real64),
dimension(1) :: temp
1655 class(errors),
pointer :: errmgr
1656 type(errors),
target :: deferr
1661 if (
present(err))
then
1668 if (
size(a, 2) /= n)
then
1669 call errmgr%report_error(
"mtx_inverse", &
1670 "The matrix must be squre to invert.", la_array_size_error)
1675 call dgetri(n, a, n, itemp, temp, -1, flag)
1676 lwork = int(temp(1), int32)
1677 if (
present(olwork))
then
1683 if (
present(work))
then
1684 if (
size(work) < lwork)
then
1686 call errmgr%report_error(
"mtx_inverse_dbl", &
1687 "Incorrectly sized input array WORK, argument 3.", &
1688 la_array_size_error)
1691 wptr => work(1:lwork)
1693 allocate(wrk(lwork), stat = istat)
1694 if (istat /= 0)
then
1696 call errmgr%report_error(
"mtx_inverse_dbl", &
1697 "Insufficient memory available.", &
1698 la_out_of_memory_error)
1705 if (
present(iwork))
then
1706 if (
size(iwork) < liwork)
then
1708 call errmgr%report_error(
"mtx_inverse_dbl", &
1709 "Incorrectly sized input array IWORK, argument 2.", &
1710 la_array_size_error)
1713 iptr => iwork(1:liwork)
1715 allocate(iwrk(liwork), stat = istat)
1716 if (istat /= 0)
then
1718 call errmgr%report_error(
"mtx_inverse_dbl", &
1719 "Insufficient memory available.", &
1720 la_out_of_memory_error)
1727 call dgetrf(n, n, a, n, iptr, flag)
1730 call dgetri(n, a, n, iptr, wptr, lwork, flag)
1734 call errmgr%report_error(
"mtx_inverse_dbl", &
1735 "The matrix is singular; therefore, the inverse could " // &
1736 "not be computed.", la_singular_matrix_error)
1741 module subroutine mtx_inverse_cmplx(a, iwork, work, olwork, err)
1743 complex(real64),
intent(inout),
dimension(:,:) :: a
1744 integer(int32),
intent(out),
target,
optional,
dimension(:) :: iwork
1745 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
1746 integer(int32),
intent(out),
optional :: olwork
1747 class(errors),
intent(inout),
optional,
target :: err
1750 integer(int32) :: n, liwork, lwork, istat, flag, itemp(1)
1751 integer(int32),
pointer,
dimension(:) :: iptr
1752 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
1753 complex(real64),
pointer,
dimension(:) :: wptr
1754 complex(real64),
allocatable,
target,
dimension(:) :: wrk
1755 complex(real64),
dimension(1) :: temp
1756 class(errors),
pointer :: errmgr
1757 type(errors),
target :: deferr
1762 if (
present(err))
then
1769 if (
size(a, 2) /= n)
then
1770 call errmgr%report_error(
"mtx_inverse_cmplx", &
1771 "The matrix must be squre to invert.", la_array_size_error)
1776 call zgetri(n, a, n, itemp, temp, -1, flag)
1777 lwork = int(temp(1), int32)
1778 if (
present(olwork))
then
1784 if (
present(work))
then
1785 if (
size(work) < lwork)
then
1787 call errmgr%report_error(
"mtx_inverse_cmplx", &
1788 "Incorrectly sized input array WORK, argument 3.", &
1789 la_array_size_error)
1792 wptr => work(1:lwork)
1794 allocate(wrk(lwork), stat = istat)
1795 if (istat /= 0)
then
1797 call errmgr%report_error(
"mtx_inverse_cmplx", &
1798 "Insufficient memory available.", &
1799 la_out_of_memory_error)
1806 if (
present(iwork))
then
1807 if (
size(iwork) < liwork)
then
1809 call errmgr%report_error(
"mtx_inverse_cmplx", &
1810 "Incorrectly sized input array IWORK, argument 2.", &
1811 la_array_size_error)
1814 iptr => iwork(1:liwork)
1816 allocate(iwrk(liwork), stat = istat)
1817 if (istat /= 0)
then
1819 call errmgr%report_error(
"mtx_inverse_cmplx", &
1820 "Insufficient memory available.", &
1821 la_out_of_memory_error)
1828 call zgetrf(n, n, a, n, iptr, flag)
1831 call zgetri(n, a, n, iptr, wptr, lwork, flag)
1835 call errmgr%report_error(
"mtx_inverse_cmplx", &
1836 "The matrix is singular; therefore, the inverse could " // &
1837 "not be computed.", la_singular_matrix_error)
1842 module subroutine mtx_pinverse_dbl(a, ainv, tol, work, olwork, err)
1844 real(real64),
intent(inout),
dimension(:,:) :: a
1845 real(real64),
intent(out),
dimension(:,:) :: ainv
1846 real(real64),
intent(in),
optional :: tol
1847 real(real64),
intent(out),
target,
dimension(:),
optional :: work
1848 integer(int32),
intent(out),
optional :: olwork
1849 class(errors),
intent(inout),
optional,
target :: err
1853 function dlamch(cmach)
result(x)
1854 use,
intrinsic :: iso_fortran_env, only : real64
1855 character,
intent(in) :: cmach
1861 real(real64),
parameter :: zero = 0.0d0
1862 real(real64),
parameter :: one = 1.0d0
1865 integer(int32) :: i, m, n, mn, lwork, istat, flag, i1, i2a, i2b, i3a, &
1867 real(real64),
pointer,
dimension(:) :: s, wptr, w
1868 real(real64),
pointer,
dimension(:,:) :: u, vt
1869 real(real64),
allocatable,
target,
dimension(:) :: wrk
1870 real(real64),
dimension(1) :: temp
1871 real(real64) :: t, tref, tolcheck
1872 class(errors),
pointer :: errmgr
1873 type(errors),
target :: deferr
1874 character(len = :),
allocatable :: errmsg
1882 i2b = i2a + n * n - 1
1886 tolcheck = dlamch(
's')
1887 if (
present(err))
then
1894 if (
size(ainv, 1) /= n .or.
size(ainv, 2) /= m)
then
1895 allocate(
character(len = 256) :: errmsg)
1896 write(errmsg, 100) &
1897 "The output matrix AINV is not sized appropriately. " // &
1898 "It is expected to be ", n,
"-by-", m,
"."
1899 call errmgr%report_error(
"mtx_pinverse", errmsg, &
1900 la_array_size_error)
1905 call dgesvd(
'S',
'A', m, n, a, m, a(1:mn,:), a, m, a, n, temp, -1, flag)
1906 lwork = int(temp(1), int32)
1907 lwork = lwork + m * mn + n * n + mn
1908 if (
present(olwork))
then
1914 if (
present(work))
then
1915 if (
size(work) < lwork)
then
1917 call errmgr%report_error(
"mtx_pinverse", &
1918 "Incorrectly sized input array WORK, argument 4.", &
1919 la_array_size_error)
1922 wptr => work(1:lwork)
1924 allocate(wrk(lwork), stat = istat)
1925 if (istat /= 0)
then
1927 call errmgr%report_error(
"mtx_pinverse", &
1928 "Insufficient memory available.", &
1929 la_out_of_memory_error)
1934 u(1:m,1:mn) => wptr(1:i1)
1935 vt(1:n,1:n) => wptr(i2a:i2b)
1940 call dgesvd(
'S',
'A', m, n, a, m, s, u, m, vt, n, w,
size(w), flag)
1944 allocate(
character(len = 256) :: errmsg)
1945 write(errmsg, 101) flag,
" superdiagonals could not " // &
1946 "converge to zero as part of the QR iteration process."
1947 call errmgr%report_warning(
"mtx_pinverse", errmsg, &
1948 la_convergence_error)
1954 tref = max(m, n) * epsilon(t) * s(1)
1955 if (
present(tol))
then
1961 if (t < tolcheck)
then
1979 call recip_mult_array(s(i), vt(i,1:n))
1984 call mtx_mult(.true., .true., one, vt(1:mn,:), u, zero, ainv)
1987100
format(a, i0, a, i0, a)
1992 module subroutine mtx_pinverse_cmplx(a, ainv, tol, work, olwork, rwork, err)
1994 complex(real64),
intent(inout),
dimension(:,:) :: a
1995 complex(real64),
intent(out),
dimension(:,:) :: ainv
1996 real(real64),
intent(in),
optional :: tol
1997 complex(real64),
intent(out),
target,
dimension(:),
optional :: work
1998 integer(int32),
intent(out),
optional :: olwork
1999 real(real64),
intent(out),
target,
dimension(:),
optional :: rwork
2000 class(errors),
intent(inout),
optional,
target :: err
2004 function dlamch(cmach)
result(x)
2005 use,
intrinsic :: iso_fortran_env, only : real64
2006 character,
intent(in) :: cmach
2012 complex(real64),
parameter :: zero = (0.0d0, 0.0d0)
2013 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
2016 integer(int32) :: i, m, n, mn, lwork, istat, flag, i1, i2a, i2b, i3, &
2018 real(real64),
pointer,
dimension(:) :: s, rwptr, rw
2019 real(real64),
allocatable,
target,
dimension(:) :: rwrk
2020 complex(real64),
pointer,
dimension(:) :: wptr, w
2021 complex(real64),
pointer,
dimension(:,:) :: u, vt
2022 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2023 complex(real64) :: temp(1), val
2024 real(real64) :: t, tref, tolcheck, rtemp(1)
2025 class(errors),
pointer :: errmgr
2026 type(errors),
target :: deferr
2027 character(len = :),
allocatable :: errmsg
2036 i2b = i2a + n * n - 1
2038 tolcheck = dlamch(
's')
2039 if (
present(err))
then
2046 if (
size(ainv, 1) /= n .or.
size(ainv, 2) /= m)
then
2047 allocate(
character(len = 256) :: errmsg)
2048 write(errmsg, 100) &
2049 "The output matrix AINV is not sized appropriately. " // &
2050 "It is expected to be ", n,
"-by-", m,
"."
2051 call errmgr%report_error(
"mtx_pinverse_cmplx", errmsg, &
2052 la_array_size_error)
2057 call zgesvd(
'S',
'A', m, n, a, m, rtemp, a, m, a, n, temp, -1, &
2059 lwork = int(temp(1), int32)
2060 lwork = lwork + m * mn + n * n
2061 if (
present(olwork))
then
2067 if (
present(work))
then
2068 if (
size(work) < lwork)
then
2070 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2071 "Incorrectly sized input array WORK, argument 4.", &
2072 la_array_size_error)
2075 wptr => work(1:lwork)
2077 allocate(wrk(lwork), stat = istat)
2078 if (istat /= 0)
then
2080 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2081 "Insufficient memory available.", &
2082 la_out_of_memory_error)
2088 if (
present(rwork))
then
2089 if (
size(rwork) < lrwork)
then
2091 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2092 "Incorrectly sized input array RWORK, argument 6.", &
2093 la_array_size_error)
2096 rwptr => rwork(1:lrwork)
2098 allocate(rwrk(lrwork), stat = istat)
2099 if (istat /= 0)
then
2101 call errmgr%report_error(
"mtx_pinverse_cmplx", &
2102 "Insufficient memory available.", &
2103 la_out_of_memory_error)
2108 u(1:m,1:mn) => wptr(1:i1)
2109 vt(1:n,1:n) => wptr(i2a:i2b)
2112 rw => rwptr(mn+1:lrwork)
2115 call zgesvd(
'S',
'A', m, n, a, m, s, u, m, vt, n, w,
size(w), rw, flag)
2119 allocate(
character(len = 256) :: errmsg)
2120 write(errmsg, 101) flag,
" superdiagonals could not " // &
2121 "converge to zero as part of the QR iteration process."
2122 call errmgr%report_warning(
"mtx_pinverse_cmplx", errmsg, &
2123 la_convergence_error)
2129 tref = max(m, n) * epsilon(t) * s(1)
2130 if (
present(tol))
then
2136 if (t < tolcheck)
then
2155 vt(i,1:n) = conjg(vt(i,1:n)) / s(i)
2168 val = val + vt(k,i) * conjg(u(j,k))
2175100
format(a, i0, a, i0, a)
2182 module subroutine solve_least_squares_mtx(a, b, work, olwork, err)
2184 real(real64),
intent(inout),
dimension(:,:) :: a, b
2185 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2186 integer(int32),
intent(out),
optional :: olwork
2187 class(errors),
intent(inout),
optional,
target :: err
2190 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag
2191 real(real64),
pointer,
dimension(:) :: wptr
2192 real(real64),
allocatable,
target,
dimension(:) :: wrk
2193 real(real64),
dimension(1) :: temp
2194 class(errors),
pointer :: errmgr
2195 type(errors),
target :: deferr
2202 if (
present(err))
then
2209 if (
size(b, 1) /= maxmn)
then
2210 call errmgr%report_error(
"solve_least_squares_mtx", &
2211 "Input 2 is not sized correctly.", la_array_size_error)
2216 call dgels(
'N', m, n, nrhs, a, m, b, maxmn, temp, -1, flag)
2217 lwork = int(temp(1), int32)
2218 if (
present(olwork))
then
2224 if (
present(work))
then
2225 if (
size(work) < lwork)
then
2227 call errmgr%report_error(
"solve_least_squares_mtx", &
2228 "Incorrectly sized input array WORK, argument 3.", &
2229 la_array_size_error)
2232 wptr => work(1:lwork)
2234 allocate(wrk(lwork), stat = istat)
2235 if (istat /= 0)
then
2237 call errmgr%report_error(
"solve_least_squares_mtx", &
2238 "Insufficient memory available.", &
2239 la_out_of_memory_error)
2246 call dgels(
'N', m, n, nrhs, a, m, b, maxmn, wptr, lwork, flag)
2248 call errmgr%report_error(
"solve_least_squares_mtx", &
2249 "The supplied matrix is not of full rank; therefore, " // &
2250 "the solution could not be computed via this routine. " // &
2251 "Try a routine that utilizes column pivoting.", &
2252 la_invalid_operation_error)
2257 module subroutine solve_least_squares_mtx_cmplx(a, b, work, olwork, err)
2259 complex(real64),
intent(inout),
dimension(:,:) :: a, b
2260 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2261 integer(int32),
intent(out),
optional :: olwork
2262 class(errors),
intent(inout),
optional,
target :: err
2265 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag
2266 complex(real64),
pointer,
dimension(:) :: wptr
2267 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2268 complex(real64),
dimension(1) :: temp
2269 class(errors),
pointer :: errmgr
2270 type(errors),
target :: deferr
2277 if (
present(err))
then
2284 if (
size(b, 1) /= maxmn)
then
2285 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2286 "Input 2 is not sized correctly.", la_array_size_error)
2291 call zgels(
'N', m, n, nrhs, a, m, b, maxmn, temp, -1, flag)
2292 lwork = int(temp(1), int32)
2293 if (
present(olwork))
then
2299 if (
present(work))
then
2300 if (
size(work) < lwork)
then
2302 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2303 "Incorrectly sized input array WORK, argument 3.", &
2304 la_array_size_error)
2307 wptr => work(1:lwork)
2309 allocate(wrk(lwork), stat = istat)
2310 if (istat /= 0)
then
2312 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2313 "Insufficient memory available.", &
2314 la_out_of_memory_error)
2321 call zgels(
'N', m, n, nrhs, a, m, b, maxmn, wptr, lwork, flag)
2323 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2324 "The supplied matrix is not of full rank; therefore, " // &
2325 "the solution could not be computed via this routine. " // &
2326 "Try a routine that utilizes column pivoting.", &
2327 la_invalid_operation_error)
2332 module subroutine solve_least_squares_vec(a, b, work, olwork, err)
2334 real(real64),
intent(inout),
dimension(:,:) :: a
2335 real(real64),
intent(inout),
dimension(:) :: b
2336 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2337 integer(int32),
intent(out),
optional :: olwork
2338 class(errors),
intent(inout),
optional,
target :: err
2341 integer(int32) :: m, n, maxmn, lwork, istat, flag
2342 real(real64),
pointer,
dimension(:) :: wptr
2343 real(real64),
allocatable,
target,
dimension(:) :: wrk
2344 real(real64),
dimension(1) :: temp
2345 class(errors),
pointer :: errmgr
2346 type(errors),
target :: deferr
2352 if (
present(err))
then
2359 if (
size(b) /= maxmn)
then
2360 call errmgr%report_error(
"solve_least_squares_vec", &
2361 "Input 2 is not sized correctly.", la_array_size_error)
2366 call dgels(
'N', m, n, 1, a, m, b, maxmn, temp, -1, flag)
2367 lwork = int(temp(1), int32)
2368 if (
present(olwork))
then
2374 if (
present(work))
then
2375 if (
size(work) < lwork)
then
2377 call errmgr%report_error(
"solve_least_squares_vec", &
2378 "Incorrectly sized input array WORK, argument 3.", &
2379 la_array_size_error)
2382 wptr => work(1:lwork)
2384 allocate(wrk(lwork), stat = istat)
2385 if (istat /= 0)
then
2387 call errmgr%report_error(
"solve_least_squares_vec", &
2388 "Insufficient memory available.", &
2389 la_out_of_memory_error)
2396 call dgels(
'N', m, n, 1, a, m, b, maxmn, wptr, lwork, flag)
2398 call errmgr%report_error(
"solve_least_squares_mtx", &
2399 "The supplied matrix is not of full rank; therefore, " // &
2400 "the solution could not be computed via this routine. " // &
2401 "Try a routine that utilizes column pivoting.", &
2402 la_invalid_operation_error)
2407 module subroutine solve_least_squares_vec_cmplx(a, b, work, olwork, err)
2409 complex(real64),
intent(inout),
dimension(:,:) :: a
2410 complex(real64),
intent(inout),
dimension(:) :: b
2411 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2412 integer(int32),
intent(out),
optional :: olwork
2413 class(errors),
intent(inout),
optional,
target :: err
2416 integer(int32) :: m, n, maxmn, lwork, istat, flag
2417 complex(real64),
pointer,
dimension(:) :: wptr
2418 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2419 complex(real64),
dimension(1) :: temp
2420 class(errors),
pointer :: errmgr
2421 type(errors),
target :: deferr
2427 if (
present(err))
then
2434 if (
size(b) /= maxmn)
then
2435 call errmgr%report_error(
"solve_least_squares_vec_cmplx", &
2436 "Input 2 is not sized correctly.", la_array_size_error)
2441 call zgels(
'N', m, n, 1, a, m, b, maxmn, temp, -1, flag)
2442 lwork = int(temp(1), int32)
2443 if (
present(olwork))
then
2449 if (
present(work))
then
2450 if (
size(work) < lwork)
then
2452 call errmgr%report_error(
"solve_least_squares_vec_cmplx", &
2453 "Incorrectly sized input array WORK, argument 3.", &
2454 la_array_size_error)
2457 wptr => work(1:lwork)
2459 allocate(wrk(lwork), stat = istat)
2460 if (istat /= 0)
then
2462 call errmgr%report_error(
"solve_least_squares_vec_cmplx", &
2463 "Insufficient memory available.", &
2464 la_out_of_memory_error)
2471 call zgels(
'N', m, n, 1, a, m, b, maxmn, wptr, lwork, flag)
2473 call errmgr%report_error(
"solve_least_squares_mtx_cmplx", &
2474 "The supplied matrix is not of full rank; therefore, " // &
2475 "the solution could not be computed via this routine. " // &
2476 "Try a routine that utilizes column pivoting.", &
2477 la_invalid_operation_error)
2482 module subroutine solve_least_squares_mtx_pvt(a, b, ipvt, arnk, work, olwork, err)
2484 real(real64),
intent(inout),
dimension(:,:) :: a, b
2485 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2486 integer(int32),
intent(out),
optional :: arnk
2487 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2488 integer(int32),
intent(out),
optional :: olwork
2489 class(errors),
intent(inout),
optional,
target :: err
2492 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag, rnk
2493 real(real64),
pointer,
dimension(:) :: wptr
2494 real(real64),
allocatable,
target,
dimension(:) :: wrk
2495 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2496 integer(int32),
pointer,
dimension(:) :: iptr
2497 real(real64),
dimension(1) :: temp
2498 integer(int32),
dimension(1) :: itemp
2500 class(errors),
pointer :: errmgr
2501 type(errors),
target :: deferr
2502 character(len = :),
allocatable :: errmsg
2510 if (
present(arnk)) arnk = 0
2511 if (
present(err))
then
2519 if (
size(b, 1) /= maxmn)
then
2523 allocate(
character(len = 256) :: errmsg)
2524 write(errmsg, 100)
"Input number ", flag, &
2525 " is not sized correctly."
2526 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2527 trim(errmsg), la_array_size_error)
2532 call dgelsy(m, n, nrhs, a, m, b, maxmn, itemp, rc, rnk, temp, -1, flag)
2533 lwork = int(temp(1), int32)
2534 if (
present(olwork))
then
2540 if (
present(ipvt))
then
2541 if (
size(ipvt) < n)
then
2543 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2544 "Incorrectly sized pivot array, argument 3.", &
2545 la_array_size_error)
2550 allocate(iwrk(n), stat = istat)
2551 if (istat /= 0)
then
2553 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2554 "Insufficient memory available.", &
2555 la_out_of_memory_error)
2562 if (
present(work))
then
2563 if (
size(work) < lwork)
then
2565 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2566 "Incorrectly sized input array WORK, argument 5.", &
2567 la_array_size_error)
2570 wptr => work(1:lwork)
2572 allocate(wrk(lwork), stat = istat)
2573 if (istat /= 0)
then
2575 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2576 "Insufficient memory available.", &
2577 la_out_of_memory_error)
2584 call dgelsy(m, n, nrhs, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2586 if (
present(arnk)) arnk = rnk
2593 module subroutine solve_least_squares_mtx_pvt_cmplx(a, b, ipvt, arnk, &
2594 work, olwork, rwork, err)
2596 complex(real64),
intent(inout),
dimension(:,:) :: a, b
2597 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2598 integer(int32),
intent(out),
optional :: arnk
2599 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2600 integer(int32),
intent(out),
optional :: olwork
2601 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
2602 class(errors),
intent(inout),
optional,
target :: err
2605 integer(int32) :: m, n, maxmn, nrhs, lwork, istat, flag, rnk, lrwork
2606 complex(real64),
pointer,
dimension(:) :: wptr
2607 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2608 real(real64),
pointer,
dimension(:) :: rwptr
2609 real(real64),
allocatable,
target,
dimension(:) :: rwrk
2610 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2611 integer(int32),
pointer,
dimension(:) :: iptr
2612 complex(real64),
dimension(1) :: temp
2613 real(real64),
dimension(1) :: rtemp
2614 integer(int32),
dimension(1) :: itemp
2616 class(errors),
pointer :: errmgr
2617 type(errors),
target :: deferr
2618 character(len = :),
allocatable :: errmsg
2627 if (
present(arnk)) arnk = 0
2628 if (
present(err))
then
2636 if (
size(b, 1) /= maxmn)
then
2640 allocate(
character(len = 256) :: errmsg)
2641 write(errmsg, 100)
"Input number ", flag, &
2642 " is not sized correctly."
2643 call errmgr%report_error(
"solve_least_squares_mtx_pvt_cmplx", &
2644 trim(errmsg), la_array_size_error)
2649 call zgelsy(m, n, nrhs, a, m, b, maxmn, itemp, rc, rnk, temp, -1, &
2651 lwork = int(temp(1), int32)
2652 if (
present(olwork))
then
2658 if (
present(ipvt))
then
2659 if (
size(ipvt) < n)
then
2661 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2662 "Incorrectly sized pivot array, argument 3.", &
2663 la_array_size_error)
2668 allocate(iwrk(n), 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)
2680 if (
present(work))
then
2681 if (
size(work) < lwork)
then
2683 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2684 "Incorrectly sized input array WORK, argument 5.", &
2685 la_array_size_error)
2688 wptr => work(1:lwork)
2690 allocate(wrk(lwork), stat = istat)
2691 if (istat /= 0)
then
2693 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2694 "Insufficient memory available.", &
2695 la_out_of_memory_error)
2701 if (
present(rwork))
then
2702 if (
size(rwork) < lrwork)
then
2704 call errmgr%report_error(
"solve_least_squares_mtx_pvt_cmplx", &
2705 "Incorrectly sized input array RWORK, argument 7.", &
2706 la_array_size_error)
2709 rwptr => rwork(1:lrwork)
2711 allocate(rwrk(lrwork), stat = istat)
2712 if (istat /= 0)
then
2714 call errmgr%report_error(
"solve_least_squares_mtx_pvt_cmplx", &
2715 "Insufficient memory available.", &
2716 la_out_of_memory_error)
2723 call zgelsy(m, n, nrhs, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2725 if (
present(arnk)) arnk = rnk
2732 module subroutine solve_least_squares_vec_pvt(a, b, ipvt, arnk, work, olwork, err)
2734 real(real64),
intent(inout),
dimension(:,:) :: a
2735 real(real64),
intent(inout),
dimension(:) :: b
2736 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2737 integer(int32),
intent(out),
optional :: arnk
2738 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
2739 integer(int32),
intent(out),
optional :: olwork
2740 class(errors),
intent(inout),
optional,
target :: err
2743 integer(int32) :: m, n, maxmn, lwork, istat, flag, rnk
2744 real(real64),
pointer,
dimension(:) :: wptr
2745 real(real64),
allocatable,
target,
dimension(:) :: wrk
2746 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2747 integer(int32),
pointer,
dimension(:) :: iptr
2748 real(real64),
dimension(1) :: temp
2749 integer(int32),
dimension(1) :: itemp
2751 class(errors),
pointer :: errmgr
2752 type(errors),
target :: deferr
2753 character(len = :),
allocatable :: errmsg
2760 if (
present(arnk)) arnk = 0
2761 if (
present(err))
then
2769 if (
size(b, 1) /= maxmn)
then
2773 allocate(
character(len = 256) :: errmsg)
2774 write(errmsg, 100)
"Input number ", flag, &
2775 " is not sized correctly."
2776 call errmgr%report_error(
"solve_least_squares_vec_pvt", &
2777 trim(errmsg), la_array_size_error)
2782 call dgelsy(m, n, 1, a, m, b, maxmn, itemp, rc, rnk, temp, -1, flag)
2783 lwork = int(temp(1), int32)
2784 if (
present(olwork))
then
2790 if (
present(ipvt))
then
2791 if (
size(ipvt) < n)
then
2793 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2794 "Incorrectly sized pivot array, argument 3.", &
2795 la_array_size_error)
2800 allocate(iwrk(n), stat = istat)
2801 if (istat /= 0)
then
2803 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2804 "Insufficient memory available.", &
2805 la_out_of_memory_error)
2812 if (
present(work))
then
2813 if (
size(work) < lwork)
then
2815 call errmgr%report_error(
"solve_least_squares_vec_pvt", &
2816 "Incorrectly sized input array WORK, argument 5.", &
2817 la_array_size_error)
2820 wptr => work(1:lwork)
2822 allocate(wrk(lwork), stat = istat)
2823 if (istat /= 0)
then
2825 call errmgr%report_error(
"solve_least_squares_vec_pvt", &
2826 "Insufficient memory available.", &
2827 la_out_of_memory_error)
2834 call dgelsy(m, n, 1, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, flag)
2835 if (
present(arnk)) arnk = rnk
2842 module subroutine solve_least_squares_vec_pvt_cmplx(a, b, ipvt, arnk, &
2843 work, olwork, rwork, err)
2845 complex(real64),
intent(inout),
dimension(:,:) :: a
2846 complex(real64),
intent(inout),
dimension(:) :: b
2847 integer(int32),
intent(inout),
target,
optional,
dimension(:) :: ipvt
2848 integer(int32),
intent(out),
optional :: arnk
2849 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
2850 integer(int32),
intent(out),
optional :: olwork
2851 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork
2852 class(errors),
intent(inout),
optional,
target :: err
2855 integer(int32) :: m, n, maxmn, lwork, istat, flag, rnk
2856 complex(real64),
pointer,
dimension(:) :: wptr
2857 complex(real64),
allocatable,
target,
dimension(:) :: wrk
2858 real(real64),
pointer,
dimension(:) :: rwptr
2859 real(real64),
allocatable,
target,
dimension(:) :: rwrk
2860 integer(int32),
allocatable,
target,
dimension(:) :: iwrk
2861 integer(int32),
pointer,
dimension(:) :: iptr
2862 complex(real64),
dimension(1) :: temp
2863 real(real64),
dimension(1) :: rtemp
2864 integer(int32),
dimension(1) :: itemp
2866 class(errors),
pointer :: errmgr
2867 type(errors),
target :: deferr
2868 character(len = :),
allocatable :: errmsg
2876 if (
present(arnk)) arnk = 0
2877 if (
present(err))
then
2885 if (
size(b, 1) /= maxmn)
then
2889 allocate(
character(len = 256) :: errmsg)
2890 write(errmsg, 100)
"Input number ", flag, &
2891 " is not sized correctly."
2892 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2893 trim(errmsg), la_array_size_error)
2898 call zgelsy(m, n, 1, a, m, b, maxmn, itemp, rc, rnk, temp, -1, rtemp, &
2900 lwork = int(temp(1), int32)
2901 if (
present(olwork))
then
2907 if (
present(ipvt))
then
2908 if (
size(ipvt) < n)
then
2910 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2911 "Incorrectly sized pivot array, argument 3.", &
2912 la_array_size_error)
2917 allocate(iwrk(n), stat = istat)
2918 if (istat /= 0)
then
2920 call errmgr%report_error(
"solve_least_squares_mtx_pvt", &
2921 "Insufficient memory available.", &
2922 la_out_of_memory_error)
2929 if (
present(work))
then
2930 if (
size(work) < lwork)
then
2932 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2933 "Incorrectly sized input array WORK, argument 5.", &
2934 la_array_size_error)
2937 wptr => work(1:lwork)
2939 allocate(wrk(lwork), stat = istat)
2940 if (istat /= 0)
then
2942 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2943 "Insufficient memory available.", &
2944 la_out_of_memory_error)
2950 if (
present(rwork))
then
2951 if (
size(rwork) < lrwork)
then
2953 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2954 "Incorrectly sized input array RWORK, argument 7.", &
2955 la_array_size_error)
2958 rwptr => rwork(1:lrwork)
2960 allocate(rwrk(lrwork), stat = istat)
2961 if (istat /= 0)
then
2963 call errmgr%report_error(
"solve_least_squares_vec_pvt_cmplx", &
2964 "Insufficient memory available.", &
2965 la_out_of_memory_error)
2972 call zgelsy(m, n, 1, a, m, b, maxmn, iptr, rc, rnk, wptr, lwork, &
2974 if (
present(arnk)) arnk = rnk
2981 module subroutine solve_least_squares_mtx_svd(a, b, s, arnk, work, olwork, err)
2983 real(real64),
intent(inout),
dimension(:,:) :: a, b
2984 integer(int32),
intent(out),
optional :: arnk
2985 real(real64),
intent(out),
target,
optional,
dimension(:) :: work, s
2986 integer(int32),
intent(out),
optional :: olwork
2987 class(errors),
intent(inout),
optional,
target :: err
2990 integer(int32) :: m, n, nrhs, mn, maxmn, istat, flag, lwork, rnk
2991 real(real64),
pointer,
dimension(:) :: wptr, sptr
2992 real(real64),
allocatable,
target,
dimension(:) :: wrk, sing
2993 real(real64),
dimension(1) :: temp
2994 real(real64) :: rcond
2995 class(errors),
pointer :: errmgr
2996 type(errors),
target :: deferr
2997 character(len = :),
allocatable :: errmsg
3005 rcond = epsilon(rcond)
3006 if (
present(arnk)) arnk = 0
3007 if (
present(err))
then
3015 if (
size(b, 1) /= maxmn)
then
3020 allocate(
character(len = 256) :: errmsg)
3021 write(errmsg, 100)
"Input number ", flag, &
3022 " is not sized correctly."
3023 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3024 trim(errmsg), la_array_size_error)
3029 call dgelss(m, n, nrhs, a, m, b, maxmn, temp, rcond, rnk, temp, -1, &
3031 lwork = int(temp(1), int32)
3032 if (
present(olwork))
then
3038 if (
present(s))
then
3039 if (
size(s) < mn)
then
3041 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3042 "Incorrectly sized input array S, argument 3.", &
3043 la_array_size_error)
3048 allocate(sing(mn), stat = istat)
3049 if (istat /= 0)
then
3051 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3052 "Insufficient memory available.", &
3053 la_out_of_memory_error)
3059 if (
present(work))
then
3060 if (
size(work) < lwork)
then
3062 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3063 "Incorrectly sized input array WORK, argument 5.", &
3064 la_array_size_error)
3067 wptr => work(1:lwork)
3069 allocate(wrk(lwork), stat = istat)
3070 if (istat /= 0)
then
3072 call errmgr%report_error(
"solve_least_squares_mtx_svd", &
3073 "Insufficient memory available.", &
3074 la_out_of_memory_error)
3081 call dgelss(m, n, nrhs, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3083 if (
present(arnk)) arnk = rnk
3085 allocate(
character(len = 256) :: errmsg)
3086 write(errmsg, 101) flag,
" superdiagonals could not " // &
3087 "converge to zero as part of the QR iteration process."
3088 call errmgr%report_warning(
"solve_least_squares_mtx_svd", errmsg, &
3089 la_convergence_error)
3098 module subroutine solve_least_squares_mtx_svd_cmplx(a, b, s, arnk, work, &
3101 complex(real64),
intent(inout),
dimension(:,:) :: a, b
3102 integer(int32),
intent(out),
optional :: arnk
3103 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3104 real(real64),
intent(out),
target,
optional,
dimension(:) :: s, rwork
3105 integer(int32),
intent(out),
optional :: olwork
3106 class(errors),
intent(inout),
optional,
target :: err
3109 integer(int32) :: m, n, nrhs, mn, maxmn, istat, flag, lwork, rnk, lrwork
3110 complex(real64),
pointer,
dimension(:) :: wptr
3111 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3112 real(real64),
pointer,
dimension(:) :: rwptr, sptr
3113 real(real64),
allocatable,
target,
dimension(:) :: rwrk, sing
3114 complex(real64),
dimension(1) :: temp
3115 real(real64),
dimension(1) :: rtemp
3116 real(real64) :: rcond
3117 class(errors),
pointer :: errmgr
3118 type(errors),
target :: deferr
3119 character(len = :),
allocatable :: errmsg
3128 rcond = epsilon(rcond)
3129 if (
present(arnk)) arnk = 0
3130 if (
present(err))
then
3138 if (
size(b, 1) /= maxmn)
then
3143 allocate(
character(len = 256) :: errmsg)
3144 write(errmsg, 100)
"Input number ", flag, &
3145 " is not sized correctly."
3146 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3147 trim(errmsg), la_array_size_error)
3152 call zgelss(m, n, nrhs, a, m, b, maxmn, rtemp, rcond, rnk, temp, -1, &
3154 lwork = int(temp(1), int32)
3155 if (
present(olwork))
then
3161 if (
present(s))
then
3162 if (
size(s) < mn)
then
3164 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3165 "Incorrectly sized input array S, argument 3.", &
3166 la_array_size_error)
3171 allocate(sing(mn), stat = istat)
3172 if (istat /= 0)
then
3174 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3175 "Insufficient memory available.", &
3176 la_out_of_memory_error)
3182 if (
present(work))
then
3183 if (
size(work) < lwork)
then
3185 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3186 "Incorrectly sized input array WORK, argument 5.", &
3187 la_array_size_error)
3190 wptr => work(1:lwork)
3192 allocate(wrk(lwork), stat = istat)
3193 if (istat /= 0)
then
3195 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3196 "Insufficient memory available.", &
3197 la_out_of_memory_error)
3203 if (
present(rwork))
then
3204 if (
size(rwork) < lrwork)
then
3206 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3207 "Incorrectly sized input array RWORK, argument 7.", &
3208 la_array_size_error)
3211 rwptr => rwork(1:lrwork)
3213 allocate(rwrk(lrwork), stat = istat)
3214 if (istat /= 0)
then
3216 call errmgr%report_error(
"solve_least_squares_mtx_svd_cmplx", &
3217 "Insufficient memory available.", &
3218 la_out_of_memory_error)
3225 call zgelss(m, n, nrhs, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3227 if (
present(arnk)) arnk = rnk
3229 allocate(
character(len = 256) :: errmsg)
3230 write(errmsg, 101) flag,
" superdiagonals could not " // &
3231 "converge to zero as part of the QR iteration process."
3232 call errmgr%report_warning(
"solve_least_squares_mtx_svd_cmplx", &
3233 errmsg, la_convergence_error)
3242 module subroutine solve_least_squares_vec_svd(a, b, s, arnk, work, olwork, err)
3244 real(real64),
intent(inout),
dimension(:,:) :: a
3245 real(real64),
intent(inout),
dimension(:) :: b
3246 integer(int32),
intent(out),
optional :: arnk
3247 real(real64),
intent(out),
target,
optional,
dimension(:) :: work, s
3248 integer(int32),
intent(out),
optional :: olwork
3249 class(errors),
intent(inout),
optional,
target :: err
3252 integer(int32) :: m, n, mn, maxmn, istat, flag, lwork, rnk
3253 real(real64),
pointer,
dimension(:) :: wptr, sptr
3254 real(real64),
allocatable,
target,
dimension(:) :: wrk, sing
3255 real(real64),
dimension(1) :: temp
3256 real(real64) :: rcond
3257 class(errors),
pointer :: errmgr
3258 type(errors),
target :: deferr
3259 character(len = :),
allocatable :: errmsg
3266 rcond = epsilon(rcond)
3267 if (
present(arnk)) arnk = 0
3268 if (
present(err))
then
3276 if (
size(b) /= maxmn)
then
3281 allocate(
character(len = 256) :: errmsg)
3282 write(errmsg, 100)
"Input number ", flag, &
3283 " is not sized correctly."
3284 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3285 trim(errmsg), la_array_size_error)
3290 call dgelss(m, n, 1, a, m, b, maxmn, temp, rcond, rnk, temp, -1, flag)
3291 lwork = int(temp(1), int32)
3292 if (
present(olwork))
then
3298 if (
present(s))
then
3299 if (
size(s) < mn)
then
3301 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3302 "Incorrectly sized input array S, argument 3.", &
3303 la_array_size_error)
3308 allocate(sing(mn), stat = istat)
3309 if (istat /= 0)
then
3311 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3312 "Insufficient memory available.", &
3313 la_out_of_memory_error)
3319 if (
present(work))
then
3320 if (
size(work) < lwork)
then
3322 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3323 "Incorrectly sized input array WORK, argument 5.", &
3324 la_array_size_error)
3327 wptr => work(1:lwork)
3329 allocate(wrk(lwork), stat = istat)
3330 if (istat /= 0)
then
3332 call errmgr%report_error(
"solve_least_squares_vec_svd", &
3333 "Insufficient memory available.", &
3334 la_out_of_memory_error)
3341 call dgelss(m, n, 1, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3343 if (
present(arnk)) arnk = rnk
3345 allocate(
character(len = 256) :: errmsg)
3346 write(errmsg, 101) flag,
" superdiagonals could not " // &
3347 "converge to zero as part of the QR iteration process."
3348 call errmgr%report_warning(
"solve_least_squares_vec_svd", errmsg, &
3349 la_convergence_error)
3358 module subroutine solve_least_squares_vec_svd_cmplx(a, b, s, arnk, work, &
3361 complex(real64),
intent(inout),
dimension(:,:) :: a
3362 complex(real64),
intent(inout),
dimension(:) :: b
3363 integer(int32),
intent(out),
optional :: arnk
3364 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3365 real(real64),
intent(out),
target,
optional,
dimension(:) :: rwork, s
3366 integer(int32),
intent(out),
optional :: olwork
3367 class(errors),
intent(inout),
optional,
target :: err
3370 integer(int32) :: m, n, mn, maxmn, istat, flag, lwork, rnk, lrwork
3371 real(real64),
pointer,
dimension(:) :: rwptr, sptr
3372 real(real64),
allocatable,
target,
dimension(:) :: rwrk, sing
3373 complex(real64),
pointer,
dimension(:) :: wptr
3374 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3375 complex(real64),
dimension(1) :: temp
3376 real(real64),
dimension(1) :: rtemp
3377 real(real64) :: rcond
3378 class(errors),
pointer :: errmgr
3379 type(errors),
target :: deferr
3380 character(len = :),
allocatable :: errmsg
3388 rcond = epsilon(rcond)
3389 if (
present(arnk)) arnk = 0
3390 if (
present(err))
then
3398 if (
size(b) /= maxmn)
then
3403 allocate(
character(len = 256) :: errmsg)
3404 write(errmsg, 100)
"Input number ", flag, &
3405 " is not sized correctly."
3406 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3407 trim(errmsg), la_array_size_error)
3412 call zgelss(m, n, 1, a, m, b, maxmn, rtemp, rcond, rnk, temp, -1, &
3414 lwork = int(temp(1), int32)
3415 if (
present(olwork))
then
3421 if (
present(s))
then
3422 if (
size(s) < mn)
then
3424 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3425 "Incorrectly sized input array S, argument 3.", &
3426 la_array_size_error)
3431 allocate(sing(mn), stat = istat)
3432 if (istat /= 0)
then
3434 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3435 "Insufficient memory available.", &
3436 la_out_of_memory_error)
3442 if (
present(work))
then
3443 if (
size(work) < lwork)
then
3445 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3446 "Incorrectly sized input array WORK, argument 5.", &
3447 la_array_size_error)
3450 wptr => work(1:lwork)
3452 allocate(wrk(lwork), stat = istat)
3453 if (istat /= 0)
then
3455 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3456 "Insufficient memory available.", &
3457 la_out_of_memory_error)
3463 if (
present(rwork))
then
3464 if (
size(rwork) < lrwork)
then
3466 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3467 "Incorrectly sized input array RWORK, argument 7.", &
3468 la_array_size_error)
3471 rwptr => rwork(1:lrwork)
3473 allocate(rwrk(lrwork), stat = istat)
3474 if (istat /= 0)
then
3476 call errmgr%report_error(
"solve_least_squares_vec_svd_cmplx", &
3477 "Insufficient memory available.", &
3478 la_out_of_memory_error)
3485 call zgelss(m, n, 1, a, m, b, maxmn, sptr, rcond, rnk, wptr, lwork, &
3487 if (
present(arnk)) arnk = rnk
3489 allocate(
character(len = 256) :: errmsg)
3490 write(errmsg, 101) flag,
" superdiagonals could not " // &
3491 "converge to zero as part of the QR iteration process."
3492 call errmgr%report_warning(
"solve_least_squares_vec_svd_cmplx", &
3493 errmsg, la_convergence_error)
3504 module subroutine solve_lq_mtx(a, tau, b, work, olwork, err)
3506 real(real64),
intent(in),
dimension(:,:) :: a
3507 real(real64),
intent(in),
dimension(:) :: tau
3508 real(real64),
intent(inout),
dimension(:,:) :: b
3509 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
3510 integer(int32),
intent(out),
optional :: olwork
3511 class(errors),
intent(inout),
optional,
target :: err
3514 real(real64),
parameter :: one = 1.0d0
3517 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
3518 real(real64),
pointer,
dimension(:) :: wptr
3519 real(real64),
allocatable,
target,
dimension(:) :: wrk
3520 class(errors),
pointer :: errmgr
3521 type(errors),
target :: deferr
3522 character(len = :),
allocatable :: errmsg
3529 if (
present(err))
then
3539 else if (
size(tau) /= k)
then
3541 else if (
size(b, 1) /= n)
then
3547 allocate(
character(len = 256) :: errmsg)
3548 write(errmsg, 100)
"Input number ", flag, &
3549 " is not sized correctly."
3550 call errmgr%report_error(
"solve_lq_mtx", trim(errmsg), &
3551 la_array_size_error)
3556 call mult_lq(.true., .true., a, tau, b, olwork = lwork)
3558 if (
present(olwork))
then
3564 if (
present(work))
then
3565 if (
size(work) < lwork)
then
3567 call errmgr%report_error(
"solve_lq_mtx", &
3568 "Incorrectly sized input array WORK, argument 4.", &
3569 la_array_size_error)
3572 wptr => work(1:lwork)
3574 allocate(wrk(lwork), stat = istat)
3575 if (istat /= 0)
then
3577 call errmgr%report_error(
"solve_lq_mtx", &
3578 "Insufficient memory available.", &
3579 la_out_of_memory_error)
3587 call solve_triangular_system(.true., .false., .false., .true., one, &
3588 a(1:m,1:m), b(1:m,:), errmgr)
3589 if (errmgr%has_error_occurred())
return
3592 call mult_lq(.true., .true., a, tau, b, work = wptr, err = errmgr)
3593 if (errmgr%has_error_occurred())
return
3600 module subroutine solve_lq_mtx_cmplx(a, tau, b, work, olwork, err)
3602 complex(real64),
intent(in),
dimension(:,:) :: a
3603 complex(real64),
intent(in),
dimension(:) :: tau
3604 complex(real64),
intent(inout),
dimension(:,:) :: b
3605 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3606 integer(int32),
intent(out),
optional :: olwork
3607 class(errors),
intent(inout),
optional,
target :: err
3610 complex(real64),
parameter :: one = (1.0d0, 0.0d0)
3613 integer(int32) :: m, n, nrhs, k, lwork, flag, istat
3614 complex(real64),
pointer,
dimension(:) :: wptr
3615 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3616 class(errors),
pointer :: errmgr
3617 type(errors),
target :: deferr
3618 character(len = :),
allocatable :: errmsg
3625 if (
present(err))
then
3635 else if (
size(tau) /= k)
then
3637 else if (
size(b, 1) /= n)
then
3643 allocate(
character(len = 256) :: errmsg)
3644 write(errmsg, 100)
"Input number ", flag, &
3645 " is not sized correctly."
3646 call errmgr%report_error(
"solve_lq_mtx_cmplx", trim(errmsg), &
3647 la_array_size_error)
3652 call mult_lq(.true., .true., a, tau, b, olwork = lwork)
3654 if (
present(olwork))
then
3660 if (
present(work))
then
3661 if (
size(work) < lwork)
then
3663 call errmgr%report_error(
"solve_lq_mtx_cmplx", &
3664 "Incorrectly sized input array WORK, argument 4.", &
3665 la_array_size_error)
3668 wptr => work(1:lwork)
3670 allocate(wrk(lwork), stat = istat)
3671 if (istat /= 0)
then
3673 call errmgr%report_error(
"solve_lq_mtx_cmplx", &
3674 "Insufficient memory available.", &
3675 la_out_of_memory_error)
3683 call solve_triangular_system(.true., .false., .false., .true., one, &
3684 a(1:m,1:m), b(1:m,:), errmgr)
3685 if (errmgr%has_error_occurred())
return
3688 call mult_lq(.true., .true., a, tau, b, work = wptr, err = errmgr)
3689 if (errmgr%has_error_occurred())
return
3696 module subroutine solve_lq_vec(a, tau, b, work, olwork, err)
3698 real(real64),
intent(in),
dimension(:,:) :: a
3699 real(real64),
intent(in),
dimension(:) :: tau
3700 real(real64),
intent(inout),
dimension(:) :: b
3701 real(real64),
intent(out),
target,
optional,
dimension(:) :: work
3702 integer(int32),
intent(out),
optional :: olwork
3703 class(errors),
intent(inout),
optional,
target :: err
3706 integer(int32) :: m, n, k, lwork, flag, istat
3707 real(real64),
pointer,
dimension(:) :: wptr
3708 real(real64),
allocatable,
target,
dimension(:) :: wrk
3709 class(errors),
pointer :: errmgr
3710 type(errors),
target :: deferr
3711 character(len = :),
allocatable :: errmsg
3717 if (
present(err))
then
3727 else if (
size(tau) /= k)
then
3729 else if (
size(b) /= n)
then
3735 allocate(
character(len = 256) :: errmsg)
3736 write(errmsg, 100)
"Input number ", flag, &
3737 " is not sized correctly."
3738 call errmgr%report_error(
"solve_lq_vec", trim(errmsg), &
3739 la_array_size_error)
3744 call mult_lq(.true., a, tau, b, olwork = lwork)
3746 if (
present(olwork))
then
3752 if (
present(work))
then
3753 if (
size(work) < lwork)
then
3755 call errmgr%report_error(
"solve_lq_vec", &
3756 "Incorrectly sized input array WORK, argument 4.", &
3757 la_array_size_error)
3760 wptr => work(1:lwork)
3762 allocate(wrk(lwork), stat = istat)
3763 if (istat /= 0)
then
3765 call errmgr%report_error(
"solve_lq_vec", &
3766 "Insufficient memory available.", &
3767 la_out_of_memory_error)
3775 call solve_triangular_system(.false., .false., .true., a(1:m,1:m), &
3777 if (errmgr%has_error_occurred())
return
3780 call mult_lq(.true., a, tau, b, work = wptr, err = errmgr)
3781 if (errmgr%has_error_occurred())
return
3788 module subroutine solve_lq_vec_cmplx(a, tau, b, work, olwork, err)
3790 complex(real64),
intent(in),
dimension(:,:) :: a
3791 complex(real64),
intent(in),
dimension(:) :: tau
3792 complex(real64),
intent(inout),
dimension(:) :: b
3793 complex(real64),
intent(out),
target,
optional,
dimension(:) :: work
3794 integer(int32),
intent(out),
optional :: olwork
3795 class(errors),
intent(inout),
optional,
target :: err
3798 integer(int32) :: m, n, k, lwork, flag, istat
3799 complex(real64),
pointer,
dimension(:) :: wptr
3800 complex(real64),
allocatable,
target,
dimension(:) :: wrk
3801 class(errors),
pointer :: errmgr
3802 type(errors),
target :: deferr
3803 character(len = :),
allocatable :: errmsg
3809 if (
present(err))
then
3819 else if (
size(tau) /= k)
then
3821 else if (
size(b) /= n)
then
3827 allocate(
character(len = 256) :: errmsg)
3828 write(errmsg, 100)
"Input number ", flag, &
3829 " is not sized correctly."
3830 call errmgr%report_error(
"solve_lq_vec_cmplx", trim(errmsg), &
3831 la_array_size_error)
3836 call mult_lq(.true., a, tau, b, olwork = lwork)
3838 if (
present(olwork))
then
3844 if (
present(work))
then
3845 if (
size(work) < lwork)
then
3847 call errmgr%report_error(
"solve_lq_vec_cmplx", &
3848 "Incorrectly sized input array WORK, argument 4.", &
3849 la_array_size_error)
3852 wptr => work(1:lwork)
3854 allocate(wrk(lwork), stat = istat)
3855 if (istat /= 0)
then
3857 call errmgr%report_error(
"solve_lq_vec_cmplx", &
3858 "Insufficient memory available.", &
3859 la_out_of_memory_error)
3867 call solve_triangular_system(.false., .false., .true., a(1:m,1:m), &
3869 if (errmgr%has_error_occurred())
return
3872 call mult_lq(.true., a, tau, b, work = wptr, err = errmgr)
3873 if (errmgr%has_error_occurred())
return
Provides a set of common linear algebra routines.