!-*-f90-*- ! ! Interfaces: Linear Algebra support ! ! LU ! function gsl_linalg_lu_decomp (a, p, signum) & bind(c, name='gsl_linalg_LU_decomp') import :: c_ptr, c_int type(c_ptr), value :: a, p integer(c_int) :: signum integer(c_int) :: gsl_linalg_lu_decomp end function gsl_linalg_lu_decomp function gsl_linalg_complex_lu_decomp (a, p, signum) & bind(c, name='gsl_linalg_complex_LU_decomp') import :: c_ptr, c_int type(c_ptr), value :: a, p integer(c_int) :: signum integer(c_int) :: gsl_linalg_complex_lu_decomp end function gsl_linalg_complex_lu_decomp function gsl_linalg_lu_solve (lu, p, b, x) & bind(c, name='gsl_linalg_LU_solve') import :: c_ptr, c_int type(c_ptr), value :: lu, p, b, x integer(c_int) :: gsl_linalg_lu_solve end function gsl_linalg_lu_solve function gsl_linalg_complex_lu_solve (lu, p, b, x) & bind(c, name='gsl_linalg_complex_LU_solve') import :: c_ptr, c_int type(c_ptr), value :: lu, p, b, x integer(c_int) :: gsl_linalg_complex_lu_solve end function gsl_linalg_complex_lu_solve function gsl_linalg_lu_svx (lu, p, x) & bind(c, name='gsl_linalg_LU_svx') import :: c_ptr, c_int type(c_ptr), value :: lu, p, x integer(c_int) :: gsl_linalg_lu_svx end function gsl_linalg_lu_svx function gsl_linalg_complex_lu_svx (lu, p, x) & bind(c, name='gsl_linalg_complex_LU_svx') import :: c_ptr, c_int type(c_ptr), value :: lu, p, x integer(c_int) :: gsl_linalg_complex_lu_svx end function gsl_linalg_complex_lu_svx function gsl_linalg_lu_refine (a, lu, p, b, x, residual) & bind(c, name='gsl_linalg_LU_refine') import :: c_ptr, c_int type(c_ptr), value :: a, lu, p, b, x, residual integer(c_int) :: gsl_linalg_lu_refine end function gsl_linalg_lu_refine function gsl_linalg_complex_lu_refine (a, lu, p, b, x, residual) & bind(c, name='gsl_linalg_complex_LU_refine') import :: c_ptr, c_int type(c_ptr), value :: a, lu, p, b, x, residual integer(c_int) :: gsl_linalg_complex_lu_refine end function gsl_linalg_complex_lu_refine function gsl_linalg_lu_invert (lu, p, inv) & bind(c, name='gsl_linalg_LU_invert') import :: c_ptr, c_int type(c_ptr), value :: lu, p, inv integer(c_int) :: gsl_linalg_lu_invert end function gsl_linalg_lu_invert function gsl_linalg_complex_lu_invert (lu, p, inv) & bind(c, name='gsl_linalg_complex_LU_invert') import :: c_ptr, c_int type(c_ptr), value :: lu, p, inv integer(c_int) :: gsl_linalg_complex_lu_invert end function gsl_linalg_complex_lu_invert function gsl_linalg_lu_invx (lu, p) & bind(c, name='gsl_linalg_LU_invx') import :: c_ptr, c_int type(c_ptr), value :: lu, p integer(c_int) :: gsl_linalg_lu_invx end function gsl_linalg_lu_invx function gsl_linalg_complex_lu_invx (lu, p) & bind(c, name='gsl_linalg_complex_LU_invx') import :: c_ptr, c_int type(c_ptr), value :: lu, p integer(c_int) :: gsl_linalg_complex_lu_invx end function gsl_linalg_complex_lu_invx function gsl_linalg_lu_det(lu, signum) bind(c, name='gsl_linalg_LU_det') import :: c_ptr, c_double, c_int type(c_ptr), value :: lu integer(c_int), value :: signum real(c_double) :: gsl_linalg_lu_det end function gsl_linalg_lu_det function gsl_linalg_complex_lu_det(lu, signum) & bind(c, name='gsl_linalg_complex_LU_det') import :: c_ptr, c_int, gsl_complex type(c_ptr), value :: lu integer(c_int), value :: signum type(gsl_complex) :: gsl_linalg_complex_lu_det end function gsl_linalg_complex_lu_det function gsl_linalg_lu_lndet(lu) bind(c, name='gsl_linalg_LU_lndet') import :: c_ptr, c_double type(c_ptr), value :: lu real(c_double) :: gsl_linalg_lu_lndet end function gsl_linalg_lu_lndet function gsl_linalg_complex_lu_lndet(lu) & bind(c, name='gsl_linalg_complex_LU_lndet') import :: c_ptr, c_double type(c_ptr), value :: lu real(c_double) :: gsl_linalg_complex_lu_lndet end function gsl_linalg_complex_lu_lndet function gsl_linalg_lu_sgndet(lu, signum) bind(c, name='gsl_linalg_LU_sgndet') import :: c_ptr, c_double, c_int type(c_ptr), value :: lu integer(c_int), value :: signum integer(c_int) :: gsl_linalg_lu_sgndet end function gsl_linalg_lu_sgndet function gsl_linalg_complex_lu_sgndet(lu, signum) & bind(c, name='gsl_linalg_complex_LU_sgndet') import :: c_ptr, c_int, gsl_complex type(c_ptr), value :: lu integer(c_int), value :: signum type(gsl_complex) :: gsl_linalg_complex_lu_sgndet end function gsl_linalg_complex_lu_sgndet ! ! QR ! function gsl_linalg_qr_decomp (a, tau) bind(c, name='gsl_linalg_QR_decomp') import :: c_ptr, c_int type(c_ptr), value :: a, tau integer(c_int) :: gsl_linalg_qr_decomp end function gsl_linalg_qr_decomp function gsl_linalg_qr_decomp_r (a, t) bind(c, name='gsl_linalg_QR_decomp_r') import :: c_ptr, c_int type(c_ptr), value :: a, t integer(c_int) :: gsl_linalg_qr_decomp_r end function gsl_linalg_qr_decomp_r function gsl_linalg_qr_solve (qr, tau, b, x) & bind(c, name='gsl_linalg_QR_solve') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, b, x integer(c_int) :: gsl_linalg_qr_solve end function gsl_linalg_qr_solve function gsl_linalg_qr_solve_r (qr, t, b, x) & bind(c, name='gsl_linalg_QR_solve_r') import :: c_ptr, c_int type(c_ptr), value :: qr, t, b, x integer(c_int) :: gsl_linalg_qr_solve_r end function gsl_linalg_qr_solve_r function gsl_linalg_qr_svx (qr, tau, x) & bind(c, name='gsl_linalg_QR_svx') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, x integer(c_int) :: gsl_linalg_qr_svx end function gsl_linalg_qr_svx function gsl_linalg_qr_lssolve (qr, tau, b, x, residual) & bind(c, name='gsl_linalg_QR_lssolve') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, b, x, residual integer(c_int) :: gsl_linalg_qr_lssolve end function gsl_linalg_qr_lssolve function gsl_linalg_qr_lssolve_r (qr, t, b, x, work) & bind(c, name='gsl_linalg_QR_lssolve_r') import :: c_ptr, c_int type(c_ptr), value :: qr, t, b, x, work integer(c_int) :: gsl_linalg_qr_lssolve_r end function gsl_linalg_qr_lssolve_r function gsl_linalg_qr_qtvec (qr, tau, v) & bind(c, name='gsl_linalg_QR_QTvec') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, v integer(c_int) :: gsl_linalg_qr_qtvec end function gsl_linalg_qr_qtvec function gsl_linalg_qr_qtvec_r (qr, t, v, work) & bind(c, name='gsl_linalg_QR_QTvec_r') import :: c_ptr, c_int type(c_ptr), value :: qr, t, v,work integer(c_int) :: gsl_linalg_qr_qtvec_r end function gsl_linalg_qr_qtvec_r function gsl_linalg_qr_qvec (qr, tau, v) & bind(c, name='gsl_linalg_QR_Qvec') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, v integer(c_int) :: gsl_linalg_qr_qvec end function gsl_linalg_qr_qvec function gsl_linalg_qr_qtmat (qr, tau, a) & bind(c, name='gsl_linalg_QR_QTmat') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, a integer(c_int) :: gsl_linalg_qr_qtmat end function gsl_linalg_qr_qtmat function gsl_linalg_qr_qtmat_r (qr, t, a, work) & bind(c, name='gsl_linalg_QR_QTmat_r') import :: c_ptr, c_int type(c_ptr), value :: qr, t, a, work integer(c_int) :: gsl_linalg_qr_qtmat_r end function gsl_linalg_qr_qtmat_r function gsl_linalg_qr_rsolve (qr, b, x) & bind(c, name='gsl_linalg_QR_Rsolve') import :: c_ptr, c_int type(c_ptr), value :: qr, b, x integer(c_int) :: gsl_linalg_qr_rsolve end function gsl_linalg_qr_rsolve function gsl_linalg_qr_rsvx (qr, x) & bind(c, name='gsl_linalg_QR_Rsvx') import :: c_ptr, c_int type(c_ptr), value :: qr, x integer(c_int) :: gsl_linalg_qr_rsvx end function gsl_linalg_qr_rsvx function gsl_linalg_qr_unpack (qr, tau, q, r) & bind(c, name='gsl_linalg_QR_unpack') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, q, r integer(c_int) :: gsl_linalg_qr_unpack end function gsl_linalg_qr_unpack function gsl_linalg_qr_unpack_r (qr, t, q, r) & bind(c, name='gsl_linalg_QR_unpack_r') import :: c_ptr, c_int type(c_ptr), value :: qr, t, q, r integer(c_int) :: gsl_linalg_qr_unpack_r end function gsl_linalg_qr_unpack_r function gsl_linalg_qr_qrsolve (q, r, b, x) & bind(c, name='gsl_linalg_QR_QRsolve') import :: c_ptr, c_int type(c_ptr), value :: q, r, b, x integer(c_int) :: gsl_linalg_qr_qrsolve end function gsl_linalg_qr_qrsolve function gsl_linalg_qr_update (q, r, w, v) & bind(c, name='gsl_linalg_QR_update') import :: c_ptr, c_int type(c_ptr), value :: q, r, w, v integer(c_int) :: gsl_linalg_qr_update end function gsl_linalg_qr_update function gsl_linalg_r_solve (r, b, x) & bind(c, name='gsl_linalg_R_solve') import :: c_ptr, c_int type(c_ptr), value :: r, b, x integer(c_int) :: gsl_linalg_r_solve end function gsl_linalg_r_solve function gsl_linalg_r_svx (r, x) & bind(c, name='gsl_linalg_R_svx') import :: c_ptr, c_int type(c_ptr), value :: r, x integer(c_int) :: gsl_linalg_r_svx end function gsl_linalg_r_svx function gsl_linalg_qrpt_decomp (a, tau, p, signum, norm) & bind(c, name='gsl_linalg_QRPT_decomp') import :: c_ptr, c_int type(c_ptr), value :: a, tau, p, norm integer(c_int), intent(out) :: signum integer(c_int) :: gsl_linalg_qrpt_decomp end function gsl_linalg_qrpt_decomp function gsl_linalg_qrpt_decomp2 (a, q, r, tau, p, signum, norm) & bind(c, name='gsl_linalg_QRPT_decomp2') import :: c_ptr, c_int type(c_ptr), value :: a, q, r, tau, p, norm integer(c_int), intent(out) :: signum integer(c_int) :: gsl_linalg_qrpt_decomp2 end function gsl_linalg_qrpt_decomp2 function gsl_linalg_qrpt_solve (qr, tau, p, b, x) & bind(c, name='gsl_linalg_QRPT_solve') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, p, b, x integer(c_int) :: gsl_linalg_qrpt_solve end function gsl_linalg_qrpt_solve function gsl_linalg_qrpt_svx (qr, tau, p, x) & bind(c, name='gsl_linalg_QRPT_svx') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, p, x integer(c_int) :: gsl_linalg_qrpt_svx end function gsl_linalg_qrpt_svx function gsl_linalg_qrpt_lssolve (qr, tau, p, b, x, r) & bind(c, name='gsl_linalg_QRPT_lssolve') import :: c_ptr, c_int type(c_ptr), value :: qr, tau, p, b, x, r integer(c_int) :: gsl_linalg_qrpt_lssolve end function gsl_linalg_qrpt_lssolve function gsl_linalg_qrpt_lssolve2 (qr, tau, p, b, rank, x, r) & bind(c, name='gsl_linalg_QRPT_lssolve2') import :: c_ptr, c_int, c_size_t type(c_ptr), value :: qr, tau, p, b, x, r integer(c_size_t), value :: rank integer(c_int) :: gsl_linalg_qrpt_lssolve2 end function gsl_linalg_qrpt_lssolve2 function gsl_linalg_qrpt_qrsolve (q, r, p, b, x) & bind(c, name='gsl_linalg_QRPT_QRsolve') import :: c_ptr, c_int type(c_ptr), value :: q, r, p, b, x integer(c_int) :: gsl_linalg_qrpt_qrsolve end function gsl_linalg_qrpt_qrsolve function gsl_linalg_qrpt_update (q, r, p, w, v) & bind(c, name='gsl_linalg_QRPT_update') import :: c_ptr, c_int type(c_ptr), value :: q, r, p, w, v integer(c_int) :: gsl_linalg_qrpt_update end function gsl_linalg_qrpt_update function gsl_linalg_qrpt_rsolve (qr, p, b, x) & bind(c, name='gsl_linalg_QRPT_Rsolve') import :: c_ptr, c_int type(c_ptr), value :: qr, p, b, x integer(c_int) :: gsl_linalg_qrpt_rsolve end function gsl_linalg_qrpt_rsolve function gsl_linalg_qrpt_rsvx (qr, p, x) & bind(c, name='gsl_linalg_QRPT_Rsvx') import :: c_ptr, c_int type(c_ptr), value :: qr, p, x integer(c_int) :: gsl_linalg_qrpt_rsvx end function gsl_linalg_qrpt_rsvx function gsl_linalg_qrpt_rank (qr, tol) & bind(c, name='gsl_linalg_QRPT_rank') import :: c_ptr, c_size_t, c_double type(c_ptr), value :: qr real(c_double), value :: tol integer(c_size_t) :: gsl_linalg_qrpt_rank end function gsl_linalg_qrpt_rank function gsl_linalg_qrpt_rcond (qr, rcond, wk) & bind(c, name='gsl_linalg_QRPT_rcond') import :: c_ptr, c_int, c_double type(c_ptr), value :: qr, wk real(c_double) :: rcond integer(c_int) :: gsl_linalg_qrpt_rcond end function gsl_linalg_qrpt_rcond ! ! LQ decomposition ! integer(c_int) function gsl_linalg_lq_decomp(a, tau) & bind(c, name='gsl_linalg_LQ_decomp') import :: c_int, c_ptr type(c_ptr), value :: a, tau end function gsl_linalg_lq_decomp integer(c_int) function gsl_linalg_lq_lssolve(lq, tau, b, x, residual) & bind(c, name='gsl_linalg_LQ_lssolve') import :: c_int, c_ptr type(c_ptr), value :: lq, tau, b, x, residual end function gsl_linalg_lq_lssolve integer(c_int) function gsl_linalg_lq_unpack(lq, tau, q, l) & bind(c, name='gsl_linalg_LQ_unpack') import :: c_int, c_ptr type(c_ptr), value :: lq, tau, q, l end function gsl_linalg_lq_unpack integer(c_int) function gsl_linalg_lq_qtvec(lq, tau, v) & bind(c, name='gsl_linalg_LQ_QTvec') import :: c_int, c_ptr type(c_ptr), value :: lq, tau, v end function gsl_linalg_lq_qtvec ! ! Complete orthogonal decomposition (COD) ! function gsl_linalg_cod_decomp (a, tau_q, tau_z, p, rank, work) & bind(c, name='gsl_linalg_COD_decomp') import :: c_ptr, c_int, c_size_t type(c_ptr), value :: a, tau_q, tau_z, p, work integer(c_size_t) :: rank integer(c_int) :: gsl_linalg_cod_decomp end function gsl_linalg_cod_decomp function gsl_linalg_cod_decomp_e (a, tau_q, tau_z, p, tol, rank, work) & bind(c, name='gsl_linalg_COD_decomp_e') import :: c_ptr, c_int, c_size_t, c_double type(c_ptr), value :: a, tau_q, tau_z, p, work real(c_double), value :: tol integer(c_size_t) :: rank integer(c_int) :: gsl_linalg_cod_decomp_e end function gsl_linalg_cod_decomp_e function gsl_linalg_cod_lssolve (qrzt, tau_q, tau_z, p, rank, b, x, residual) & bind(c, name='gsl_linalg_COD_lssolve') import :: c_ptr, c_int, c_size_t type(c_ptr), value :: qrzt, tau_q, tau_z, p, b, x, residual integer(c_size_t), value :: rank integer(c_int) :: gsl_linalg_cod_lssolve end function gsl_linalg_cod_lssolve function gsl_linalg_cod_lssolve2 (lambda, qrzt, tau_q, tau_z, p, rank, b, & x, residual, s, work) bind(c, name='gsl_linalg_COD_lssolve2') import :: c_ptr, c_int, c_size_t, c_double real(c_double), value :: lambda type(c_ptr), value :: qrzt, tau_q, tau_z, p, b, x, residual, s, work integer(c_size_t), value :: rank integer(c_int) :: gsl_linalg_cod_lssolve2 end function gsl_linalg_cod_lssolve2 function gsl_linalg_cod_unpack (qrzt, tau_q, tau_z, p, rank, q, r, z) & bind(c, name='gsl_linalg_COD_unpack') import :: c_ptr, c_int, c_size_t type(c_ptr), value :: qrzt, tau_q, tau_z, p, q, r, z integer(c_size_t), value :: rank integer(c_int) :: gsl_linalg_cod_unpack end function gsl_linalg_cod_unpack function gsl_linalg_cod_matz (qrzt, tau_z, rank, a, work) & bind(c, name='gsl_linalg_COD_matZ') import :: c_ptr, c_int, c_size_t type(c_ptr), value :: qrzt, tau_z, a, work integer(c_size_t), value :: rank integer(c_int) :: gsl_linalg_cod_matz end function gsl_linalg_cod_matz ! ! SVD ! function gsl_linalg_sv_decomp (a, v, s, work) & bind(c, name='gsl_linalg_SV_decomp') import :: c_ptr, c_int type(c_ptr), value :: a, v, s, work integer(c_int) :: signum integer(c_int) :: gsl_linalg_sv_decomp end function gsl_linalg_sv_decomp function gsl_linalg_sv_decomp_mod (a, x, v, s, work) & bind(c, name='gsl_linalg_SV_decomp_mod') import :: c_ptr, c_int type(c_ptr), value :: a, x, v, s, work integer(c_int) :: signum integer(c_int) :: gsl_linalg_sv_decomp_mod end function gsl_linalg_sv_decomp_mod function gsl_linalg_sv_decomp_jacobi (a, v, s) & bind(c, name='gsl_linalg_SV_decomp_jacobi') import :: c_ptr, c_int type(c_ptr), value :: a, v, s integer(c_int) :: signum integer(c_int) :: gsl_linalg_sv_decomp_jacobi end function gsl_linalg_sv_decomp_jacobi function gsl_linalg_sv_solve (u, v, s, b, x) & bind(c, name='gsl_linalg_SV_solve') import :: c_ptr, c_int type(c_ptr), value :: u, v, s, b, x integer(c_int) :: signum integer(c_int) :: gsl_linalg_sv_solve end function gsl_linalg_sv_solve function gsl_linalg_sv_leverage(u, h) & bind(c, name='gsl_linalg_SV_leverage') import :: c_ptr, c_int type(c_ptr), value :: u, h integer(c_int) :: gsl_linalg_sv_leverage end function gsl_linalg_sv_leverage ! ! Cholesky ! function gsl_linalg_cholesky_decomp1 (a) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a integer(c_int) :: gsl_linalg_cholesky_decomp1 end function gsl_linalg_cholesky_decomp1 function gsl_linalg_cholesky_decomp (a) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a integer(c_int) :: gsl_linalg_cholesky_decomp end function gsl_linalg_cholesky_decomp function gsl_linalg_complex_cholesky_decomp (a) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a integer(c_int) :: gsl_linalg_complex_cholesky_decomp end function gsl_linalg_complex_cholesky_decomp function gsl_linalg_cholesky_solve (chol, b, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: chol, b, x integer(c_int) :: gsl_linalg_cholesky_solve end function gsl_linalg_cholesky_solve function gsl_linalg_complex_cholesky_solve (chol, b, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: chol, b, x integer(c_int) :: gsl_linalg_complex_cholesky_solve end function gsl_linalg_complex_cholesky_solve function gsl_linalg_cholesky_svx (chol, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: chol, x integer(c_int) :: gsl_linalg_cholesky_svx end function gsl_linalg_cholesky_svx function gsl_linalg_complex_cholesky_svx (chol, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: chol, x integer(c_int) :: gsl_linalg_complex_cholesky_svx end function gsl_linalg_complex_cholesky_svx function gsl_linalg_cholesky_invert (chol) bind(c) import :: c_ptr, c_int integer(c_int) :: gsl_linalg_cholesky_invert type(c_ptr), value :: chol end function gsl_linalg_cholesky_invert function gsl_linalg_complex_cholesky_invert (chol) bind(c) import :: c_ptr, c_int integer(c_int) :: gsl_linalg_complex_cholesky_invert type(c_ptr), value :: chol end function gsl_linalg_complex_cholesky_invert function gsl_linalg_cholesky_decomp2 (a,s) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a,s integer(c_int) :: gsl_linalg_cholesky_decomp2 end function gsl_linalg_cholesky_decomp2 function gsl_linalg_cholesky_solve2 (chol, s, b, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: chol, s, b, x integer(c_int) :: gsl_linalg_cholesky_solve2 end function gsl_linalg_cholesky_solve2 function gsl_linalg_cholesky_svx2 (chol, s, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: chol, s, x integer(c_int) :: gsl_linalg_cholesky_svx2 end function gsl_linalg_cholesky_svx2 function gsl_linalg_cholesky_scale (a,s) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a,s integer(c_int) :: gsl_linalg_cholesky_scale end function gsl_linalg_cholesky_scale function gsl_linalg_cholesky_scale_apply (a,s) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a,s integer(c_int) :: gsl_linalg_cholesky_scale_apply end function gsl_linalg_cholesky_scale_apply function gsl_linalg_cholesky_rcond (chol, rcond, work) & bind(c) import :: c_ptr, c_int, c_double type(c_ptr), value :: chol, work real(c_double), intent(inout) :: rcond integer(c_int) :: gsl_linalg_cholesky_rcond end function gsl_linalg_cholesky_rcond ! ! pivoted Cholesky ! function gsl_linalg_pcholesky_decomp (a, p) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, p integer(c_int) :: gsl_linalg_pcholesky_decomp end function gsl_linalg_pcholesky_decomp function gsl_linalg_pcholesky_solve (ldlt, p, b, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: ldlt, p, b, x integer(c_int) :: gsl_linalg_pcholesky_solve end function gsl_linalg_pcholesky_solve function gsl_linalg_pcholesky_svx (ldlt, p, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: ldlt, p, x integer(c_int) :: gsl_linalg_pcholesky_svx end function gsl_linalg_pcholesky_svx function gsl_linalg_pcholesky_invert (ldlt, p, ainv) bind(c) import :: c_ptr, c_int integer(c_int) :: gsl_linalg_pcholesky_invert type(c_ptr), value :: ldlt, p, ainv end function gsl_linalg_pcholesky_invert function gsl_linalg_pcholesky_decomp2 (a, p,s) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, p,s integer(c_int) :: gsl_linalg_pcholesky_decomp2 end function gsl_linalg_pcholesky_decomp2 function gsl_linalg_pcholesky_solve2 (ldlt, p, s, b, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: ldlt, p, s, b, x integer(c_int) :: gsl_linalg_pcholesky_solve2 end function gsl_linalg_pcholesky_solve2 function gsl_linalg_pcholesky_svx2 (ldlt, p, s, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: ldlt, p, s, x integer(c_int) :: gsl_linalg_pcholesky_svx2 end function gsl_linalg_pcholesky_svx2 function gsl_linalg_pcholesky_rcond (ldlt, p, rcond, work) & bind(c) import :: c_ptr, c_int, c_double type(c_ptr), value :: ldlt, p, work real(c_double), intent(inout) :: rcond integer(c_int) :: gsl_linalg_pcholesky_rcond end function gsl_linalg_pcholesky_rcond ! ! modified Cholesky ! function gsl_linalg_mcholesky_decomp (a, p, e) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, p, e integer(c_int) :: gsl_linalg_mcholesky_decomp end function gsl_linalg_mcholesky_decomp function gsl_linalg_mcholesky_solve (ldlt, p, b, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: ldlt, p, b, x integer(c_int) :: gsl_linalg_mcholesky_solve end function gsl_linalg_mcholesky_solve function gsl_linalg_mcholesky_svx (ldlt, p, x) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: ldlt, p, x integer(c_int) :: gsl_linalg_mcholesky_svx end function gsl_linalg_mcholesky_svx function gsl_linalg_mcholesky_invert (ldlt, p, ainv) bind(c) import :: c_ptr, c_int integer(c_int) :: gsl_linalg_mcholesky_invert type(c_ptr), value :: ldlt, p, ainv end function gsl_linalg_mcholesky_invert function gsl_linalg_mcholesky_rcond (ldlt, p, rcond, work) & bind(c) import :: c_ptr, c_int, c_double type(c_ptr), value :: ldlt, p, work real(c_double), intent(inout) :: rcond integer(c_int) :: gsl_linalg_mcholesky_rcond end function gsl_linalg_mcholesky_rcond ! ! LDLT ! function gsl_linalg_ldlt_decomp (a) bind(c) import :: c_ptr, c_int type(c_ptr), value :: a integer(c_int) :: gsl_linalg_ldlt_decomp end function gsl_linalg_ldlt_decomp integer(c_int) function gsl_linalg_ldlt_solve(ldlt, b, x) bind(c) import :: c_int, c_ptr type(c_ptr), value :: ldlt, b, x end function gsl_linalg_ldlt_solve integer(c_int) function gsl_linalg_ldlt_svx(ldlt, x) bind(c) import :: c_int, c_ptr type(c_ptr), value :: ldlt, x end function gsl_linalg_ldlt_svx integer(c_int) function gsl_linalg_ldlt_invert(ldlt, ainv) bind(c) import :: c_int, c_ptr type(c_ptr), value :: ldlt, ainv end function gsl_linalg_ldlt_invert integer(c_int) function gsl_linalg_ldlt_rcond(ldlt, rcond, w) bind(c) import :: c_int, c_ptr, c_double type(c_ptr), value :: ldlt, w real(c_double) :: rcond end function gsl_linalg_ldlt_rcond ! ! Tridiagonal Decomposition of Real Symmetric Matrices ! function gsl_linalg_symmtd_decomp (a, tau) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau integer(c_int) :: gsl_linalg_symmtd_decomp end function gsl_linalg_symmtd_decomp function gsl_linalg_symmtd_unpack (a, tau, q, diag, subdiag) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau, q, diag, subdiag integer(c_int) :: gsl_linalg_symmtd_unpack end function gsl_linalg_symmtd_unpack function gsl_linalg_symmtd_unpack_t (a, diag, subdiag) & bind(c, name='gsl_linalg_symmtd_unpack_T') import :: c_ptr, c_int type(c_ptr), value :: a, diag, subdiag integer(c_int) :: gsl_linalg_symmtd_unpack_t end function gsl_linalg_symmtd_unpack_t ! ! Tridiagonal Decomposition of Hermitian Matrices ! function gsl_linalg_hermtd_decomp (a, tau) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau integer(c_int) :: gsl_linalg_hermtd_decomp end function gsl_linalg_hermtd_decomp function gsl_linalg_hermtd_unpack (a, tau, q, diag, subdiag) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau, q, diag, subdiag integer(c_int) :: gsl_linalg_hermtd_unpack end function gsl_linalg_hermtd_unpack function gsl_linalg_hermtd_unpack_t (a, diag, subdiag) & bind(c, name='gsl_linalg_hermtd_unpack_T') import :: c_ptr, c_int type(c_ptr), value :: a, diag, subdiag integer(c_int) :: gsl_linalg_hermtd_unpack_t end function gsl_linalg_hermtd_unpack_t ! ! Hessenberg ! function gsl_linalg_hessenberg_decomp (a, tau) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau integer(c_int) :: gsl_linalg_hessenberg_decomp end function gsl_linalg_hessenberg_decomp function gsl_linalg_hessenberg_unpack (h, tau, u) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: h, tau, u integer(c_int) :: gsl_linalg_hessenberg_unpack end function gsl_linalg_hessenberg_unpack function gsl_linalg_hessenberg_unpack_accum (h, tau, v) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: h, tau, v integer(c_int) :: gsl_linalg_hessenberg_unpack_accum end function gsl_linalg_hessenberg_unpack_accum function gsl_linalg_hessenberg_set_zero (h) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: h integer(c_int) :: gsl_linalg_hessenberg_set_zero end function gsl_linalg_hessenberg_set_zero function gsl_linalg_hesstri_decomp (a, b, u, v, work) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, b, u, v, work integer(c_int) :: gsl_linalg_hesstri_decomp end function gsl_linalg_hesstri_decomp ! ! Bidiag ! function gsl_linalg_bidiag_decomp (a, tau_u, tau_v) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau_u, tau_v integer(c_int) :: gsl_linalg_bidiag_decomp end function gsl_linalg_bidiag_decomp function gsl_linalg_bidiag_unpack (a, tau_u, u, tau_v, v, diag, & superdiag) bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau_u, tau_v, u, v, diag, superdiag integer(c_int) :: gsl_linalg_bidiag_unpack end function gsl_linalg_bidiag_unpack function gsl_linalg_bidiag_unpack2 (a, tau_u, tau_v, v) & bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, tau_u, tau_v, v integer(c_int) :: gsl_linalg_bidiag_unpack2 end function gsl_linalg_bidiag_unpack2 function gsl_linalg_bidiag_unpack_b (a, diag, superdiag) & bind(c, name='gsl_linalg_bidiag_unpack_B') import :: c_ptr, c_int type(c_ptr), value :: a, diag, superdiag integer(c_int) :: gsl_linalg_bidiag_unpack_b end function gsl_linalg_bidiag_unpack_b ! ! Householder ! function gsl_linalg_householder_transform (v) bind(c) import :: c_ptr, c_double type(c_ptr), value :: v real(c_double) :: gsl_linalg_householder_transform end function gsl_linalg_householder_transform function gsl_linalg_complex_householder_transform (v) bind(c) import :: c_ptr, gsl_complex type(c_ptr), value :: v type(gsl_complex) :: gsl_linalg_complex_householder_transform end function gsl_linalg_complex_householder_transform function gsl_linalg_householder_hm (tau, v, a) bind(c) import :: c_ptr, c_double, c_int real(c_double), value :: tau type(c_ptr), value :: v, a integer(c_int) :: gsl_linalg_householder_hm end function gsl_linalg_householder_hm function gsl_linalg_complex_householder_hm (tau, v, a) bind(c) import :: c_ptr, c_int, gsl_complex type(gsl_complex), value :: tau type(c_ptr), value :: v, a integer(c_int) :: gsl_linalg_complex_householder_hm end function gsl_linalg_complex_householder_hm function gsl_linalg_householder_mh (tau, v, a) bind(c) import :: c_ptr, c_double, c_int real(c_double), value :: tau type(c_ptr), value :: v, a integer(c_int) :: gsl_linalg_householder_mh end function gsl_linalg_householder_mh function gsl_linalg_complex_householder_mh (tau, v, a) bind(c) import :: c_ptr, c_int, gsl_complex type(gsl_complex), value :: tau type(c_ptr), value :: v, a integer(c_int) :: gsl_linalg_complex_householder_mh end function gsl_linalg_complex_householder_mh function gsl_linalg_householder_hv (tau, v, w) bind(c) import :: c_ptr, c_double, c_int real(c_double), value :: tau type(c_ptr), value :: v, w integer(c_int) :: gsl_linalg_householder_hv end function gsl_linalg_householder_hv function gsl_linalg_complex_householder_hv (tau, v, w) bind(c) import :: c_ptr, c_int, gsl_complex type(gsl_complex), value :: tau type(c_ptr), value :: v, w integer(c_int) :: gsl_linalg_complex_householder_hv end function gsl_linalg_complex_householder_hv function gsl_linalg_hh_solve (a, b, x) & bind(c, name='gsl_linalg_HH_solve') import :: c_ptr, c_int type(c_ptr), value :: a, b, x integer(c_int) :: gsl_linalg_hh_solve end function gsl_linalg_hh_solve function gsl_linalg_hh_svx (a, x) & bind(c, name='gsl_linalg_HH_svx') import :: c_ptr, c_int type(c_ptr), value :: a, x integer(c_int) :: gsl_linalg_hh_svx end function gsl_linalg_hh_svx ! ! Tridiagonal ! function gsl_linalg_solve_tridiag(diag, e, f, b, x) bind(c) import :: c_ptr, c_int type(c_ptr), value :: diag, e, f, b, x integer(c_int) :: gsl_linalg_solve_tridiag end function gsl_linalg_solve_tridiag function gsl_linalg_solve_symm_tridiag(diag, e, b, x) bind(c) import :: c_ptr, c_int type(c_ptr), value :: diag, e, b, x integer(c_int) :: gsl_linalg_solve_symm_tridiag end function gsl_linalg_solve_symm_tridiag function gsl_linalg_solve_cyc_tridiag(diag, e, f, b, x) bind(c) import :: c_ptr, c_int type(c_ptr), value :: diag, e, f, b, x integer(c_int) :: gsl_linalg_solve_cyc_tridiag end function gsl_linalg_solve_cyc_tridiag function gsl_linalg_solve_symm_cyc_tridiag(diag, e, b, x) bind(c) import :: c_ptr, c_int type(c_ptr), value :: diag, e, b, x integer(c_int) :: gsl_linalg_solve_symm_cyc_tridiag end function gsl_linalg_solve_symm_cyc_tridiag function gsl_linalg_qr_matq(QR, tau, A) & bind(c, name='gsl_linalg_QR_matQ') import :: c_ptr, c_int type(c_ptr), value :: QR, tau, A integer(c_int) :: gsl_linalg_qr_matq end function gsl_linalg_qr_matq subroutine gsl_linalg_givens(a, b, c, s) bind(c) import :: c_double real(c_double), value :: a, b real(c_double):: c, s end subroutine gsl_linalg_givens subroutine gsl_linalg_givens_gv(v, i, j, c, s) bind(c) import :: c_ptr, c_size_t, c_double type(c_ptr), value :: v integer(c_size_t), value :: i, j real(c_double), value :: c, s end subroutine gsl_linalg_givens_gv ! ! Triangular matrices ! integer(c_int) function gsl_linalg_tri_invert(uplo, diag, t) bind(c) import :: c_int, c_ptr integer(c_int), value :: uplo, diag type(c_ptr), value :: t end function gsl_linalg_tri_invert integer(c_int) function gsl_linalg_complex_tri_invert(uplo, diag, t) bind(c) import :: c_int, c_ptr integer(c_int), value :: uplo, diag type(c_ptr), value :: t end function gsl_linalg_complex_tri_invert integer(c_int) function gsl_linalg_tri_ltl(l) bind(c, name='gsl_linalg_tri_LTL') import :: c_int, c_ptr type(c_ptr), value :: l end function gsl_linalg_tri_ltl integer(c_int) function gsl_linalg_complex_tri_lhl(l) bind(c, name='gsl_linalg_complex_tri_LHL') import :: c_int, c_ptr type(c_ptr), value :: l end function gsl_linalg_complex_tri_lhl integer(c_int) function gsl_linalg_tri_ul(lu) bind(c, name='gsl_linalg_tri_UL') import :: c_int, c_ptr type(c_ptr), value :: lu end function gsl_linalg_tri_ul integer(c_int) function gsl_linalg_complex_tri_ul(lu) bind(c, name='gsl_linalg_complex_tri_UL') import :: c_int, c_ptr type(c_ptr), value :: lu end function gsl_linalg_complex_tri_ul integer(c_int) function gsl_linalg_tri_rcond(uplo, a, rcond, work) bind(c) import :: c_int, c_ptr, c_double integer(c_int), value :: uplo real(c_double), intent(inout) :: rcond type(c_ptr), value :: a, work end function gsl_linalg_tri_rcond ! ! Triangular matrices (legacy) ! integer(c_int) function gsl_linalg_tri_upper_invert(t) bind(c) import :: c_int, c_ptr type(c_ptr), value :: t end function integer(c_int) function gsl_linalg_tri_lower_invert(t) bind(c) import :: c_int, c_ptr type(c_ptr), value :: t end function integer(c_int) function gsl_linalg_tri_upper_unit_invert(t) bind(c) import :: c_int, c_ptr type(c_ptr), value :: t end function integer(c_int) function gsl_linalg_tri_lower_unit_invert(t) bind(c) import :: c_int, c_ptr type(c_ptr), value :: t end function integer(c_int) function gsl_linalg_tri_upper_rcond(t, rcond, work) bind(c) import :: c_int, c_double, c_ptr type(c_ptr), value :: t, work real(c_double) :: rcond end function integer(c_int) function gsl_linalg_tri_lower_rcond(t, rcond, work) bind(c) import :: c_int, c_double, c_ptr type(c_ptr), value :: t, work real(c_double) :: rcond end function ! ! Banded systems ! integer(c_int) function gsl_linalg_cholesky_band_decomp(a) bind(c) import :: c_int, c_ptr type(c_ptr), value :: a end function gsl_linalg_cholesky_band_decomp integer(c_int) function gsl_linalg_cholesky_band_solve(llt, b, x) bind(c) import :: c_int, c_ptr type(c_ptr), value :: llt, b, x end function gsl_linalg_cholesky_band_solve integer(c_int) function gsl_linalg_cholesky_band_svx(llt, x) bind(c) import :: c_int, c_ptr type(c_ptr), value :: llt, x end function gsl_linalg_cholesky_band_svx integer(c_int) function gsl_linalg_cholesky_band_invert(llt, ainv) bind(c) import :: c_int, c_ptr type(c_ptr), value :: llt, ainv end function gsl_linalg_cholesky_band_invert integer(c_int) function gsl_linalg_cholesky_band_unpack(llt, l) bind(c) import :: c_int, c_ptr type(c_ptr), value :: llt, l end function gsl_linalg_cholesky_band_unpack integer(c_int) function gsl_linalg_cholesky_band_rcond(llt, rcond, w) bind(c) import :: c_int, c_ptr, c_double type(c_ptr), value :: llt, w real(c_double) :: rcond end function gsl_linalg_cholesky_band_rcond integer(c_int) function gsl_linalg_ldlt_band_decomp(a) bind(c) import :: c_int, c_ptr type(c_ptr), value :: a end function gsl_linalg_ldlt_band_decomp integer(c_int) function gsl_linalg_ldlt_band_solve(ldlt, b, x) bind(c) import :: c_int, c_ptr type(c_ptr), value :: ldlt, b, x end function gsl_linalg_ldlt_band_solve integer(c_int) function gsl_linalg_ldlt_band_svx(ldlt, x) bind(c) import :: c_int, c_ptr type(c_ptr), value :: ldlt, x end function gsl_linalg_ldlt_band_svx integer(c_int) function gsl_linalg_ldlt_band_invert(ldlt, ainv) bind(c) import :: c_int, c_ptr type(c_ptr), value :: ldlt, ainv end function gsl_linalg_ldlt_band_invert integer(c_int) function gsl_linalg_ldlt_band_unpack(ldlt, l, d) bind(c) import :: c_int, c_ptr type(c_ptr), value :: ldlt, l, d end function gsl_linalg_ldlt_band_unpack integer(c_int) function gsl_linalg_ldlt_band_rcond(ldlt, rcond, w) bind(c) import :: c_int, c_ptr, c_double type(c_ptr), value :: ldlt, w real(c_double) :: rcond end function gsl_linalg_ldlt_band_rcond ! ! Balancing ! function gsl_linalg_balance_matrix (a, d) bind(c) import :: c_ptr, c_int type(c_ptr), value :: a, d integer(c_int) :: gsl_linalg_balance_matrix end function gsl_linalg_balance_matrix