minpack_hybrd1 Subroutine

public subroutine minpack_hybrd1(fcn, n, x, fvec, tol, info, Wa, Lwa, udata) bind(c)

Arguments

TypeIntentOptionalAttributesName
procedure(minpack_func) :: fcn
integer(kind=c_int), value:: n
real(kind=c_double), intent(inout) :: x(n)
real(kind=c_double), intent(out) :: fvec(n)
real(kind=c_double), value:: tol
integer(kind=c_int), intent(out) :: info
real(kind=c_double), intent(inout) :: Wa(Lwa)
integer(kind=c_int), value:: Lwa
type(c_ptr), value:: udata

Calls

proc~~minpack_hybrd1~~CallsGraph proc~minpack_hybrd1 minpack_hybrd1 proc~hybrd1 hybrd1 proc~minpack_hybrd1->proc~hybrd1 proc~hybrd hybrd proc~hybrd1->proc~hybrd proc~dogleg dogleg proc~hybrd->proc~dogleg proc~r1updt r1updt proc~hybrd->proc~r1updt proc~qrfac qrfac proc~hybrd->proc~qrfac proc~qform qform proc~hybrd->proc~qform proc~enorm enorm proc~hybrd->proc~enorm proc~r1mpyq r1mpyq proc~hybrd->proc~r1mpyq proc~fdjac1 fdjac1 proc~hybrd->proc~fdjac1 proc~dogleg->proc~enorm proc~qrfac->proc~enorm

Contents

Source Code


Source Code

    subroutine minpack_hybrd1(fcn, n, x, Fvec, Tol, Info, Wa, Lwa, udata) &
          & bind(c)
        procedure(minpack_func) :: fcn
        integer(c_int), value :: n
        integer(c_int), intent(out) :: info
        real(c_double), value :: tol
        real(c_double), intent(inout) :: x(n)
        real(c_double), intent(out) :: fvec(n)
        integer(c_int), value :: Lwa
        real(c_double), intent(inout) :: Wa(Lwa)
        type(c_ptr), value :: udata

        call hybrd1(wrap_fcn, n, x, fvec, tol, info, Wa, Lwa)

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

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