lmder Subroutine

public subroutine lmder(fcn, m, n, x, Fvec, Fjac, Ldfjac, Ftol, Xtol, Gtol, Maxfev, Diag, Mode, Factor, Nprint, Info, Nfev, Njev, Ipvt, Qtf, Wa1, Wa2, Wa3, Wa4)

the purpose of lmder is to minimize the sum of the squares of m nonlinear functions in n variables by a modification of the levenberg-marquardt algorithm. the user must provide a subroutine which calculates the functions and the jacobian.

Arguments

TypeIntentOptionalAttributesName
procedure(fcn_lmder) :: fcn

the 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) :: Ftol

a nonnegative input variable. termination occurs when both the actual and predicted relative reductions in the sum of squares are at most ftol. therefore, ftol measures the relative error desired in the sum of squares.

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

a nonnegative input variable. termination occurs when the relative error between two consecutive iterates is at most xtol. therefore, xtol measures the relative error desired in the approximate solution.

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

a nonnegative input variable. termination occurs when the cosine of the angle between fvec and any column of the jacobian is at most gtol in absolute value. therefore, gtol measures the orthogonality desired between the function vector and the columns of the jacobian.

integer, intent(in) :: Maxfev

a positive integer input variable. termination occurs when the number of calls to fcn with iflag = 1 has reached maxfev.

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

an array of length n. if mode = 1 (see below), diag is internally set. if mode = 2, diag must contain positive entries that serve as multiplicative scale factors for the variables.

integer, intent(in) :: Mode

an integer input variable. if mode = 1, the variables will be scaled internally. if mode = 2, the scaling is specified by the input diag. other values of mode are equivalent to mode = 1.

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

a positive input variable used in determining the initial step bound. this bound is set to the product of factor and the euclidean norm of diag*x if nonzero, or else to factor itself. in most cases factor should lie in the interval (.1,100.).100. is a generally recommended value.

integer, intent(in) :: Nprint

an integer input variable that enables controlled printing of iterates if it is positive. in this case, fcn is called with iflag = 0 at the beginning of the first iteration and every nprint iterations thereafter and immediately prior to return, with x, fvec, and fjac available for printing. fvec and fjac should not be altered. if nprint is not positive, no special calls of fcn with iflag = 0 are made.

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 both actual and predicted relative reductions in the sum of squares are at most ftol.
  • info = 2 relative error between two consecutive iterates is at most xtol.
  • info = 3 conditions for info = 1 and info = 2 both hold.
  • info = 4 the cosine of the angle between fvec and any column of the jacobian is at most gtol in absolute value.
  • info = 5 number of calls to fcn with iflag = 1 has reached maxfev.
  • info = 6 ftol is too small. no further reduction in the sum of squares is possible.
  • info = 7 xtol is too small. no further improvement in the approximate solution x is possible.
  • info = 8 gtol is too small. fvec is orthogonal to the columns of the jacobian to machine precision.
integer, intent(out) :: Nfev

an integer output variable set to the number of calls to fcn with iflag = 1.

integer, intent(out) :: Njev

an integer output variable set to the number of calls to fcn with iflag = 2.

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(out) :: Qtf(n)

an output array of length n which contains the first n elements of the vector (q transpose)*fvec.

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

work array of length n.

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

work array of length n.

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

work array of length n.

real(kind=wp), intent(inout) :: Wa4(m)

work array of length m.


Calls

proc~~lmder~~CallsGraph proc~lmder 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~~lmder~~CalledByGraph proc~lmder lmder proc~minpack_lmder minpack_lmder proc~minpack_lmder->proc~lmder proc~lmder1 lmder1 proc~lmder1->proc~lmder proc~minpack_lmder1 minpack_lmder1 proc~minpack_lmder1->proc~lmder1

Contents

Source Code


