the purpose of hybrj is to find a zero of a system of n nonlinear functions in n variables by a modification of the powell hybrid method. 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) | :: | Xtol | a nonnegative input variable. termination occurs when the relative error between two consecutive iterates is at most xtol. |
||
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 and fvec 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:
|
||
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. |
||
real(kind=wp), | intent(out) | :: | r(Lr) | an output array of length lr which contains the upper triangular matrix produced by the qr factorization of the final approximate jacobian, stored rowwise. |
||
integer, | intent(in) | :: | Lr | a positive integer input variable not less than (n*(n+1))/2. |
||
real(kind=wp), | intent(out) | :: | Qtf(n) | an output array of length n which contains 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(n) | work array of length n. |
subroutine hybrj(fcn, n, x, Fvec, Fjac, Ldfjac, Xtol, Maxfev, Diag, Mode, &
Factor, Nprint, Info, Nfev, Njev, r, Lr, Qtf, Wa1, Wa2, &
Wa3, Wa4)
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(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 and fvec 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*** relative error between two consecutive iterates
!! is at most xtol.
!! * ***info = 2*** number of calls to fcn with iflag = 1 has
!! reached maxfev.
!! * ***info = 3*** xtol is too small. no further improvement in
!! the approximate solution x is possible.
!! * ***info = 4*** iteration is not making good progress, as
!! measured by the improvement from the last
!! five jacobian evaluations.
!! * ***info = 5*** iteration is not making good progress, as
!! measured by the improvement from the last
!! ten iterations.
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(in) :: Lr !! a positive integer input variable not less than
!! (n*(n+1))/2.
real(wp), intent(in) :: Xtol !! a nonnegative input variable. termination
!! occurs when the relative error between two consecutive
!! iterates is at most xtol.
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(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) :: 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) :: r(Lr) !! an output array of length lr which contains the
!! upper triangular matrix produced by the qr factorization
!! of the final approximate jacobian, stored rowwise.
real(wp), intent(out) :: Qtf(n) !! an output array of length n which contains
!! 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(n) !! work array of length n.
integer :: i, iflag, iter, j, jm1, l, ncfail, ncsuc, nslow1, nslow2
integer :: iwa(1)
logical :: jeval, sing
real(wp) :: actred, delta, fnorm, fnorm1, pnorm, prered, ratio, sum, temp, xnorm
real(wp), parameter :: p1 = 1.0e-1_wp
real(wp), parameter :: p5 = 5.0e-1_wp
real(wp), parameter :: p001 = 1.0e-3_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 .or. Ldfjac < n .or. Xtol < zero .or. Maxfev <= 0 .or. &
Factor <= zero .or. Lr < (n*(n + 1))/2) exit main
if (Mode == 2) then
do j = 1, n
if (Diag(j) <= zero) exit main
end do
end if
! evaluate the function at the starting point
! and calculate its norm.
iflag = 1
call fcn(n, x, Fvec, Fjac, Ldfjac, iflag)
Nfev = 1
if (iflag < 0) exit main
fnorm = enorm(n, Fvec)
! initialize iteration counter and monitors.
iter = 1
ncsuc = 0
ncfail = 0
nslow1 = 0
nslow2 = 0
! beginning of the outer loop.
outer : do
jeval = .true.
! calculate the jacobian matrix.
iflag = 2
call fcn(n, x, Fvec, Fjac, Ldfjac, iflag)
Njev = Njev + 1
if (iflag < 0) exit main
! compute the qr factorization of the jacobian.
call qrfac(n, n, Fjac, Ldfjac, .false., iwa, 1, 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 in qtf.
do i = 1, n
Qtf(i) = Fvec(i)
end do
do j = 1, n
if (Fjac(j, j) /= zero) then
sum = zero
do i = j, n
sum = sum + Fjac(i, j)*Qtf(i)
end do
temp = -sum/Fjac(j, j)
do i = j, n
Qtf(i) = Qtf(i) + Fjac(i, j)*temp
end do
end if
end do
! copy the triangular factor of the qr factorization into r.
sing = .false.
do j = 1, n
l = j
jm1 = j - 1
if (jm1 >= 1) then
do i = 1, jm1
r(l) = Fjac(i, j)
l = l + n - i
end do
end if
r(l) = Wa1(j)
if (Wa1(j) == zero) sing = .true.
end do
! accumulate the orthogonal factor in fjac.
call qform(n, n, Fjac, Ldfjac, Wa1)
! 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
! if requested, call fcn to enable printing of iterates.
if (Nprint > 0) then
iflag = 0
if (mod(iter - 1, Nprint) == 0) &
call fcn(n, x, Fvec, Fjac, Ldfjac, iflag)
if (iflag < 0) exit main
end if
! determine the direction p.
call dogleg(n, r, Lr, Diag, Qtf, delta, Wa1, Wa2, Wa3)
! 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(n, Wa2, Wa4, Fjac, Ldfjac, iflag)
Nfev = Nfev + 1
if (iflag < 0) exit main
fnorm1 = enorm(n, Wa4)
! compute the scaled actual reduction.
actred = -one
if (fnorm1 < fnorm) actred = one - (fnorm1/fnorm)**2
! compute the scaled predicted reduction.
l = 1
do i = 1, n
sum = zero
do j = i, n
sum = sum + r(l)*Wa1(j)
l = l + 1
end do
Wa3(i) = Qtf(i) + sum
end do
temp = enorm(n, Wa3)
prered = zero
if (temp < fnorm) prered = one - (temp/fnorm)**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 >= p1) then
ncfail = 0
ncsuc = ncsuc + 1
if (ratio >= p5 .or. ncsuc > 1) delta = max(delta, pnorm/p5)
if (abs(ratio - one) <= p1) delta = pnorm/p5
else
ncsuc = 0
ncfail = ncfail + 1
delta = p5*delta
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)
Fvec(j) = Wa4(j)
end do
xnorm = enorm(n, Wa2)
fnorm = fnorm1
iter = iter + 1
end if
! determine the progress of the iteration.
nslow1 = nslow1 + 1
if (actred >= p001) nslow1 = 0
if (jeval) nslow2 = nslow2 + 1
if (actred >= p1) nslow2 = 0
! test for convergence.
if (delta <= Xtol*xnorm .or. fnorm == zero) Info = 1
if (Info /= 0) exit main
! tests for termination and stringent tolerances.
if (Nfev >= Maxfev) Info = 2
if (p1*max(p1*delta, pnorm) <= epsmch*xnorm) Info = 3
if (nslow2 == 5) Info = 4
if (nslow1 == 10) Info = 5
if (Info /= 0) exit main
! criterion for recalculating jacobian.
if (ncfail == 2) cycle outer
! calculate the rank one modification to the jacobian
! and update qtf if necessary.
do j = 1, n
sum = zero
do i = 1, n
sum = sum + Fjac(i, j)*Wa4(i)
end do
Wa2(j) = (sum - Wa3(j))/pnorm
Wa1(j) = Diag(j)*((Diag(j)*Wa1(j))/pnorm)
if (ratio >= p0001) Qtf(j) = sum
end do
! compute the qr factorization of the updated jacobian.
call r1updt(n, n, r, Lr, Wa1, Wa2, Wa3, sing)
call r1mpyq(n, n, Fjac, Ldfjac, Wa2, Wa3)
call r1mpyq(1, n, Qtf, 1, Wa2, Wa3)
jeval = .false.
end do inner ! end of the inner loop.
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(n, x, Fvec, Fjac, Ldfjac, iflag)
end subroutine hybrj