!-*-f90-*- ! ! Interfaces: large linear least squares systems ! function gsl_multilarge_linear_alloc(T, p) bind(c) import :: c_ptr, c_size_t type(c_ptr), value :: T integer(c_size_t), value :: p type(c_ptr) :: gsl_multilarge_linear_alloc end function gsl_multilarge_linear_alloc subroutine gsl_multilarge_linear_free(w) bind(c) import :: c_ptr type(c_ptr), value :: w end subroutine gsl_multilarge_linear_free function gsl_multilarge_linear_name(w) bind(c) import :: c_ptr type(c_ptr), value :: w type(c_ptr) :: gsl_multilarge_linear_name end function gsl_multilarge_linear_name function gsl_multilarge_linear_reset(w) bind(c) import :: c_ptr, c_int type(c_ptr), value :: w integer(c_int) :: gsl_multilarge_linear_reset end function gsl_multilarge_linear_reset function gsl_multilarge_linear_accumulate(X, y, w) bind(c) import :: c_ptr, c_int type(c_ptr), value :: X, y, w integer(c_int) :: gsl_multilarge_linear_accumulate end function gsl_multilarge_linear_accumulate function gsl_multilarge_linear_solve(lambda, c, rnorm, snorm, w) bind(c) import :: c_int, c_ptr, c_double real(c_double), value :: lambda type(c_ptr), value :: c, w real(c_double) :: rnorm, snorm integer(c_int) :: gsl_multilarge_linear_solve end function gsl_multilarge_linear_solve function gsl_multilarge_linear_rcond(rcond, w) bind(c) import :: c_ptr, c_int, c_double real(c_double) :: rcond type(c_ptr), value :: w integer(c_int) :: gsl_multilarge_linear_rcond end function gsl_multilarge_linear_rcond function gsl_multilarge_linear_lcurve(reg_param, rho, eta, w) bind(c) import :: c_ptr, c_int type(c_ptr), value :: reg_param, rho, eta, w integer(c_int) :: gsl_multilarge_linear_lcurve end function gsl_multilarge_linear_lcurve function gsl_multilarge_linear_wstdform1(L, X, w, y, Xs, ys, work) bind(c) import :: c_ptr, c_int type(c_ptr), value :: L, X, w, y, Xs, ys, work integer(c_int) :: gsl_multilarge_linear_wstdform1 end function gsl_multilarge_linear_wstdform1 function gsl_multilarge_linear_stdform1(L, X, y, Xs, ys, work) bind(c) import :: c_ptr, c_int type(c_ptr), value :: L, X, y, Xs, ys, work integer(c_int) :: gsl_multilarge_linear_stdform1 end function gsl_multilarge_linear_stdform1 function gsl_multilarge_linear_l_decomp(L, tau) & bind(c, name='gsl_multilarge_linear_L_decomp') import :: c_ptr, c_int type(c_ptr), value :: L, tau integer(c_int) :: gsl_multilarge_linear_l_decomp end function gsl_multilarge_linear_l_decomp function gsl_multilarge_linear_wstdform2(LQR, Ltau, X, w, y, Xs, ys, work) bind(c) import :: c_ptr, c_int type(c_ptr), value :: LQR, Ltau, X, w, y, Xs, ys, work integer(c_int) :: gsl_multilarge_linear_wstdform2 end function gsl_multilarge_linear_wstdform2 function gsl_multilarge_linear_stdform2(LQR, Ltau, X, y, Xs, ys, work) bind(c) import :: c_ptr, c_int type(c_ptr), value :: LQR, Ltau, X, y, Xs, ys, work integer(c_int) :: gsl_multilarge_linear_stdform2 end function gsl_multilarge_linear_stdform2 function gsl_multilarge_linear_genform1(L, cs, c, work) bind(c) import :: c_ptr, c_int type(c_ptr), value :: L, cs, c, work integer(c_int) :: gsl_multilarge_linear_genform1 end function gsl_multilarge_linear_genform1 function gsl_multilarge_linear_genform2(LQR, Ltau, cs, c, work) bind(c) import :: c_ptr, c_int type(c_ptr), value :: LQR, Ltau, cs, c, work integer(c_int) :: gsl_multilarge_linear_genform2 end function gsl_multilarge_linear_genform2 function fgsl_aux_multilarge_linear_alloc(i) bind(c) import :: c_ptr, c_int integer(c_int), value :: i type(c_ptr) :: fgsl_aux_multilarge_linear_alloc end function fgsl_aux_multilarge_linear_alloc