Source Code

    subroutine lmder(fcn, m, n, x, Fvec, Fjac, Ldfjac, Ftol, Xtol, Gtol, Maxfev, &
                     Diag, Mode, Factor, Nprint, Info, Nfev, Njev, Ipvt, Qtf, &
                     Wa1, Wa2, Wa3, Wa4)

        implicit none

        procedure(fcn_lmder) :: fcn !! the 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(in) :: Maxfev !! a positive integer input variable. termination
                                 !! occurs when the number of calls to fcn with iflag = 1
                                 !! has reached maxfev.
        integer, intent(in) :: Mode !! an integer input variable. if mode = 1, the
                                !! variables will be scaled internally. if mode = 2,
                                !! the scaling is specified by the input diag. other
                                !! values of mode are equivalent to mode = 1.
        integer, intent(in) :: Nprint !! an integer input variable that enables controlled
                                 !! printing of iterates if it is positive. in this case,
                                 !! fcn is called with iflag = 0 at the beginning of the first
                                 !! iteration and every nprint iterations thereafter and
                                 !! immediately prior to return, with x, fvec, and fjac
                                 !! available for printing. fvec and fjac should not be
                                 !! altered. if nprint is not positive, no special calls
                                 !! of fcn with iflag = 0 are made.
        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***  both actual and predicted relative reductions
                                !!    in the sum of squares are at most ftol.
                                !!  * ***info = 2***  relative error between two consecutive iterates
                                !!    is at most xtol.
                                !!  * ***info = 3***  conditions for info = 1 and info = 2 both hold.
                                !!  * ***info = 4***  the cosine of the angle between fvec and any
                                !!    column of the jacobian is at most gtol in
                                !!    absolute value.
                                !!  * ***info = 5***  number of calls to fcn with iflag = 1 has
                                !!    reached maxfev.
                                !!  * ***info = 6***  ftol is too small. no further reduction in
                                !!    the sum of squares is possible.
                                !!  * ***info = 7***  xtol is too small. no further improvement in
                                !!    the approximate solution x is possible.
                                !!  * ***info = 8***  gtol is too small. fvec is orthogonal to the
                                !!    columns of the jacobian to machine precision.
        integer, intent(out) :: Nfev !! an integer output variable set to the number of
                                !! calls to fcn with iflag = 1.
        integer, intent(out) :: Njev !! an integer output variable set to the number of
                                !! calls to fcn with iflag = 2.
        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) :: Ftol !! a nonnegative input variable. termination
                                !! occurs when both the actual and predicted relative
                                !! reductions in the sum of squares are at most ftol.
                                !! therefore, ftol measures the relative error desired
                                !! in the sum of squares.
        real(wp), intent(in) :: Xtol !! a nonnegative input variable. termination
                                !! occurs when the relative error between two consecutive
                                !! iterates is at most xtol. therefore, xtol measures the
                                !! relative error desired in the approximate solution.
        real(wp), intent(in) :: Gtol !! a nonnegative input variable. termination
                                !! occurs when the cosine of the angle between fvec and
                                !! any column of the jacobian is at most gtol in absolute
                                !! value. therefore, gtol measures the orthogonality
                                !! desired between the function vector and the columns
                                !! of the jacobian.
        real(wp), intent(in) :: Factor !! a positive input variable used in determining the
                                  !! initial step bound. this bound is set to the product of
                                  !! factor and the euclidean norm of diag*x if nonzero, or else
                                  !! to factor itself. in most cases factor should lie in the
                                  !! interval (.1,100.).100. is a generally recommended value.
        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) :: Diag(n) !! an array of length n. if mode = 1 (see
                                      !! below), diag is internally set. if mode = 2, diag
                                      !! must contain positive entries that serve as
                                      !! multiplicative scale factors for the variables.
        real(wp), intent(out) :: Qtf(n) !! an output array of length n which contains
                                   !! the first n elements of the vector (q transpose)*fvec.
        real(wp), intent(inout) :: Wa1(n) !! work array of length n.
        real(wp), intent(inout) :: Wa2(n) !! work array of length n.
        real(wp), intent(inout) :: Wa3(n) !! work array of length n.
        real(wp), intent(inout) :: Wa4(m) !! work array of length m.

        integer :: i, iflag, iter, j, l
        real(wp) :: actred, delta, dirder, fnorm, fnorm1, gnorm, par, &
                    pnorm, prered, ratio, sum, temp, temp1, temp2, xnorm

        real(wp), parameter :: p1 = 1.0e-1_wp
        real(wp), parameter :: p5 = 5.0e-1_wp
        real(wp), parameter :: p25 = 2.5e-1_wp
        real(wp), parameter :: p75 = 7.5e-1_wp
        real(wp), parameter :: p0001 = 1.0e-4_wp

        Info = 0
        iflag = 0
        Nfev = 0
        Njev = 0

        main : block

            ! check the input parameters for errors.

            if (n > 0 .and. m >= n .and. Ldfjac >= m .and. Ftol >= zero .and. &
                Xtol >= zero .and. Gtol >= zero .and. Maxfev > 0 .and. &
                Factor > zero) then
                if (Mode == 2) then
                    do j = 1, n
                        if (Diag(j) <= zero) exit main
                    end do
                end if
            else
                exit main
            end if

            ! evaluate the function at the starting point
            ! and calculate its norm.

            iflag = 1
            call fcn(m, n, x, Fvec, Fjac, Ldfjac, iflag)
            Nfev = 1
            if (iflag < 0) exit main
            fnorm = enorm(m, Fvec)

            ! initialize levenberg-marquardt parameter and iteration counter.

            par = zero
            iter = 1

            ! beginning of the outer loop.

            outer : do

                ! calculate the jacobian matrix.

                iflag = 2
                call fcn(m, n, x, Fvec, Fjac, Ldfjac, iflag)
                Njev = Njev + 1
                if (iflag < 0) exit main

                ! if requested, call fcn to enable printing of iterates.

                if (Nprint > 0) then
                    iflag = 0
                    if (mod(iter - 1, Nprint) == 0) &
                        call fcn(m, n, x, Fvec, Fjac, Ldfjac, iflag)
                    if (iflag < 0) exit main
                end if

                ! compute the qr factorization of the jacobian.

                call qrfac(m, n, Fjac, Ldfjac, .true., Ipvt, n, Wa1, Wa2, Wa3)

                ! on the first iteration and if mode is 1, scale according
                ! to the norms of the columns of the initial jacobian.

                if (iter == 1) then
                    if (Mode /= 2) then
                        do j = 1, n
                            Diag(j) = Wa2(j)
                            if (Wa2(j) == zero) Diag(j) = one
                        end do
                    end if

                    ! on the first iteration, calculate the norm of the scaled x
                    ! and initialize the step bound delta.

                    do j = 1, n
                        Wa3(j) = Diag(j)*x(j)
                    end do
                    xnorm = enorm(n, Wa3)
                    delta = Factor*xnorm
                    if (delta == zero) delta = Factor
                end if

                ! form (q transpose)*fvec and store the first n components in
                ! qtf.

                do i = 1, m
                    Wa4(i) = Fvec(i)
                end do
                do j = 1, n
                    if (Fjac(j, j) /= zero) then
                        sum = zero
                        do i = j, m
                            sum = sum + Fjac(i, j)*Wa4(i)
                        end do
                        temp = -sum/Fjac(j, j)
                        do i = j, m
                            Wa4(i) = Wa4(i) + Fjac(i, j)*temp
                        end do
                    end if
                    Fjac(j, j) = Wa1(j)
                    Qtf(j) = Wa4(j)
                end do

                ! compute the norm of the scaled gradient.

                gnorm = zero
                if (fnorm /= zero) then
                    do j = 1, n
                        l = Ipvt(j)
                        if (Wa2(l) /= zero) then
                            sum = zero
                            do i = 1, j
                                sum = sum + Fjac(i, j)*(Qtf(i)/fnorm)
                            end do
                            gnorm = max(gnorm, abs(sum/Wa2(l)))
                        end if
                    end do
                end if

                ! test for convergence of the gradient norm.

                if (gnorm <= Gtol) Info = 4
                if (Info /= 0) exit main

                ! rescale if necessary.

                if (Mode /= 2) then
                    do j = 1, n
                        Diag(j) = max(Diag(j), Wa2(j))
                    end do
                end if

                ! beginning of the inner loop.
                inner : do

                    ! determine the levenberg-marquardt parameter.

                    call lmpar(n, Fjac, Ldfjac, Ipvt, Diag, Qtf, delta, par, Wa1, Wa2, Wa3, Wa4)

                    ! store the direction p and x + p. calculate the norm of p.

                    do j = 1, n
                        Wa1(j) = -Wa1(j)
                        Wa2(j) = x(j) + Wa1(j)
                        Wa3(j) = Diag(j)*Wa1(j)
                    end do
                    pnorm = enorm(n, Wa3)

                    ! on the first iteration, adjust the initial step bound.

                    if (iter == 1) delta = min(delta, pnorm)

                    ! evaluate the function at x + p and calculate its norm.

                    iflag = 1
                    call fcn(m, n, Wa2, Wa4, Fjac, Ldfjac, iflag)
                    Nfev = Nfev + 1
                    if (iflag < 0) exit main
                    fnorm1 = enorm(m, Wa4)

                    ! compute the scaled actual reduction.

                    actred = -one
                    if (p1*fnorm1 < fnorm) actred = one - (fnorm1/fnorm)**2

                    ! compute the scaled predicted reduction and
                    ! the scaled directional derivative.

                    do j = 1, n
                        Wa3(j) = zero
                        l = Ipvt(j)
                        temp = Wa1(l)
                        do i = 1, j
                            Wa3(i) = Wa3(i) + Fjac(i, j)*temp
                        end do
                    end do
                    temp1 = enorm(n, Wa3)/fnorm
                    temp2 = (sqrt(par)*pnorm)/fnorm
                    prered = temp1**2 + temp2**2/p5
                    dirder = -(temp1**2 + temp2**2)

                    ! compute the ratio of the actual to the predicted
                    ! reduction.

                    ratio = zero
                    if (prered /= zero) ratio = actred/prered

                    ! update the step bound.

                    if (ratio <= p25) then
                        if (actred >= zero) temp = p5
                        if (actred < zero) temp = p5*dirder/(dirder + p5*actred)
                        if (p1*fnorm1 >= fnorm .or. temp < p1) temp = p1
                        delta = temp*min(delta, pnorm/p1)
                        par = par/temp
                    elseif (par == zero .or. ratio >= p75) then
                        delta = pnorm/p5
                        par = p5*par
                    end if

                    ! test for successful iteration.

                    if (ratio >= p0001) then
                        ! successful iteration. update x, fvec, and their norms.
                        do j = 1, n
                            x(j) = Wa2(j)
                            Wa2(j) = Diag(j)*x(j)
                        end do
                        do i = 1, m
                            Fvec(i) = Wa4(i)
                        end do
                        xnorm = enorm(n, Wa2)
                        fnorm = fnorm1
                        iter = iter + 1
                    end if

                    ! tests for convergence.
                    if (abs(actred) <= Ftol .and. prered <= Ftol .and. p5*ratio <= one) Info = 1
                    if (delta <= Xtol*xnorm) Info = 2
                    if (abs(actred) <= Ftol .and. prered <= Ftol .and. p5*ratio <= one .and. Info == 2) Info = 3
                    if (Info /= 0) exit main

                    ! tests for termination and stringent tolerances.
                    if (Nfev >= Maxfev) Info = 5
                    if (abs(actred) <= epsmch .and. prered <= epsmch .and. p5*ratio <= one) Info = 6
                    if (delta <= epsmch*xnorm) Info = 7
                    if (gnorm <= epsmch) Info = 8
                    if (Info /= 0) exit main

                    if (ratio >= p0001) exit inner

                end do inner ! end of the inner loop. repeat if iteration unsuccessful.

            end do outer ! end of the outer loop

        end block main

        ! termination, either normal or user imposed.

        if (iflag < 0) Info = iflag
        iflag = 0
        if (Nprint > 0) call fcn(m, n, x, Fvec, Fjac, Ldfjac, iflag)

    end subroutine lmder