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.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(fcn_hybrj) | :: | fcn | the user-supplied subroutine which calculates the functions and the jacobian |
|||
integer, | intent(in) | :: | n | a positive integer input variable set to the number of functions and variables. |
||
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(n) | an output array of length n which contains the functions evaluated at the output x. |
||
real(kind=wp), | intent(out) | :: | Fjac(Ldfjac,n) | an output n by n array which contains the orthogonal matrix q produced by the qr factorization of the final approximate jacobian. |
||
integer, | intent(in) | :: | Ldfjac | a positive integer input variable not less than n 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 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:
|
||
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 (n*(n+13))/2. |
subroutine hybrj1(fcn, n, x, Fvec, Fjac, Ldfjac, Tol, Info, Wa, Lwa)
implicit none
procedure(fcn_hybrj) :: fcn !! the user-supplied subroutine which
!! calculates the functions and the jacobian
integer, intent(in) :: n !! a positive integer input variable set to the number
!! of functions and variables.
integer, intent(in) :: Ldfjac !! a positive integer input variable not less than n
!! 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
!! 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.
integer, intent(in) :: Lwa !! a positive integer input variable not less than
!! (n*(n+13))/2.
real(wp), intent(in) :: Tol !! a nonnegative input variable. termination occurs
!! when the algorithm estimates 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(n) !! an output array of length n which contains
!! the functions evaluated at the output x.
real(wp), intent(out) :: Fjac(Ldfjac, n) !! an output n by n array which contains the
!! orthogonal matrix q produced by the qr factorization
!! of the final approximate jacobian.
real(wp), intent(inout) :: Wa(Lwa) !! a work array of length 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
end do
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
end if
end subroutine hybrj1