lmder1 Subroutine

public subroutine lmder1(fcn, m, n, x, Fvec, Fjac, Ldfjac, Tol, Info, Ipvt, Wa, Lwa)

the purpose of lmder1 is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm. this is done by using the more general least-squares solver lmder. the user must provide a subroutine which calculates the functions and the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_lmder) :: fcn

user-supplied subroutine which calculates the functions and the jacobian.

integer, intent(in) :: m

a positive integer input variable set to the number of functions.

integer, intent(in) :: n

a positive integer input variable set to the number of variables. n must not exceed m.

real(kind=wp), intent(inout) :: x(n)

an array of length n. on input x must contain an initial estimate of the solution vector. on output x contains the final estimate of the solution vector.

real(kind=wp), intent(out) :: Fvec(m)

an output array of length m which contains the functions evaluated at the output x.

real(kind=wp), intent(out) :: Fjac(Ldfjac,n)

an output m by n array. the upper n by n submatrix of fjac contains an upper triangular matrix r with diagonal elements of nonincreasing magnitude such that

        t     t           t
       p *(jac *jac)*p = r *r,

where p is a permutation matrix and jac is the final calculated jacobian. column j of p is column ipvt(j) (see below) of the identity matrix. the lower trapezoidal part of fjac contains information generated during the computation of r.

integer, intent(in) :: Ldfjac

a positive integer input variable not less than m which specifies the leading dimension of the array fjac.

real(kind=wp), intent(in) :: Tol

a nonnegative input variable. termination occurs when the algorithm estimates either that the relative error in the sum of squares is at most tol or that the relative error between x and the solution is at most tol.

integer, intent(out) :: Info

an integer output variable. if the user has terminated execution, info is set to the (negative) value of iflag. see description of fcn. otherwise, info is set as follows.

  • info = 0 improper input parameters.
  • info = 1 algorithm estimates that the relative error in the sum of squares is at most tol.
  • info = 2 algorithm estimates that the relative error between x and the solution is at most tol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 fvec is orthogonal to the columns of the jacobian to machine precision.
  • info = 5 number of calls to fcn with iflag = 1 has reached 100*(n+1).
  • info = 6 tol is too small. no further reduction in the sum of squares is possible.
  • info = 7 tol is too small. no further improvement in the approximate solution x is possible.
integer, intent(out) :: Ipvt(n)

an integer output array of length n. ipvt defines a permutation matrix p such that jacp = qr, where jac is the final calculated jacobian, q is orthogonal (not stored), and r is upper triangular with diagonal elements of nonincreasing magnitude. column j of p is column ipvt(j) of the identity matrix.

real(kind=wp), intent(inout) :: Wa(Lwa)

a work array of length lwa.

integer, intent(in) :: Lwa

a positive integer input variable not less than 5*n+m.


Calls

proc~~lmder1~~CallsGraph proc~lmder1 lmder1 proc~lmder lmder proc~lmder1->proc~lmder proc~enorm enorm proc~lmder->proc~enorm proc~lmpar lmpar proc~lmder->proc~lmpar proc~qrfac qrfac proc~lmder->proc~qrfac proc~lmpar->proc~enorm proc~qrsolv qrsolv proc~lmpar->proc~qrsolv proc~qrfac->proc~enorm

Called by

proc~~lmder1~~CalledByGraph proc~lmder1 lmder1 proc~minpack_lmder1 minpack_lmder1 proc~minpack_lmder1->proc~lmder1

Contents

Source Code


