hybrj1 Subroutine

public subroutine hybrj1(fcn, n, x, fvec, fjac, ldfjac, tol, info, wa, lwa)

the purpose of hybrj1 is to find a zero of a system of n nonlinear functions in n variables by a modification of the powell hybrid method. this is done by using the more general nonlinear equation solver hybrj. the user must provide a subroutine which calculates the functions and the jacobian.

the subroutine statement is

subroutine hybrj1(fcn,n,x,fvec,fjac,ldfjac,tol,info,wa,lwa)

where

fcn is the name of the user-supplied subroutine which calculates the functions and the jacobian. fcn must be declared in an external statement in the user calling program, and should be written as follows.

subroutine fcn(n,x,fvec,fjac,ldfjac,iflag)
integer n,ldfjac,iflag
real(wp) x(n),fvec(n),fjac(ldfjac,n)
----------
if iflag = 1 calculate the functions at x and
return this vector in fvec. do not alter fjac.
if iflag = 2 calculate the jacobian at x and
return this matrix in fjac. do not alter fvec.
---------
return
end

the value of iflag should not be changed by fcn unless
the user wants to terminate execution of hybrj1.
in this case set iflag to a negative integer.

n is a positive integer input variable set to the number of functions and variables.

x is 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.

fvec is an output array of length n which contains the functions evaluated at the output x.

fjac is an output n by n array which contains the orthogonal matrix q produced by the qr factorization of the final approximate jacobian.

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

tol is a nonnegative input variable. termination occurs when the algorithm estimates that the relative error between x and the solution is at most tol.

info is 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
           between x and the solution is at most tol.

info = 2   number of calls to fcn with iflag = 1 has
           reached 100*(n+1).

info = 3   tol is too small. no further improvement in
           the approximate solution x is possible.

info = 4   iteration is not making good progress.

wa is a work array of length lwa.

lwa is a positive integer input variable not less than (n*(n+13))/2.

Arguments

Type IntentOptional Attributes Name
procedure(fcn_hybrj) :: fcn
integer :: n
real(kind=wp) :: x(n)
real(kind=wp) :: fvec(n)
real(kind=wp) :: fjac(ldfjac,n)
integer :: ldfjac
real(kind=wp) :: tol
integer :: info
real(kind=wp) :: wa(lwa)
integer :: lwa

Calls

proc~~hybrj1~~CallsGraph proc~hybrj1 minpack_module::hybrj1 proc~hybrj minpack_module::hybrj proc~hybrj1->proc~hybrj proc~dogleg minpack_module::dogleg proc~hybrj->proc~dogleg proc~dpmpar minpack_module::dpmpar proc~hybrj->proc~dpmpar proc~enorm minpack_module::enorm proc~hybrj->proc~enorm proc~qform minpack_module::qform proc~hybrj->proc~qform proc~qrfac minpack_module::qrfac proc~hybrj->proc~qrfac proc~r1mpyq minpack_module::r1mpyq proc~hybrj->proc~r1mpyq proc~r1updt minpack_module::r1updt proc~hybrj->proc~r1updt proc~dogleg->proc~dpmpar proc~dogleg->proc~enorm proc~qrfac->proc~dpmpar proc~qrfac->proc~enorm proc~r1updt->proc~dpmpar

Source Code

    subroutine hybrj1(fcn,n,x,fvec,fjac,ldfjac,tol,info,wa,lwa)

    implicit none

    procedure(fcn_hybrj) :: fcn
    integer n , ldfjac , info , lwa
    real(wp) tol
    real(wp) x(n) , fvec(n) , fjac(ldfjac,n) , wa(lwa)

      integer j , lr , maxfev , mode , nfev , njev , nprint
      real(wp) xtol

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

      info = 0

      ! check the input parameters for errors.

      if ( n>0 .and. ldfjac>=n .and. tol>=zero .and. lwa>=(n*(n+13))/2 ) then

         ! call hybrj.

         maxfev = 100*(n+1)
         xtol = tol
         mode = 2
         do j = 1 , n
            wa(j) = one
         enddo
         nprint = 0
         lr = (n*(n+1))/2
         call hybrj(fcn,n,x,fvec,fjac,ldfjac,xtol,maxfev,wa(1),mode,    &
                    factor,nprint,info,nfev,njev,wa(6*n+1),lr,wa(n+1),  &
                    wa(2*n+1),wa(3*n+1),wa(4*n+1),wa(5*n+1))
         if ( info==5 ) info = 4

      endif

      end subroutine hybrj1