linalg  1.5.0
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
linalg_c_binding.f90
1 ! linalg_c_binding.f90
2 
8  use, intrinsic :: iso_c_binding
9  use, intrinsic :: iso_fortran_env, only : int32, real64
10  use linalg_core
11  use ferror, only : errors
12  use ferror_c_binding, only : errorhandler, get_errorhandler
13 contains
14 ! ******************************************************************************
15 ! LINALG_CORE ROUTINES
16 ! ------------------------------------------------------------------------------
41  subroutine mtx_mult_c(transa, transb, m, n, k, alpha, a, lda, b, ldb, &
42  beta, c) bind(c, name = "mtx_mult")
43  ! Arguments
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)
49 
50  ! Process
51  character :: ta, tb
52  if (transa) then
53  ta = 'T'
54  else
55  ta = 'N'
56  end if
57  if (transb) then
58  tb = 'T'
59  else
60  tb = 'N'
61  end if
62  call dgemm(ta, tb, m, n, k, alpha, a, lda, b, ldb, beta, c, m)
63  end subroutine
64 
65 ! ------------------------------------------------------------------------------
90  subroutine cmtx_mult_c(transa, transb, m, n, k, alpha, a, lda, b, ldb, &
91  beta, c) bind(c, name = "cmtx_mult")
92  ! Arguments
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)
98 
99  ! Local Variables
100  complex(real64) :: calpha, cbeta
101 
102  ! Process
103  character :: ta, tb
104  if (transa) then
105  ta = 'T'
106  else
107  ta = 'N'
108  end if
109  if (transb) then
110  tb = 'T'
111  else
112  tb = 'N'
113  end if
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)
117  end subroutine
118 
119 ! ------------------------------------------------------------------------------
132  subroutine diag_mtx_mult_left_c(m, n, k, alpha, a, b, beta, c) &
133  bind(c, name = "diag_mtx_mult")
134  ! Arguments
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)
139 
140  ! Process
141  call diag_mtx_mult(.true., .false., alpha, a, b, beta, c)
142  end subroutine
143 
144 ! ------------------------------------------------------------------------------
157  subroutine disg_mtx_mult_right_c(m, n, k, alpha, a, b, beta, c) &
158  bind(c, name = "diag_mtx_rmult")
159  ! Arguments
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)
164 
165  ! Process
166  call diag_mtx_mult(.false., .false., alpha, b, a, beta, c)
167  end subroutine
168 
169 ! ------------------------------------------------------------------------------
182  subroutine diag_mtx_mult_cmplx_left_c(m, n, k, alpha, a, b, beta, c) &
183  bind(c, name = "diag_mtx_mult_cmplx")
184  ! Arguments
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)
190 
191  ! Process
192  call diag_mtx_mult(.true., .false., alpha, a, b, beta, c)
193  end subroutine
194 
195 ! ------------------------------------------------------------------------------
208  subroutine diag_mtx_mult_cmplx_right_c(m, n, k, alpha, a, b, beta, c) &
209  bind(c, name = "diag_mtx_rmult_cmplx")
210  ! Arguments
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)
216 
217  ! Process
218  call diag_mtx_mult(.false., .false., alpha, b, a, beta, c)
219  end subroutine
220 
221 ! ------------------------------------------------------------------------------
234  subroutine diag_cmtx_mult_left_c(m, n, k, alpha, a, b, beta, c) &
235  bind(c, name = "diag_cmtx_mult")
236  ! Arguments
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)
241 
242  ! Process
243  call diag_mtx_mult(.true., .false., alpha, a, b, beta, c)
244  end subroutine
245 
246 ! ------------------------------------------------------------------------------
259  subroutine diag_cmtx_mult_right_c(m, n, k, alpha, a, b, beta, c) &
260  bind(c, name = "diag_cmtx_rmult")
261  ! Arguments
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)
266 
267  ! Process
268  call diag_mtx_mult(.false., .false., alpha, b, a, beta, c)
269  end subroutine
270 
271 ! ------------------------------------------------------------------------------
288  subroutine rank1_update_c(m, n, alpha, x, y, a) &
289  bind(c, name = "rank1_update")
290  ! Arguments
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)
295 
296  ! Process
297  call rank1_update(alpha, x, y, a)
298  end subroutine
299 
300 ! ------------------------------------------------------------------------------
309  pure function trace_c(m, n, x) result(y) bind(C, name = "trace")
310  ! Arguments
311  integer(int32), intent(in), value :: m, n
312  real(real64), intent(in) :: x(m,n)
313  real(real64) :: y
314 
315  ! Process
316  y = trace(x)
317  end function
318 
319 ! ------------------------------------------------------------------------------
336  function mtx_rank_c(m, n, a, err) result(rnk) bind(C, name = "mtx_rank")
337  ! Arguments
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
342 
343  ! Local Variables
344  type(errors), pointer :: eptr
345 
346  ! Process
347  call get_errorhandler(err, eptr)
348  if (associated(eptr)) then
349  rnk = mtx_rank(a, err = eptr)
350  else
351  rnk = mtx_rank(a)
352  end if
353  end function
354 
355 ! ------------------------------------------------------------------------------
369  function det_c(n, a, err) result(x) bind(C, name = "det")
370  ! Arguments
371  integer(int32), intent(in), value :: n
372  real(real64), intent(inout) :: a(n,n)
373  type(errorhandler), intent(inout) :: err
374  real(real64) :: x
375 
376  ! Local Variables
377  type(errors), pointer :: eptr
378 
379  ! Process
380  call get_errorhandler(err, eptr)
381  if (associated(eptr)) then
382  x = det(a, err = eptr)
383  else
384  x = det(a)
385  end if
386  end function
387 
388 ! ------------------------------------------------------------------------------
394  subroutine swap_c(n, x, y) bind(C, name = "swap")
395  ! Arguments
396  integer(int32), intent(in), value :: n
397  real(real64), intent(inout) :: x(n), y(n)
398 
399  ! Process
400  call swap(x, y)
401  end subroutine
402 
403 ! ------------------------------------------------------------------------------
404  !@ brief Computes the triangular matrix operation:
405  !! B = alpha * A**T * A + beta * B, or B = alpha * A * A**T + beta * B,
406  !! where A is a triangular matrix.
407  !!
408  !! @param[in] upper Set to true if matrix A is upper triangular, and
409  !! B = alpha * A**T * A + beta * B is to be calculated; else, set to false
410  !! if A is lower triangular, and B = alpha * A * A**T + beta * B is to
411  !! be computed.
412  !! @param[in] n The size of the matrix.
413  !! @param[in] alpha A scalar multiplier.
414  !! @param[in] a The N-by-N triangular matrix. Notice, if @p upper is true
415  !! only the upper triangular portion of this matrix is referenced; else,
416  !! if @p upper is false, only the lower triangular portion of this matrix
417  !! is referenced.
418  !! @param[in] beta A scalar multiplier.
419  !! @param[in,out] b On input, the N-by-N matrix B. On output, the N-by-N
420  !! solution matrix.
421  !! @param[in,out] err The errorhandler object. If no error handling is
422  !! desired, simply pass NULL, and errors will be dealt with by the default
423  !! internal error handler. Possible errors that may be encountered are as
424  !! follows.
425  !! - LA_ARRAY_SIZE_ERROR: Occurs if any of the input arrays are not sized
426  !! appropriately.
427  subroutine tri_mtx_mult_c(upper, n, alpha, a, beta, b, err) &
428  bind(c, name = "tri_mtx_mult")
429  ! Arguments
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
436 
437  ! Local Variables
438  type(errors), pointer :: eptr
439 
440  ! Process
441  call get_errorhandler(err, eptr)
442  if (associated(eptr)) then
443  call tri_mtx_mult(logical(upper), alpha, a, beta, b, eptr)
444  else
445  call tri_mtx_mult(logical(upper), alpha, a, beta, b)
446  end if
447  end subroutine
448 
449 ! ******************************************************************************
450 ! LINALG_FACTOR ROUTINES
451 ! ------------------------------------------------------------------------------
470  subroutine lu_factor_c(m, n, a, ipvt, err) bind(C, name = "lu_factor")
471  ! Arguments
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
476 
477  ! Local Variables
478  type(errors), pointer :: eptr
479 
480  ! Process
481  call get_errorhandler(err, eptr)
482  if (associated(eptr)) then
483  call lu_factor(a, ipvt, eptr)
484  else
485  call lu_factor(a, ipvt)
486  end if
487  end subroutine
488 
489 ! ------------------------------------------------------------------------------
516  subroutine form_lu_c(n, lu, ipvt, u, p) bind(C, name = "form_lu")
517  ! Arguments
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)
522 
523  ! Process
524  call form_lu(lu, ipvt, u, p)
525  end subroutine
526 
527 ! ------------------------------------------------------------------------------
548  subroutine qr_factor_c(m, n, a, tau, err) bind(C, name = "qr_factor")
549  ! Arguments
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
554 
555  ! Local Variables
556  type(errors), pointer :: eptr
557 
558  ! Process
559  call get_errorhandler(err, eptr)
560  if (associated(eptr)) then
561  call qr_factor(a, tau, err = eptr)
562  else
563  call qr_factor(a, tau)
564  end if
565  end subroutine
566 
567 ! ------------------------------------------------------------------------------
592  subroutine qr_factor_pivot_c(m, n, a, tau, jpvt, err) &
593  bind(c, name = "qr_factor_pivot")
594  ! Arguments
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
600 
601  ! Local Variables
602  type(errors), pointer :: eptr
603 
604  ! Process
605  call get_errorhandler(err, eptr)
606  if (associated(eptr)) then
607  call qr_factor(a, tau, jpvt, err = eptr)
608  else
609  call qr_factor(a, tau, jpvt)
610  end if
611  end subroutine
612 
613 ! ------------------------------------------------------------------------------
638  subroutine form_qr_c(m, n, r, tau, q, err) bind(C, name = "form_qr")
639  ! Arguments
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
645 
646  ! Local Variables
647  type(errors), pointer :: eptr
648 
649  ! Process
650  call get_errorhandler(err, eptr)
651  if (associated(eptr)) then
652  call form_qr(r, tau, q, err = eptr)
653  else
654  call form_qr(r, tau, q)
655  end if
656  end subroutine
657 
658 ! ------------------------------------------------------------------------------
686  subroutine form_qr_pivot_c(m, n, r, tau, pvt, q, p, err) &
687  bind(c, name = "form_qr_pivot")
688  ! Arguments
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
695 
696  ! Local Variables
697  type(errors), pointer :: eptr
698 
699  ! Process
700  call get_errorhandler(err, eptr)
701  if (associated(eptr)) then
702  call form_qr(r, tau, pvt, q, p, err = eptr)
703  else
704  call form_qr(r, tau, pvt, q, p)
705  end if
706  end subroutine
707 
708 ! ------------------------------------------------------------------------------
731  subroutine mult_qr_c(trans, m, n, q, tau, c, err) &
732  bind(c, name = "mult_qr")
733  ! Arguments
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
739 
740  ! Local Variables
741  type(errors), pointer :: eptr
742 
743  ! Process
744  call get_errorhandler(err, eptr)
745  if (associated(eptr)) then
746  call mult_qr(.true., logical(trans), q, tau, c, err = eptr)
747  else
748  call mult_qr(.true., logical(trans), q, tau, c)
749  end if
750  end subroutine
751 
752 ! ------------------------------------------------------------------------------
772  subroutine qr_rank1_update_c(m, n, q, r, u, v, err) &
773  bind(c, name = "qr_rank1_update")
774  ! Arguments
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
778 
779  ! Local Variables
780  type(errors), pointer :: eptr
781 
782  ! Process
783  call get_errorhandler(err, eptr)
784  if (associated(eptr)) then
785  call qr_rank1_update(q, r, u, v, err = eptr)
786  else
787  call qr_rank1_update(q, r, u, v)
788  endif
789  end subroutine
790 
791 ! ------------------------------------------------------------------------------
808  subroutine cholesky_factor_c(n, a, upper, err) &
809  bind(c, name = "cholesky_factor")
810  ! Arguments
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
815 
816  ! Local Variables
817  type(errors), pointer :: eptr
818 
819  ! Process
820  call get_errorhandler(err, eptr)
821  if (associated(eptr)) then
822  call cholesky_factor(a, logical(upper), eptr)
823  else
824  call cholesky_factor(a, logical(upper))
825  end if
826  end subroutine
827 
828 ! ------------------------------------------------------------------------------
843  subroutine cholesky_rank1_update_c(n, r, u, err) &
844  bind(c, name = "cholesky_rank1_update")
845  ! Arguments
846  integer(int32), intent(in), value :: n
847  real(real64), intent(inout) :: r(n,n), u(n)
848  type(errorhandler), intent(inout) :: err
849 
850  ! Local Variables
851  type(errors), pointer :: eptr
852 
853  ! Process
854  call get_errorhandler(err, eptr)
855  if (associated(eptr)) then
856  call cholesky_rank1_update(r, u, err = eptr)
857  else
858  call cholesky_rank1_update(r, u)
859  end if
860  end subroutine
861 
862 ! ------------------------------------------------------------------------------
882  subroutine cholesky_rank1_downdate_c(n, r, u, err) &
883  bind(c, name = "cholesky_rank1_downdate")
884  ! Arguments
885  integer(int32), intent(in), value :: n
886  real(real64), intent(inout) :: r(n,n), u(n)
887  type(errorhandler), intent(inout) :: err
888 
889  ! Local Variables
890  type(errors), pointer :: eptr
891 
892  ! Process
893  call get_errorhandler(err, eptr)
894  if (associated(eptr)) then
895  call cholesky_rank1_downdate(r, u, err = eptr)
896  else
897  call cholesky_rank1_downdate(r, u)
898  end if
899  end subroutine
900 
901 ! ------------------------------------------------------------------------------
922  subroutine rz_factor_c(m, n, a, tau, err) bind(C, name = "rz_factor")
923  ! Arguments
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
928 
929  ! Local Variables
930  type(errors), pointer :: eptr
931 
932  ! Process
933  call get_errorhandler(err, eptr)
934  if (associated(eptr)) then
935  call rz_factor(a, tau, err = eptr)
936  else
937  call rz_factor(a, tau)
938  end if
939  end subroutine
940 
941 ! ------------------------------------------------------------------------------
963  subroutine mult_rz_c(trans, m, n, l, a, tau, c, err) &
964  bind(c, name = "mult_rz")
965  ! Arguments
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
971 
972  ! Local Variables
973  type(errors), pointer :: eptr
974 
975  ! Process
976  call get_errorhandler(err, eptr)
977  if (associated(eptr)) then
978  call mult_rz(.true., logical(trans), l, a, tau, c, err = eptr)
979  else
980  call mult_rz(.true., logical(trans), l, a, tau, c)
981  end if
982  end subroutine
983 
984 ! ------------------------------------------------------------------------------
1011  subroutine svd_c(m, n, a, s, u, vt, err) bind(C, name = "svd")
1012  ! Arguments
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
1017 
1018  ! Local Variables
1019  type(errors), pointer :: eptr
1020 
1021  ! Process
1022  call get_errorhandler(err, eptr)
1023  if (associated(eptr)) then
1024  call svd(a, s, u, vt, err = eptr)
1025  else
1026  call svd(a, s, u, vt)
1027  end if
1028  end subroutine
1029 
1030 ! ******************************************************************************
1031 ! LINALG_SOLVE ROUTINES
1032 ! ------------------------------------------------------------------------------
1050  subroutine solve_tri_mtx_c(upper, trans, nounit, n, nrhs, alpha, a, b) &
1051  bind(c, name = "solve_triangular_system")
1052  ! Arguments
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)
1058 
1059  ! Process
1060  call solve_triangular_system(.true., logical(upper), logical(trans), &
1061  logical(nounit), alpha, a, b)
1062  end subroutine
1063 
1064 ! ------------------------------------------------------------------------------
1074  subroutine solve_lu_c(n, nrhs, a, ipvt, b) bind(C, name = "solve_lu")
1075  ! Arguments
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)
1080 
1081  ! Process
1082  call solve_lu(a, ipvt, b)
1083  end subroutine
1084 
1085 ! ------------------------------------------------------------------------------
1106  subroutine solve_qr_c(m, n, nrhs, a, tau, b, err) &
1107  bind(c, name = "solve_qr")
1108  ! Arguments
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
1113 
1114  ! Local Variables
1115  type(errors), pointer :: eptr
1116 
1117  ! Process
1118  call get_errorhandler(err, eptr)
1119  if (associated(eptr)) then
1120  call solve_qr(a, tau, b, err = eptr)
1121  else
1122  call solve_qr(a, tau, b)
1123  end if
1124  end subroutine
1125 
1126 ! ------------------------------------------------------------------------------
1151  subroutine solve_qr_pivot_c(m, n, nrhs, a, tau, jpvt, b, err) &
1152  bind(c, name = "solve_qr_pivot")
1153  ! Arguments
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
1159 
1160  ! Local Variables
1161  type(errors), pointer :: eptr
1162 
1163  ! Process
1164  call get_errorhandler(err, eptr)
1165  if (associated(eptr)) then
1166  call solve_qr(a, tau, jpvt, b, err = eptr)
1167  else
1168  call solve_qr(a, tau, jpvt, b)
1169  end if
1170  end subroutine
1171 
1172 ! ------------------------------------------------------------------------------
1184  subroutine solve_cholesky_c(upper, n, nrhs, a, b) &
1185  bind(c, name = "solve_cholesky")
1186  ! Arguments
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)
1191 
1192  ! Process
1193  call solve_cholesky(logical(upper), a, b)
1194  end subroutine
1195 
1196 ! ------------------------------------------------------------------------------
1209  subroutine mtx_inverse_c(n, a, err) bind(C, name = "mtx_inverse")
1210  ! Arguments
1211  integer(int32), intent(in), value :: n
1212  real(real64), intent(inout) :: a(n,n)
1213  type(errorhandler), intent(inout) :: err
1214 
1215  ! Local Variables
1216  type(errors), pointer :: eptr
1217 
1218  ! Process
1219  call get_errorhandler(err, eptr)
1220  if (associated(eptr)) then
1221  call mtx_inverse(a, err = eptr)
1222  else
1223  call mtx_inverse(a)
1224  end if
1225  end subroutine
1226 
1227 ! ------------------------------------------------------------------------------
1245  subroutine mtx_pinverse_c(m, n, a, ainv, err) &
1246  bind(c, name = "mtx_pinverse")
1247  ! Arguments
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
1252 
1253  ! Local Variables
1254  type(errors), pointer :: eptr
1255 
1256  ! Process
1257  call get_errorhandler(err, eptr)
1258  if (associated(eptr)) then
1259  call mtx_pinverse(a, ainv, err = eptr)
1260  else
1261  call mtx_pinverse(a, ainv)
1262  end if
1263  end subroutine
1264 
1265 ! ------------------------------------------------------------------------------
1286  subroutine solve_least_squares_c(m, n, nrhs, a, b, err) &
1287  bind(c, name = "solve_least_squares")
1288  ! Arguments
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
1292 
1293  ! Local Variables
1294  type(errors), pointer :: eptr
1295 
1296  ! Process
1297  call get_errorhandler(err, eptr)
1298  if (associated(eptr)) then
1299  call solve_least_squares(a, b, err = eptr)
1300  else
1301  call solve_least_squares(a, b)
1302  end if
1303  end subroutine
1304 
1305 ! ******************************************************************************
1306 ! LINALG_EIGEN ROUTINES
1307 ! ------------------------------------------------------------------------------
1328  subroutine eigen_symm_c(n, vecs, a, vals, err) bind(C, name = "eigen_symm")
1329  ! Arguments
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
1335 
1336  ! Local Variables
1337  type(errors), pointer :: eptr
1338 
1339  ! Process
1340  call get_errorhandler(err, eptr)
1341  if (associated(eptr)) then
1342  call eigen(logical(vecs), a, vals, err = eptr)
1343  else
1344  call eigen(logical(vecs), a, vals)
1345  end if
1346  end subroutine
1347 
1348 ! ------------------------------------------------------------------------------
1366  subroutine eigen_asymm_c(n, a, vals, vecs, err) &
1367  bind(c, name = "eigen_asymm")
1368  ! Arguments
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
1373 
1374  ! Local Variables
1375  type(errors), pointer :: eptr
1376 
1377  ! Process
1378  call get_errorhandler(err, eptr)
1379  if (associated(eptr)) then
1380  call eigen(a, vals, vecs, err = eptr)
1381  else
1382  call eigen(a, vals, vecs)
1383  end if
1384  end subroutine
1385 
1386 ! ------------------------------------------------------------------------------
1415  subroutine eigen_gen_c(n, a, b, alpha, beta, vecs, err) &
1416  bind(c, name = "eigen_gen")
1417  ! Arguments
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
1423 
1424  ! Local Variables
1425  type(errors), pointer :: eptr
1426 
1427  ! Process
1428  call get_errorhandler(err, eptr)
1429  if (associated(eptr)) then
1430  call eigen(a, b, alpha, beta = beta, vecs = vecs, err = eptr)
1431  else
1432  call eigen(a, b, alpha, beta = beta, vecs = vecs)
1433  end if
1434  end subroutine
1435 
1436 ! ******************************************************************************
1437 ! SORTING ROUTINES
1438 ! ------------------------------------------------------------------------------
1454  subroutine sort_dbl_ind_c(ascend, n, x, ind) bind(C, name = "sort_dbl")
1455  ! Arguments
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
1460 
1461  ! Local Variables
1462  integer(int32), pointer, dimension(:) :: ptr
1463 
1464  ! Process
1465  if (c_associated(ind)) then
1466  call c_f_pointer(ind, ptr, [n])
1467  call sort(x, ptr, logical(ascend))
1468  else
1469  call sort(x, logical(ascend))
1470  end if
1471  end subroutine
1472 
1473 ! ------------------------------------------------------------------------------
1490  subroutine sort_cmplx_ind_c(ascend, n, x, ind) bind(C, name = "sort_cmplx")
1491  ! Arguments
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
1496 
1497  ! Local Variables
1498  integer(int32), pointer, dimension(:) :: ptr
1499 
1500  ! Process
1501  if (c_associated(ind)) then
1502  call c_f_pointer(ind, ptr, [n])
1503  call sort(x, ptr, logical(ascend))
1504  else
1505  call sort(x, logical(ascend))
1506  end if
1507  end subroutine
1508 
1509 ! ------------------------------------------------------------------------------
1521  subroutine sort_eigen_cmplx_c(ascend, n, vals, vecs) &
1522  bind(c, name = "sort_eigen_cmplx")
1523  ! Arguments
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)
1527 
1528  ! Process
1529  call sort(vals, vecs, logical(ascend))
1530  end subroutine
1531 
1532 ! ------------------------------------------------------------------------------
1544  subroutine sort_eigen_dbl_c(ascend, n, vals, vecs) &
1545  bind(c, name = "sort_eigen_dbl")
1546  ! Arguments
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)
1550 
1551  ! Process
1552  call sort(vals, vecs, logical(ascend))
1553  end subroutine
1554 
1555 ! ------------------------------------------------------------------------------
1556 end module
Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR fact...
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...
linalg_core
Definition: linalg_core.f90:18
Extracts the L and U matrices from the condensed [L\U] storage format used by the lu_factor...
Solves the overdetermined or underdetermined system (A*X = B) of M equations of N unknowns...
linalg_c_binding
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...
Definition: linalg_core.f90:72
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.
Sorts an array.
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.