minpack_lmstr1 Subroutine

private subroutine minpack_lmstr1(fcn, m, n, x, fvec, fjac, ldfjac, tol, info, ipvt, wa, lwa, udata) bind(c)

Arguments

TypeIntentOptionalAttributesName
procedure(minpack_fcn_lmstr) :: fcn
integer(kind=c_int), value:: m
integer(kind=c_int), value:: n
real(kind=c_double), intent(inout) :: x(n)
real(kind=c_double), intent(out) :: fvec(m)
real(kind=c_double), intent(out) :: fjac(ldfjac,n)
integer(kind=c_int), value:: ldfjac
real(kind=c_double), value:: tol
integer(kind=c_int), intent(out) :: info
integer(kind=c_int), intent(out) :: ipvt(n)
real(kind=c_double), intent(inout) :: wa(lwa)
integer(kind=c_int), value:: lwa
type(c_ptr), value:: udata

Calls

proc~~minpack_lmstr1~~CallsGraph proc~minpack_lmstr1 minpack_lmstr1 proc~lmstr1 lmstr1 proc~minpack_lmstr1->proc~lmstr1 proc~lmstr lmstr proc~lmstr1->proc~lmstr proc~rwupdt rwupdt proc~lmstr->proc~rwupdt proc~lmpar lmpar proc~lmstr->proc~lmpar proc~enorm enorm proc~lmstr->proc~enorm proc~qrfac qrfac proc~lmstr->proc~qrfac proc~lmpar->proc~enorm proc~qrsolv qrsolv proc~lmpar->proc~qrsolv proc~qrfac->proc~enorm

Contents

Source Code


Source Code

    subroutine minpack_lmstr1(fcn, m, n, x, fvec, fjac, ldfjac, tol, info, ipvt, wa, lwa, &
            & udata) &
            & bind(c)
        procedure(minpack_fcn_lmstr) :: fcn
        integer(c_int), value :: m
        integer(c_int), value :: n
        integer(c_int), value :: ldfjac
        integer(c_int), intent(out) :: info
        integer(c_int), value :: lwa
        integer(c_int), intent(out) :: ipvt(n)
        real(c_double), value :: tol
        real(c_double), intent(inout) :: x(n)
        real(c_double), intent(out) :: fvec(m)
        real(c_double), intent(out) :: fjac(ldfjac, n)
        real(c_double), intent(inout) :: wa(lwa)
        type(c_ptr), value :: udata

        call lmstr1(wrap_fcn, m, n, x, fvec, fjac, ldfjac, tol, info, ipvt, wa, lwa)
    contains
        subroutine wrap_fcn(m, n, x, fvec, fjrow, iflag)
            integer, intent(in) :: m
            integer, intent(in) :: n
            integer, intent(inout) :: iflag
            real(wp), intent(in) :: x(n)
            real(wp), intent(inout) :: fvec(m)
            real(wp), intent(inout) :: fjrow(n)

            call fcn(m, n, x, fvec, fjrow, iflag, udata)
        end subroutine wrap_fcn
    end subroutine minpack_lmstr1