UMING
Minimizes a function of N variables using a quasi-Newton method and a user-supplied gradient.
Required Arguments
FCN — User-supplied subroutine to evaluate the function to be minimized. The usage is CALL FCN (NXF), where
N – Length of X. (Input)
X – Vector of length N at which point the function is evaluated. (Input)
X should not be changed by FCN.
F – The computed function value at the point X. (Output)
FCN must be declared EXTERNAL in the calling program.
GRAD — User-supplied subroutine to compute the gradient at the point X. The usage is CALL GRAD (NXG), where
N – Length of X and G. (Input)
X – Vector of length N at which point the function is evaluated. (Input)
X should not be changed by GRAD .
G – The gradient evaluated at the point X. (Output)
GRAD must be declared EXTERNAL in the calling program.
X — Vector of length N containing the computed solution. (Output)
Optional Arguments
N — Dimension of the problem. (Input)
Default: N = SIZE (X,1).
XGUESS — Vector of length N containing the initial guess of the minimum. (Input)
Default: XGUESS = 0.0.
XSCALE — Vector of length N containing the diagonal scaling matrix for the variables. (Input)
XSCALE is used mainly in scaling the gradient and the distance between two points. In the absence of other information, set all entries to 1.0.
Default: XSCALE = 1.0.
FSCALE — Scalar containing the function scaling. (Input)
FSCALE is used mainly in scaling the gradient. In the absence of other information, set FSCALE to 1.0.
Default: FSCALE = 1.0.
IPARAM — Parameter vector of length 7. (Input/Output)
Set IPARAM(1) to zero for default values of IPARAM and RPARAM. See Comment 4.
Default: IPARAM = 0.
RPARAM — Parameter vector of length 7. (Input/Output)
See Comment 4.
FVALUE — Scalar containing the value of the function at the computed solution. (Output)
FORTRAN 90 Interface
Generic: CALL UMING (FCN, GRAD, X [])
Specific: The specific interface names are S_UMING and D_UMING.
FORTRAN 77 Interface
Single: CALL UMING (FCN, GRAD, N, XGUESS, XSCALE, FSCALE, IPARAM, RPARAM, X, FVALUE)
Double: The double precision name is DUMING.
Description
The routine UMING uses a quasi-Newton method to find the minimum of a function f(x) of n variables. Function values and first derivatives are required. The problem is stated as follows:
Given a starting point xc, the search direction is computed according to the formula
d = -B-1 gc
where B is a positive definite approximation of the Hessian and gc is the gradient evaluated at xc. A line search is then used to find a new point
xn = xc + λd,            λ > 0
such that
f(xn) ≤ f(xc) + αgT dα  (0, 0.5)
Finally, the optimality condition g(x) = ɛ is checked where ɛ is a gradient tolerance.
When optimality is not achieved, B is updated according to the BFGS formula
where s = xn  xc and y = gn  gc. Another search direction is then computed to begin the next iteration. For more details, see Dennis and Schnabel (1983, Appendix A).
Comments
1. Workspace may be explicitly provided, if desired, by use of U2ING/DU2ING. The reference is:
CALL U2ING (FCN, GRAD, N, XGUESS, XSCALE, FSCALE, IPARAM, RPARAM, X, FVALUE, WK)
The additional argument is
WK — Work vector of length N * (N + 8). WK contains the following information on output: The second N locations contain the last step taken. The third N locations contain the last Newton step. The fourth N locations contain an estimate of the gradient at the solution. The final N2 locations contain the Cholesky factorization of a BFGS approximation to the Hessian at the solution.
2. Informational errors
Type
Code
Description
4
2
The iterates appear to be converging to a noncritical point.
4
3
Maximum number of iterations exceeded.
4
4
Maximum number of function evaluations exceeded.
4
5
Maximum number of gradient evaluations exceeded.
4
6
Five consecutive steps have been taken with the maximum step length.
2
7
Scaled step tolerance satisfied; the current point may be an approximate local solution, or the algorithm is making very slow progress and is not near a solution, or STEPTL is too big.
3
8
The last global step failed to locate a lower point than the current X value.
3. The first stopping criterion for UMING occurs when the infinity norm of the scaled gradient is less than the given gradient tolerance (RPARAM(1)). The second stopping criterion for UMING occurs when the scaled distance between the last two steps is less than the step tolerance (RPARAM(2)).
4. If the default parameters are desired for UMING, then set IPARAM(1) to zero and call routine UMING. Otherwise, if any nondefault parameters are desired for IPARAM or RPARAM, then the following steps should be taken before calling UMING:
CALL U4INF (IPARAM, RPARAM)
Set nondefault values for desired IPARAM, RPARAM elements.
Note that the call to U4INF will set IPARAM and RPARAM to their default values so only nondefault values need to be set above.
The following is a list of the parameters and the default values:
IPARAM — Integer vector of length 7.
IPARAM(1) = Initialization flag.
IPARAM(2) = Number of good digits in the function.
Default: Machine dependent.
IPARAM(3) = Maximum number of iterations.
Default: 100.
IPARAM(4) = Maximum number of function evaluations.
Default: 400.
IPARAM(5) = Maximum number of gradient evaluations.
Default: 400.
IPARAM(6) = Hessian initialization parameter
If IPARAM(6) = 0, the Hessian is initialized to the identity matrix; otherwise, it is initialized to a diagonal matrix containing
on the diagonal where t = XGUESS, fs = FSCALE, and s = XSCALE.
Default: 0.
IPARAM(7) = Maximum number of Hessian evaluations.
Default: Not used in UMING.
RPARAM — Real vector of length 7.
RPARAM(1) = Scaled gradient tolerance.
The i-th component of the scaled gradient at x is calculated as
where g = f (x), s = XSCALE, and fs = FSCALE.
Default:
in double where ɛ is the machine precision.
RPARAM(2) = Scaled step tolerance. (STEPTL)
The i-th component of the scaled step between two points x and y is computed as
where s = XSCALE.
Default: ɛ2/3 where ɛ is the machine precision.
RPARAM(3) = Relative function tolerance.
Default: Not used in UMING.
RPARAM(4) = Absolute function tolerance.
Default: Not used in UMING.
RPARAM(5) = False convergence tolerance.
Default: Not used in UMING.
RPARAM(6) = Maximum allowable step size.
Default: 1000 max(ɛ1, ɛ2) where
ɛ2 = s2s = XSCALE, and t = XGUESS.
RPARAM(7) = Size of initial trust region radius.
Default: Not used in UMING.
If double precision is required, then DU4INF is called, and RPARAM is declared double precision.
5. Users wishing to override the default print/stop attributes associated with error messages issued by this routine are referred to “Error Handling” in the Introduction.
Example
The function
is minimized. Default values for parameters are used.
 
