minpack_lmdif Subroutine

public subroutine minpack_lmdif(fcn, m, n, x, fvec, ftol, xtol, gtol, maxfev, epsfcn, diag, mode, factor, nprint, info, nfev, fjac, ldfjac, ipvt, qtf, wa1, wa2, wa3, wa4, udata) bind(c)

Arguments

TypeIntentOptionalAttributesName
procedure(minpack_func2) :: 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), value:: ftol
real(kind=c_double), value:: xtol
real(kind=c_double), value:: gtol
integer(kind=c_int), value:: maxfev
real(kind=c_double), value:: epsfcn
real(kind=c_double), intent(inout) :: diag(n)
integer(kind=c_int), value:: mode
real(kind=c_double), value:: factor
integer(kind=c_int), value:: nprint
integer(kind=c_int), intent(out) :: info
integer(kind=c_int), intent(out) :: nfev
real(kind=c_double), intent(out) :: fjac(ldfjac,n)
integer(kind=c_int), value:: ldfjac
integer(kind=c_int), intent(out) :: ipvt(n)
real(kind=c_double), intent(out) :: qtf(n)
real(kind=c_double), intent(inout) :: wa1(n)
real(kind=c_double), intent(inout) :: wa2(n)
real(kind=c_double), intent(inout) :: wa3(n)
real(kind=c_double), intent(inout) :: wa4(m)
type(c_ptr), value:: udata

Calls

proc~~minpack_lmdif~~CallsGraph proc~minpack_lmdif minpack_lmdif proc~lmdif lmdif proc~minpack_lmdif->proc~lmdif proc~enorm enorm proc~lmdif->proc~enorm proc~lmpar lmpar proc~lmdif->proc~lmpar proc~fdjac2 fdjac2 proc~lmdif->proc~fdjac2 proc~qrfac qrfac proc~lmdif->proc~qrfac proc~lmpar->proc~enorm proc~qrsolv qrsolv proc~lmpar->proc~qrsolv proc~qrfac->proc~enorm

Contents

Source Code


Source Code

    subroutine minpack_lmdif(fcn, m, n, x, fvec, ftol, xtol, gtol, maxfev, epsfcn, diag, &
            & mode, factor, nprint, info, nfev, fjac, ldfjac, ipvt, qtf, wa1, wa2, wa3, wa4, &
            & udata) &
            & bind(c)
        procedure(minpack_func2) :: fcn
        integer(c_int), value :: m
        integer(c_int), value :: n
        integer(c_int), value :: maxfev
        integer(c_int), value :: mode
        integer(c_int), value :: nprint
        integer(c_int), intent(out) :: info
        integer(c_int), intent(out) :: nfev
        integer(c_int), value :: ldfjac
        integer(c_int), intent(out) :: ipvt(n)
        real(c_double), value :: ftol
        real(c_double), value :: xtol
        real(c_double), value :: gtol
        real(c_double), value :: epsfcn
        real(c_double), value :: factor
        real(c_double), intent(inout) :: x(n)
        real(c_double), intent(out) :: fvec(m)
        real(c_double), intent(inout) :: diag(n)
        real(c_double), intent(out) :: fjac(ldfjac, n)
        real(c_double), intent(out) :: qtf(n)
        real(c_double), intent(inout) :: wa1(n)
        real(c_double), intent(inout) :: wa2(n)
        real(c_double), intent(inout) :: wa3(n)
        real(c_double), intent(inout) :: wa4(m)
        type(c_ptr), value :: udata

        call lmdif(wrap_fcn, m, n, x, fvec, ftol, xtol, gtol, maxfev, epsfcn, diag, &
            & mode, factor, nprint, info, nfev, fjac, ldfjac, ipvt, qtf, wa1, wa2, wa3, wa4)

    contains
        subroutine wrap_fcn(m, n, x, fvec, iflag)
            integer, intent(in) :: m
            integer, intent(in) :: n
            real(wp), intent(in) :: x(n)
            real(wp), intent(out) :: fvec(m)
            integer, intent(inout) :: iflag

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