!-*-f90-*- ! ! Interfaces: nonlinear least squares fitting ! type(c_ptr) function gsl_multifit_nlinear_setup(s) BIND(C) import character(c_char) :: s end function type(c_ptr) function gsl_multilarge_nlinear_setup(s) BIND(C) import character(c_char) :: s end function function gsl_multifit_nlinear_alloc(t, params, n, p) BIND(C) import type(c_ptr), value :: t type(gsl_multifit_nlinear_parameters) :: params type(c_ptr) :: gsl_multifit_nlinear_alloc integer(c_size_t), value :: n, p end function function gsl_multilarge_nlinear_alloc(t, params, n, p) BIND(C) import type(c_ptr), value :: t type(gsl_multilarge_nlinear_parameters) :: params type(c_ptr) :: gsl_multilarge_nlinear_alloc integer(c_size_t), value :: n, p end function function gsl_multifit_nlinear_default_parameters() BIND(C) import type(gsl_multifit_nlinear_parameters) :: gsl_multifit_nlinear_default_parameters end function function gsl_multilarge_nlinear_default_parameters() BIND(C) import type(gsl_multilarge_nlinear_parameters) :: gsl_multilarge_nlinear_default_parameters end function integer(c_int) function gsl_multifit_nlinear_init(x, fdf, w) BIND(C) import type(c_ptr), value :: x, fdf, w end function integer(c_int) function gsl_multilarge_nlinear_init(x, fdf, w) BIND(C) import type(c_ptr), value :: x, fdf, w end function integer(c_int) function gsl_multifit_nlinear_winit(x, wts, fdf, w) BIND(C) import type(c_ptr), value :: x, wts, fdf, w end function integer(c_int) function gsl_multilarge_nlinear_winit(x, wts, fdf, w) BIND(C) import type(c_ptr), value :: x, wts, fdf, w end function subroutine gsl_multifit_nlinear_free(w) BIND(C) import type(c_ptr), value :: w end subroutine subroutine gsl_multilarge_nlinear_free(w) BIND(C) import type(c_ptr), value :: w end subroutine function gsl_multifit_nlinear_name(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multifit_nlinear_name end function function gsl_multilarge_nlinear_name(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multilarge_nlinear_name end function function gsl_multifit_nlinear_trs_name(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multifit_nlinear_trs_name end function function gsl_multilarge_nlinear_trs_name(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multilarge_nlinear_trs_name end function integer(c_int) function gsl_multifit_nlinear_iterate(w) BIND(C) import type(c_ptr), value :: w end function integer(c_int) function gsl_multilarge_nlinear_iterate(w) BIND(C) import type(c_ptr), value :: w end function function gsl_multifit_nlinear_position(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multifit_nlinear_position end function function gsl_multilarge_nlinear_position(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multilarge_nlinear_position end function function gsl_multifit_nlinear_residual(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multifit_nlinear_residual end function function gsl_multilarge_nlinear_residual(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multilarge_nlinear_residual end function function gsl_multifit_nlinear_jac(w) BIND(C) import type(c_ptr), value :: w type(c_ptr) :: gsl_multifit_nlinear_jac end function integer (c_int) function gsl_multifit_nlinear_niter(w) BIND(C) import type(c_ptr), value :: w end function integer (c_int) function gsl_multilarge_nlinear_niter(w) BIND(C) import type(c_ptr), value :: w end function integer (c_int) function gsl_multifit_nlinear_rcond(rcond, w) BIND(C) import real(c_double), intent(inout) :: rcond type(c_ptr), value :: w end function integer (c_int) function gsl_multilarge_nlinear_rcond(rcond, w) BIND(C) import real(c_double), intent(inout) :: rcond type(c_ptr), value :: w end function integer(c_int) function gsl_multifit_nlinear_test(xtol, gtol, ftol, info, w) BIND(C) import :: c_ptr, c_double, c_int real(c_double), value :: xtol, gtol, ftol integer(c_int), intent(inout) :: info type(c_ptr), value :: w end function integer(c_int) function gsl_multilarge_nlinear_test(xtol, gtol, ftol, info, w) BIND(C) import :: c_ptr, c_double, c_int real(c_double), value :: xtol, gtol, ftol integer(c_int), intent(inout) :: info type(c_ptr), value :: w end function integer(c_int) function gsl_multifit_nlinear_driver(maxiter, xtol, gtol, ftol, & callback, callback_params, info, w) BIND(C) import :: c_ptr, c_double, c_int, c_size_t, c_funptr integer(c_size_t), value :: maxiter real(c_double), value :: xtol, gtol, ftol type(c_funptr), value :: callback type(c_ptr), value :: callback_params, w integer(c_int), intent(inout) :: info end function integer(c_int) function gsl_multilarge_nlinear_driver(maxiter, xtol, gtol, ftol, & callback, callback_params, info, w) BIND(C) import :: c_ptr, c_double, c_int, c_size_t, c_funptr integer(c_size_t), value :: maxiter real(c_double), value :: xtol, gtol, ftol type(c_funptr), value :: callback type(c_ptr), value :: callback_params, w integer(c_int), intent(inout) :: info end function integer(c_int) function gsl_multifit_nlinear_covar(j, epsrel, covar) BIND(C) import type(c_ptr), value :: j, covar real(c_double), value :: epsrel end function integer(c_int) function gsl_multilarge_nlinear_covar(covar, w) BIND(C) import type(c_ptr), value :: covar, w end function type(c_ptr) function fgsl_multifit_nlinear_fdf_cinit(ndim, p, params, fp, dfp, fvvp) BIND(C) import integer(c_size_t), value :: ndim, p type(c_ptr), value :: params type(c_funptr), value :: fp, dfp, fvvp end function subroutine gsl_multifit_nlinear_fdf_get(fdf, fp, dfp, fvvp, & n, p, params, nevalf, nevaldf, nevalfvv) BIND(C) import type(c_ptr), value :: fdf type(c_funptr) :: fp, dfp, fvvp integer(c_size_t) :: n, p, nevalf, nevaldf, nevalfvv type(c_ptr) :: params end subroutine subroutine fgsl_multifit_nlinear_fdf_cfree(fdf) BIND(C) import type(c_ptr), value :: fdf end subroutine type(c_ptr) function gsl_multifit_nlinear_get_trs(which) BIND(C) import integer(c_int), value :: which end function type(c_ptr) function gsl_multifit_nlinear_get_scale(which) BIND(C) import integer(c_int), value :: which end function type(c_ptr) function gsl_multifit_nlinear_get_solver(which) BIND(C) import integer(c_int), value :: which end function type(c_ptr) function fgsl_multilarge_nlinear_fdf_cinit(ndim, p, params, fp, dfp, fvvp) BIND(C) import integer(c_size_t), value :: ndim, p type(c_ptr), value :: params type(c_funptr), value :: fp, dfp, fvvp end function subroutine fgsl_multilarge_nlinear_fdf_cfree(fdf) BIND(C) import type(c_ptr), value :: fdf end subroutine subroutine gsl_multilarge_nlinear_fdf_get(fdf, fp, dfp, fvvp, & n, p, params, nevalf, nevaldfu, nevaldf2, nevalfvv) BIND(C) import type(c_ptr), value :: fdf type(c_funptr) :: fp, dfp, fvvp integer(c_size_t) :: n, p, nevalf, nevaldfu, nevaldf2, nevalfvv type(c_ptr) :: params end subroutine type(c_ptr) function gsl_multilarge_nlinear_get_trs(which) BIND(C) import integer(c_int), value :: which end function type(c_ptr) function gsl_multilarge_nlinear_get_scale(which) BIND(C) import integer(c_int), value :: which end function type(c_ptr) function gsl_multilarge_nlinear_get_solver(which) BIND(C) import integer(c_int), value :: which end function