USE UMING_INT
USE UMACH_INT
 
IMPLICIT NONE
INTEGER N
PARAMETER (N=2)
!
INTEGER IPARAM(7), L, NOUT
REAL F, X(N), XGUESS(N)
EXTERNAL ROSBRK, ROSGRD
!
DATA XGUESS/-1.2E0, 1.0E0/
!
IPARAM(1) = 0
! Minimize Rosenbrock function using
! initial guesses of -1.2 and 1.0
CALL UMING (ROSBRK, ROSGRD, X, XGUESS=XGUESS, IPARAM=IPARAM, FVALUE=F)
! Print results
CALL UMACH (2, NOUT)
WRITE (NOUT,99999) X, F, (IPARAM(L),L=3,5)
!
99999 FORMAT (' The solution is ', 6X, 2F8.3, //, ' The function ', &
'value is ', F8.3, //, ' The number of iterations is ', &
10X, I3, /, ' The number of function evaluations is ', &
I3, /, ' The number of gradient evaluations is ', I3)
!
END
!
SUBROUTINE ROSBRK (N, X, F)
INTEGER N
REAL X(N), F
!
F = 1.0E2*(X(2)-X(1)*X(1))**2 + (1.0E0-X(1))**2
!
RETURN
END
!
SUBROUTINE ROSGRD (N, X, G)
INTEGER N
REAL X(N), G(N)
!
G(1) = -4.0E2*(X(2)-X(1)*X(1))*X(1) - 2.0E0*(1.0E0-X(1))
G(2) = 2.0E2*(X(2)-X(1)*X(1))
!
RETURN
END
Output
 
The solution is 1.000 1.000
 
The function value is 0.000
 
The number of iterations is 18
The number of function evaluations is 31
The number of gradient evaluations is 22