8 use,
intrinsic :: iso_c_binding
9 use,
intrinsic :: iso_fortran_env, only : int32, real64
11 use ferror
, only : errors
12 use ferror_c_binding
, only : errorhandler, get_errorhandler
41 subroutine mtx_mult_c(transa, transb, m, n, k, alpha, a, lda, b, ldb, &
42 beta, c) bind(c, name =
"mtx_mult")
44 logical(c_bool),
intent(in),
value :: transa, transb
45 integer(int32),
intent(in),
value :: m, n, k, lda, ldb
46 real(real64),
intent(in),
value :: alpha, beta
47 real(real64),
intent(in) :: a(lda,*), b(ldb,*)
48 real(real64),
intent(inout) :: c(m,n)
62 call dgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
90 subroutine cmtx_mult_c(transa, transb, m, n, k, alpha, a, lda, b, ldb, &
91 beta, c) bind(c, name =
"cmtx_mult")
93 logical(c_bool),
intent(in),
value :: transa, transb
94 integer(int32),
intent(in),
value :: m, n, k, lda, ldb
95 real(real64),
intent(in),
value :: alpha, beta
96 complex(real64),
intent(in) :: a(lda,*), b(ldb,*)
97 complex(real64),
intent(inout) :: c(m,n)
100 complex(real64) :: calpha, cbeta
114 calpha = cmplx(alpha, 0.0d0, real64)
115 cbeta = cmplx(beta, 0.0d0, real64)
116 call zgemm(ta, tb, m, n, k, calpha, a, lda, b, ldb, cbeta, c, m)
132 subroutine diag_mtx_mult_left_c(m, n, k, alpha, a, b, beta, c) &
133 bind(c, name =
"diag_mtx_mult")
135 integer(int32),
intent(in),
value :: m, n, k
136 real(real64),
intent(in),
value :: alpha, beta
137 real(real64),
intent(in) :: a(min(m,k)), b(k,n)
138 real(real64),
intent(inout) :: c(m,n)
157 subroutine disg_mtx_mult_right_c(m, n, k, alpha, a, b, beta, c) &
158 bind(c, name =
"diag_mtx_rmult")
160 integer(int32),
intent(in),
value :: m, n, k
161 real(real64),
intent(in),
value :: alpha, beta
162 real(real64),
intent(in) :: a(m, k), b(min(k, n))
163 real(real64),
intent(inout) :: c(m, n)
182 subroutine diag_mtx_mult_cmplx_left_c(m, n, k, alpha, a, b, beta, c) &
183 bind(c, name =
"diag_mtx_mult_cmplx")
185 integer(int32),
intent(in),
value :: m, n, k
186 real(real64),
intent(in),
value :: alpha, beta
187 complex(real64),
intent(in) :: a(min(m, k))
188 real(real64),
intent(in) :: b(k, n)
189 complex(real64),
intent(inout) :: c(m, n)
208 subroutine diag_mtx_mult_cmplx_right_c(m, n, k, alpha, a, b, beta, c) &
209 bind(c, name =
"diag_mtx_rmult_cmplx")
211 integer(int32),
intent(in),
value :: m, n, k
212 real(real64),
intent(in),
value :: alpha, beta
213 real(real64),
intent(in) :: a(m, k)
214 complex(real64),
intent(in) :: b(min(k, n))
215 complex(real64),
intent(inout) :: c(m, n)
234 subroutine diag_cmtx_mult_left_c(m, n, k, alpha, a, b, beta, c) &
235 bind(c, name =
"diag_cmtx_mult")
237 integer(int32),
intent(in),
value :: m, n, k
238 real(real64),
intent(in),
value :: alpha, beta
239 complex(real64),
intent(in) :: a(min(m, k)), b(k, n)
240 complex(real64),
intent(inout) :: c(m, n)
259 subroutine diag_cmtx_mult_right_c(m, n, k, alpha, a, b, beta, c) &
260 bind(c, name =
"diag_cmtx_rmult")
262 integer(int32),
intent(in),
value :: m, n, k
263 real(real64),
intent(in),
value :: alpha, beta
264 complex(real64),
intent(in) :: a(m, k), b(min(k, n))
265 complex(real64),
intent(inout) :: c(m, n)
288 subroutine rank1_update_c(m, n, alpha, x, y, a) &
289 bind(c, name =
"rank1_update")
291 real(real64),
intent(in),
value :: alpha
292 integer(int32),
intent(in),
value :: m, n
293 real(real64),
intent(in) :: x(m), y(n)
294 real(real64),
intent(inout) :: a(m,n)
309 pure function trace_c(m, n, x) result(y) bind(C, name = "trace")
311 integer(int32),
intent(in),
value :: m, n
312 real(real64),
intent(in) :: x(m,n)
336 function mtx_rank_c(m, n, a, err) result(rnk) bind(C, name = "mtx_rank")
338 integer(int32),
intent(in),
value :: m, n
339 real(real64),
intent(inout) :: a(m,n)
340 type(errorhandler),
intent(inout) :: err
341 integer(int32) :: rnk
344 type(errors),
pointer :: eptr
347 call get_errorhandler(err, eptr)
348 if (
associated(eptr))
then 369 function det_c(n, a, err) result(x) bind(C, name = "det")
371 integer(int32),
intent(in),
value :: n
372 real(real64),
intent(inout) :: a(n,n)
373 type(errorhandler),
intent(inout) :: err
377 type(errors),
pointer :: eptr
380 call get_errorhandler(err, eptr)
381 if (
associated(eptr))
then 382 x =
det(a, err = eptr)
394 subroutine swap_c(n, x, y) bind(C, name = "swap")
396 integer(int32),
intent(in),
value :: n
397 real(real64),
intent(inout) :: x(n), y(n)
427 subroutine tri_mtx_mult_c(upper, n, alpha, a, beta, b, err) &
428 bind(c, name =
"tri_mtx_mult")
430 logical(c_bool),
intent(in),
value :: upper
431 integer(int32),
intent(in),
value :: n
432 real(real64),
intent(in),
value :: alpha, beta
433 real(real64),
intent(in) :: a(n,n)
434 real(real64),
intent(inout) :: b(n,n)
435 type(errorhandler),
intent(inout) :: err
438 type(errors),
pointer :: eptr
441 call get_errorhandler(err, eptr)
442 if (
associated(eptr))
then 443 call tri_mtx_mult(
logical(upper), alpha, a, beta, b, eptr)
470 subroutine lu_factor_c(m, n, a, ipvt, err) bind(C, name = "lu_factor")
472 integer(int32),
intent(in),
value :: m, n
473 real(real64),
intent(inout) :: a(m,n)
474 integer(int32),
intent(out) :: ipvt(min(m,n))
475 type(errorhandler),
intent(inout) :: err
478 type(errors),
pointer :: eptr
481 call get_errorhandler(err, eptr)
482 if (
associated(eptr))
then 516 subroutine form_lu_c(n, lu, ipvt, u, p) bind(C, name = "form_lu")
518 integer(int32),
intent(in),
value :: n
519 real(real64),
intent(inout) :: lu(n,n)
520 integer(int32),
intent(in) :: ipvt(n)
521 real(real64),
intent(out) :: u(n,n), p(n,n)
548 subroutine qr_factor_c(m, n, a, tau, err) bind(C, name = "qr_factor")
550 integer(int32),
intent(in),
value :: m, n
551 real(real64),
intent(inout) :: a(m,n)
552 real(real64),
intent(out) :: tau(min(m,n))
553 type(errorhandler),
intent(inout) :: err
556 type(errors),
pointer :: eptr
559 call get_errorhandler(err, eptr)
560 if (
associated(eptr))
then 592 subroutine qr_factor_pivot_c(m, n, a, tau, jpvt, err) &
593 bind(c, name =
"qr_factor_pivot")
595 integer(int32),
intent(in),
value :: m, n
596 real(real64),
intent(inout) :: a(m,n)
597 real(real64),
intent(out) :: tau(min(m,n))
598 integer(int32),
intent(inout) :: jpvt(n)
599 type(errorhandler),
intent(inout) :: err
602 type(errors),
pointer :: eptr
605 call get_errorhandler(err, eptr)
606 if (
associated(eptr))
then 638 subroutine form_qr_c(m, n, r, tau, q, err) bind(C, name = "form_qr")
640 integer(int32),
intent(in),
value :: m, n
641 real(real64),
intent(inout) :: r(m,n)
642 real(real64),
intent(in) :: tau(min(m,n))
643 real(real64),
intent(out) :: q(m,m)
644 type(errorhandler),
intent(inout) :: err
647 type(errors),
pointer :: eptr
650 call get_errorhandler(err, eptr)
651 if (
associated(eptr))
then 652 call form_qr(r, tau, q, err = eptr)
686 subroutine form_qr_pivot_c(m, n, r, tau, pvt, q, p, err) &
687 bind(c, name =
"form_qr_pivot")
689 integer(int32),
intent(in),
value :: m, n
690 real(real64),
intent(inout) :: r(m,n)
691 real(real64),
intent(in) :: tau(min(m,n))
692 integer(int32),
intent(in) :: pvt(n)
693 real(real64),
intent(out) :: q(m,m), p(n,n)
694 type(errorhandler),
intent(inout) :: err
697 type(errors),
pointer :: eptr
700 call get_errorhandler(err, eptr)
701 if (
associated(eptr))
then 702 call form_qr(r, tau, pvt, q, p, err = eptr)
704 call form_qr(r, tau, pvt, q, p)
731 subroutine mult_qr_c(trans, m, n, q, tau, c, err) &
732 bind(c, name =
"mult_qr")
734 logical(c_bool),
intent(in),
value :: trans
735 integer(int32),
intent(in),
value :: m, n
736 real(real64),
intent(inout) :: q(m,m), c(m,n)
737 real(real64),
intent(in) :: tau(min(m,n))
738 type(errorhandler),
intent(inout) :: err
741 type(errors),
pointer :: eptr
744 call get_errorhandler(err, eptr)
745 if (
associated(eptr))
then 746 call mult_qr(.true.,
logical(trans), q, tau, c, err = eptr)
748 call mult_qr(.true.,
logical(trans), q, tau, c)
772 subroutine qr_rank1_update_c(m, n, q, r, u, v, err) &
773 bind(c, name =
"qr_rank1_update")
775 integer(int32),
intent(in),
value :: m, n
776 real(real64),
intent(inout) :: q(m,m), r(m,n), u(m), v(n)
777 type(errorhandler),
intent(inout) :: err
780 type(errors),
pointer :: eptr
783 call get_errorhandler(err, eptr)
784 if (
associated(eptr))
then 808 subroutine cholesky_factor_c(n, a, upper, err) &
809 bind(c, name =
"cholesky_factor")
811 integer(int32),
intent(in),
value :: n
812 real(real64),
intent(inout) :: a(n,n)
813 logical(c_bool),
intent(in),
value :: upper
814 type(errorhandler),
intent(inout) :: err
817 type(errors),
pointer :: eptr
820 call get_errorhandler(err, eptr)
821 if (
associated(eptr))
then 843 subroutine cholesky_rank1_update_c(n, r, u, err) &
844 bind(c, name =
"cholesky_rank1_update")
846 integer(int32),
intent(in),
value :: n
847 real(real64),
intent(inout) :: r(n,n), u(n)
848 type(errorhandler),
intent(inout) :: err
851 type(errors),
pointer :: eptr
854 call get_errorhandler(err, eptr)
855 if (
associated(eptr))
then 882 subroutine cholesky_rank1_downdate_c(n, r, u, err) &
883 bind(c, name =
"cholesky_rank1_downdate")
885 integer(int32),
intent(in),
value :: n
886 real(real64),
intent(inout) :: r(n,n), u(n)
887 type(errorhandler),
intent(inout) :: err
890 type(errors),
pointer :: eptr
893 call get_errorhandler(err, eptr)
894 if (
associated(eptr))
then 922 subroutine rz_factor_c(m, n, a, tau, err) bind(C, name = "rz_factor")
924 integer(int32),
intent(in),
value :: m, n
925 real(real64),
intent(inout) :: a(m,n)
926 real(real64),
intent(out) :: tau(m)
927 type(errorhandler),
intent(inout) :: err
930 type(errors),
pointer :: eptr
933 call get_errorhandler(err, eptr)
934 if (
associated(eptr))
then 963 subroutine mult_rz_c(trans, m, n, l, a, tau, c, err) &
964 bind(c, name =
"mult_rz")
966 logical(c_bool),
intent(in),
value :: trans
967 integer(int32),
intent(in),
value :: m, n, l
968 real(real64),
intent(inout) :: a(m,m), c(m,n)
969 real(real64),
intent(in) :: tau(m)
970 type(errorhandler),
intent(inout) :: err
973 type(errors),
pointer :: eptr
976 call get_errorhandler(err, eptr)
977 if (
associated(eptr))
then 978 call mult_rz(.true.,
logical(trans), l, a, tau, c, err = eptr)
980 call mult_rz(.true.,
logical(trans), l, a, tau, c)
1011 subroutine svd_c(m, n, a, s, u, vt, err) bind(C, name = "svd")
1013 integer(int32),
intent(in),
value :: m, n
1014 real(real64),
intent(inout) :: a(m,n)
1015 real(real64),
intent(out) :: s(min(m,n)), u(m,m), vt(n,n)
1016 type(errorhandler),
intent(inout) :: err
1019 type(errors),
pointer :: eptr
1022 call get_errorhandler(err, eptr)
1023 if (
associated(eptr))
then 1024 call svd(a, s, u, vt, err = eptr)
1026 call svd(a, s, u, vt)
1050 subroutine solve_tri_mtx_c(upper, trans, nounit, n, nrhs, alpha, a, b) &
1051 bind(c, name =
"solve_triangular_system")
1053 logical(c_bool),
intent(in),
value :: upper, trans, nounit
1054 integer(int32),
intent(in),
value :: n, nrhs
1055 real(real64),
intent(in),
value :: alpha
1056 real(real64),
intent(in) :: a(n,n)
1057 real(real64),
intent(inout) :: b(n,nrhs)
1061 logical(nounit), alpha, a, b)
1074 subroutine solve_lu_c(n, nrhs, a, ipvt, b) bind(C, name = "solve_lu")
1076 integer(int32),
intent(in),
value :: n, nrhs
1077 real(real64),
intent(in) :: a(n,n)
1078 integer(int32),
intent(in) :: ipvt(n)
1079 real(real64),
intent(inout) :: b(n,nrhs)
1106 subroutine solve_qr_c(m, n, nrhs, a, tau, b, err) &
1107 bind(c, name =
"solve_qr")
1109 integer(int32),
intent(in),
value :: m, n, nrhs
1110 real(real64),
intent(inout) :: a(m,n), b(m,nrhs)
1111 real(real64),
intent(in) :: tau(min(m,n))
1112 type(errorhandler),
intent(inout) :: err
1115 type(errors),
pointer :: eptr
1118 call get_errorhandler(err, eptr)
1119 if (
associated(eptr))
then 1120 call solve_qr(a, tau, b, err = eptr)
1151 subroutine solve_qr_pivot_c(m, n, nrhs, a, tau, jpvt, b, err) &
1152 bind(c, name =
"solve_qr_pivot")
1154 integer(int32),
intent(in),
value :: m, n, nrhs
1155 real(real64),
intent(inout) :: a(m,n), b(max(m,n),nrhs)
1156 real(real64),
intent(in) :: tau(min(m,n))
1157 integer(int32),
intent(in) :: jpvt(n)
1158 type(errorhandler),
intent(inout) :: err
1161 type(errors),
pointer :: eptr
1164 call get_errorhandler(err, eptr)
1165 if (
associated(eptr))
then 1166 call solve_qr(a, tau, jpvt, b, err = eptr)
1184 subroutine solve_cholesky_c(upper, n, nrhs, a, b) &
1185 bind(c, name =
"solve_cholesky")
1187 logical(c_bool),
intent(in),
value :: upper
1188 integer(int32),
intent(in),
value :: n, nrhs
1189 real(real64),
intent(in) :: a(n,n)
1190 real(real64),
intent(inout) :: b(n,nrhs)
1209 subroutine mtx_inverse_c(n, a, err) bind(C, name = "mtx_inverse")
1211 integer(int32),
intent(in),
value :: n
1212 real(real64),
intent(inout) :: a(n,n)
1213 type(errorhandler),
intent(inout) :: err
1216 type(errors),
pointer :: eptr
1219 call get_errorhandler(err, eptr)
1220 if (
associated(eptr))
then 1245 subroutine mtx_pinverse_c(m, n, a, ainv, err) &
1246 bind(c, name =
"mtx_pinverse")
1248 integer(int32),
intent(in),
value :: m, n
1249 real(real64),
intent(inout) :: a(m,n)
1250 real(real64),
intent(out) :: ainv(n,m)
1251 type(errorhandler),
intent(inout) :: err
1254 type(errors),
pointer :: eptr
1257 call get_errorhandler(err, eptr)
1258 if (
associated(eptr))
then 1286 subroutine solve_least_squares_c(m, n, nrhs, a, b, err) &
1287 bind(c, name =
"solve_least_squares")
1289 integer(int32),
intent(in),
value :: m, n, nrhs
1290 real(real64),
intent(inout) :: a(m, n), b(max(m,n), nrhs)
1291 type(errorhandler),
intent(inout) :: err
1294 type(errors),
pointer :: eptr
1297 call get_errorhandler(err, eptr)
1298 if (
associated(eptr))
then 1328 subroutine eigen_symm_c(n, vecs, a, vals, err) bind(C, name = "eigen_symm")
1330 integer(int32),
intent(in),
value :: n
1331 logical(c_bool),
intent(in),
value :: vecs
1332 real(real64),
intent(inout) :: a(n,n)
1333 real(real64),
intent(out) :: vals(n)
1334 type(errorhandler),
intent(inout) :: err
1337 type(errors),
pointer :: eptr
1340 call get_errorhandler(err, eptr)
1341 if (
associated(eptr))
then 1342 call eigen(
logical(vecs), a, vals, err = eptr)
1344 call eigen(
logical(vecs), a, vals)
1366 subroutine eigen_asymm_c(n, a, vals, vecs, err) &
1367 bind(c, name =
"eigen_asymm")
1369 integer(int32),
intent(in),
value :: n
1370 real(real64),
intent(inout) :: a(n,n)
1371 complex(real64),
intent(out) :: vals(n), vecs(n,n)
1372 type(errorhandler),
intent(inout) :: err
1375 type(errors),
pointer :: eptr
1378 call get_errorhandler(err, eptr)
1379 if (
associated(eptr))
then 1380 call eigen(a, vals, vecs, err = eptr)
1382 call eigen(a, vals, vecs)
1415 subroutine eigen_gen_c(n, a, b, alpha, beta, vecs, err) &
1416 bind(c, name =
"eigen_gen")
1418 integer(int32),
intent(in),
value :: n
1419 real(real64),
intent(inout) :: a(n, n), b(n, n)
1420 complex(real64),
intent(out) :: alpha(n), vecs(n,n)
1421 real(real64),
intent(out) :: beta(n)
1422 type(errorhandler),
intent(inout) :: err
1425 type(errors),
pointer :: eptr
1428 call get_errorhandler(err, eptr)
1429 if (
associated(eptr))
then 1430 call eigen(a, b, alpha, beta = beta, vecs = vecs, err = eptr)
1432 call eigen(a, b, alpha, beta = beta, vecs = vecs)
1454 subroutine sort_dbl_ind_c(ascend, n, x, ind) bind(C, name = "sort_dbl")
1456 logical(c_bool),
intent(in),
value :: ascend
1457 integer(int32),
intent(in),
value :: n
1458 real(real64),
intent(inout) :: x(n)
1459 type(c_ptr),
intent(in),
value :: ind
1462 integer(int32),
pointer,
dimension(:) :: ptr
1465 if (c_associated(ind))
then 1466 call c_f_pointer(ind, ptr, [n])
1467 call sort(x, ptr,
logical(ascend))
1469 call sort(x,
logical(ascend))
1490 subroutine sort_cmplx_ind_c(ascend, n, x, ind) bind(C, name = "sort_cmplx")
1492 logical(c_bool),
intent(in),
value :: ascend
1493 integer(int32),
intent(in),
value :: n
1494 complex(real64),
intent(inout) :: x(n)
1495 type(c_ptr),
intent(in),
value :: ind
1498 integer(int32),
pointer,
dimension(:) :: ptr
1501 if (c_associated(ind))
then 1502 call c_f_pointer(ind, ptr, [n])
1503 call sort(x, ptr,
logical(ascend))
1505 call sort(x,
logical(ascend))
1521 subroutine sort_eigen_cmplx_c(ascend, n, vals, vecs) &
1522 bind(c, name =
"sort_eigen_cmplx")
1524 logical(c_bool),
intent(in),
value :: ascend
1525 integer(int32),
intent(in),
value :: n
1526 complex(real64),
intent(inout) :: vals(n), vecs(n,n)
1529 call sort(vals, vecs,
logical(ascend))
1544 subroutine sort_eigen_dbl_c(ascend, n, vals, vecs) &
1545 bind(c, name =
"sort_eigen_dbl")
1547 logical(c_bool),
intent(in),
value :: ascend
1548 integer(int32),
intent(in),
value :: n
1549 real(real64),
intent(inout) :: vals(n), vecs(n,n)
1552 call sort(vals, vecs,
logical(ascend))
Computes the QR factorization of an M-by-N matrix.
Factors an upper trapezoidal matrix by means of orthogonal transformations such that A = R * Z = (R 0...
Solves the overdetermined or underdetermined system (A*X = B) of M equations of N unknowns...
Multiplies a diagonal matrix with another matrix or array.
Computes the triangular matrix operation: B = alpha * A**T * A + beta * B, or B = alpha * A * A**T + ...
Solves a system of Cholesky factored equations.
Performs the rank-1 update to matrix A such that: A = alpha * X * Y**T + A, where A is an M-by-N matr...
Computes the singular value decomposition of a matrix A. The SVD is defined as: A = U * S * V**T...
Multiplies a general matrix by the orthogonal matrix Q from a QR factorization.
Solves a system of M QR-factored equations of N unknowns.
Computes the eigenvalues, and optionally the eigenvectors, of a matrix.
Computes the inverse of a square matrix.
Swaps the contents of two arrays.
Computes the LU factorization of an M-by-N matrix.
Solves a system of LU-factored equations.
Computes the determinant of a square matrix.
Computes the Moore-Penrose pseudo-inverse of a M-by-N matrix using the singular value decomposition o...
Computes the rank 1 downdate to a Cholesky factored matrix (upper triangular).
Computes the trace of a matrix (the sum of the main diagonal elements).
Solves a triangular system of equations.
Multiplies a general matrix by the orthogonal matrix Z from an RZ factorization.
Computes the rank of a matrix.
Computes the rank 1 update to a Cholesky factored matrix (upper triangular).
Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where A = Q * R...
Computes the Cholesky factorization of a symmetric, positive definite matrix.