Source Code

    subroutine lmder1(fcn, m, n, x, Fvec, Fjac, Ldfjac, Tol, Info, Ipvt, Wa, Lwa)
        implicit none

        procedure(fcn_lmder) :: fcn !! user-supplied subroutine which
                                    !! calculates the functions and the jacobian.
        integer, intent(in) :: m !! a positive integer input variable set to the number
                                !! of functions.
        integer, intent(in) :: n !! a positive integer input variable set to the number
                                !! of variables. n must not exceed m.
        integer, intent(in) :: Ldfjac !! a positive integer input variable not less than m
                                     !! which specifies the leading dimension of the array fjac.
        integer, intent(out) :: Info !! an integer output variable. if the user has
                                    !! terminated execution, info is set to the (negative)
                                    !! value of iflag. see description of fcn. otherwise,
                                    !! info is set as follows.
                                    !!
                                    !!  * ***info = 0***  improper input parameters.
                                    !!  * ***info = 1***  algorithm estimates that the relative error
                                    !!    in the sum of squares is at most tol.
                                    !!  * ***info = 2***  algorithm estimates that the relative error
                                    !!    between x and the solution is at most tol.
                                    !!  * ***info = 3***  conditions for info = 1 and info = 2 both hold.
                                    !!  * ***info = 4***  fvec is orthogonal to the columns of the
                                    !!    jacobian to machine precision.
                                    !!  * ***info = 5***  number of calls to fcn with iflag = 1 has
                                    !!    reached 100*(n+1).
                                    !!  * ***info = 6***  tol is too small. no further reduction in
                                    !!    the sum of squares is possible.
                                    !!  * ***info = 7***  tol is too small. no further improvement in
                                    !!    the approximate solution x is possible.
        integer, intent(in) :: Lwa !! a positive integer input variable not less than 5*n+m.
        integer, intent(out) :: Ipvt(n) !! an integer output array of length n. ipvt
                                       !! defines a permutation matrix p such that jac*p = q*r,
                                       !! where jac is the final calculated jacobian, q is
                                       !! orthogonal (not stored), and r is upper triangular
                                       !! with diagonal elements of nonincreasing magnitude.
                                       !! column j of p is column ipvt(j) of the identity matrix.
        real(wp), intent(in) :: Tol !! a nonnegative input variable. termination occurs
                                   !! when the algorithm estimates either that the relative
                                   !! error in the sum of squares is at most tol or that
                                   !! the relative error between x and the solution is at
                                   !! most tol.
        real(wp), intent(inout) :: x(n) !! an array of length n. on input x must contain
                                       !! an initial estimate of the solution vector. on output x
                                       !! contains the final estimate of the solution vector.
        real(wp), intent(out) :: Fvec(m) !! an output array of length m which contains
                                        !! the functions evaluated at the output x.
        real(wp), intent(out) :: Fjac(Ldfjac, n) !! an output m by n array. the upper n by n submatrix
                                                !! of fjac contains an upper triangular matrix r with
                                                !! diagonal elements of nonincreasing magnitude such that
                                                !!```
                                                !!        t     t           t
                                                !!       p *(jac *jac)*p = r *r,
                                                !!```
                                                !! where p is a permutation matrix and jac is the final
                                                !! calculated jacobian. column j of p is column ipvt(j)
                                                !! (see below) of the identity matrix. the lower trapezoidal
                                                !! part of fjac contains information generated during
                                                !! the computation of r.
        real(wp), intent(inout) :: Wa(Lwa) !! a work array of length lwa.

        integer :: maxfev, mode, nfev, njev, nprint
        real(wp) :: ftol, gtol, xtol

        real(wp), parameter :: factor = 100.0_wp

        Info = 0

        ! check the input parameters for errors.

        if (n > 0 .and. m >= n .and. Ldfjac >= m .and. Tol >= zero .and. &
            Lwa >= 5*n + m) then
            ! call lmder.
            maxfev = 100*(n + 1)
            ftol = Tol
            xtol = Tol
            gtol = zero
            mode = 1
            nprint = 0
            call lmder(fcn, m, n, x, Fvec, Fjac, Ldfjac, ftol, xtol, gtol, maxfev,   &
                     & Wa(1), mode, factor, nprint, Info, nfev, njev, Ipvt, Wa(n + 1)&
                     & , Wa(2*n + 1), Wa(3*n + 1), Wa(4*n + 1), Wa(5*n + 1))
            if (Info == 8) Info = 4
        end if

    end subroutine lmder1