e04kd {NAGFWrappers}R Documentation

e04kd: Minimum, function of several variables, modified Newton algorithm, simple bounds, using first derivatives (comprehensive)

Description

e04kd is a comprehensive modified Newton algorithm for finding:

- an unconstrained minimum of a function of several variables;

- a minimum of a function of several variables subject to fixed upper and/or lower bounds on the variables.

First derivatives are required. The function is intended for functions which have continuous first and second derivatives (although it will usually work even if the derivatives have occasional discontinuities).

Usage

e04kd(funct, monit, eta, ibound, bl, bu, x, lh, iw, w,
n=nrow(bl),
iprint=1,
maxcal=50,
xtol=0.0,
delta=0.0,
stepmx=100000.0)

Arguments

funct

void function

funct must evaluate the function F(x) and its first derivatives ( \partial F)/( \partial x_j) at a specified point. (However, if you do not wish to calculate F or its first derivatives at a particular x, there is the option of setting a argument to cause e04kd to terminate immediately.)

monit

void function

If iprint >= 0, you must supply monit which is suitable for monitoring the minimization process. monit must not change the values of any of its arguments.

eta

double

Every iteration of e04kd involves a linear minimization (i.e., minimization of F(x + αp) with respect to α). eta specifies how accurately these linear minimizations are to be performed. The minimum with respect to α will be located more accurately for small values of eta (say, 0.01) than large values (say, 0.9).

ibound

integer

Indicates whether the problem is unconstrained or bounded. If there are bounds on the variables, ibound can be used to indicate whether the facility for dealing with bounds of special forms is to be used. It must be set to one of the following values:

ibound = 0

: If the variables are bounded and you are supplying all the l_j and u_j individually.

ibound = 1

: If the problem is unconstrained.

ibound = 2

: If the variables are bounded, but all the bounds are of the form 0 <= x_j.

ibound = 3

: If all the variables are bounded, and l_1 = l_2 = \cdots = l_n and u_1 = u_2 = \cdots = u_n.

ibound = 4

: If the problem is unconstrained. (The ibound = 4 option is provided for consistency with other functions. In e04kd it produces the same effect as ibound = 1.)

bl

double array

The fixed lower bounds l_j.

bu

double array

The fixed upper bounds u_j.

x

double array

x(j)

must be set to a guess at the jth component of the position of the minimum for j=1 . . . n.

lh

integer

iw

integer array

w

double array

n

integer: default = nrow(bl)

The number n of independent variables.

iprint

integer: default = 1

The frequency with which monit is to be called.

iprint > 0

: monit is called once every iprint iterations and just before exit from e04kd.

iprint = 0

: monit is just called at the final point.

iprint < 0

: monit is not called at all.

maxcal

integer: default = 50

The maximum permitted number of evaluations of F(x), i.e., the maximum permitted number of calls of funct with iflag set to 2. It should be borne in mind that, in addition to the calls of funct which are limited directly by maxcal, there will be calls of funct (with iflag set to 1) to evaluate only first derivatives.

xtol

double: default = 0.0

The accuracy in x to which the solution is required.

delta

double: default = 0.0

The differencing interval to be used for approximating the second derivatives of F(x). Thus, for the finite difference approximations, the first derivatives of F(x) are evaluated at points which are delta apart. If ε is the machine precision, then sqrt(ε) will usually be a suitable setting for delta. If you set delta to 0.0 (or to any positive value less than ε), e04kd will automatically use sqrt(ε) as the differencing interval.

stepmx

double: default = 100000.0

An estimate of the Euclidean distance between the solution and the starting point supplied by you. (For maximum efficiency a slight overestimate is preferable.)

Details

R interface to the NAG Fortran routine E04KDF.

Value

bl

double array

The lower bounds actually used by e04kd, e.g., If ibound = 2, bl(1) = bl(2) = \cdots = bl(n) = 0.0.

bu

double array

The upper bounds actually used by e04kd, e.g., if ibound = 2, bu(1) = bu(2) = \cdots = bu(n) = 10^6.

x

double array

The final point x^(k). Thus, if ifail =0 on exit, x(j) is the jth component of the estimated position of the minimum.

hesl

double array

During the determination of a direction p_z (see the Description in Fortran library documentation), H + E is decomposed into the product LDL^T, where L is a unit lower triangular matrix and D is a diagonal matrix. (The matrices H, E, L and D are all of dimension n_z, where n_z is the number of variables free from their bounds. H consists of those rows and columns of the full estimated second derivative matrix which relate to free variables. E is chosen so that H + E is positive definite.)

hesd

double array

During the determination of a direction p_z (see the Description in Fortran library documentation), H + E is decomposed into the product LDL^T, where L is a unit lower triangular matrix and D is a diagonal matrix. (The matrices H, E, L and D are all of dimension n_z, where n_z is the number of variables free from their bounds. H consists of those rows and columns of the full estimated second derivative matrix which relate to free variables. E is chosen so that H + E is positive definite.)

istate

integer array

Information about which variables are currently on their bounds and which are free. If istate(j) is:

- equal to - 1, x_j is fixed on its upper bound;

- equal to - 2, x_j is fixed on its lower bound;

- equal to - 3, x_j is effectively a constant (i.e., l_j = u_j);

- positive, istate(j) gives the position of x_j in the sequence of free variables.

f

double

The function value at the final point given in x.

g

double array

The first derivative vector corresponding to the final point given in x. The components of g corresponding to free variables should normally be close to zero.

Author(s)

NAG

References

http://www.nag.co.uk/numeric/FL/nagdoc_fl23/pdf/E04/e04kdf.pdf

Examples

e04kd_funct = function(iflag, n, xc, fc, gc) {
    
    gc <- as.matrix(mat.or.vec(n, 1))
    fc <- 0
    
    if (iflag != 1) {
        
        fc <- (xc[1] + 10 * xc[2])^2 + 5 * (xc[3] - xc[4])^2 + 
            (xc[2] - 2 * xc[3])^4 + 10 * (xc[1] - xc[4])^4
    }
    gc[1] <- 2 * (xc[1] + 10 * xc[2]) + 40 * (xc[1] - xc[4])^3
    gc[2] <- 20 * (xc[1] + 10 * xc[2]) + 4 * (xc[2] - 2 * xc[3])^3
    gc[3] <- 10 * (xc[3] - xc[4]) - 8 * (xc[2] - 2 * xc[3])^3
    gc[4] <- 10 * (xc[4] - xc[3]) - 40 * (xc[1] - xc[4])^3
    list(IFLAG = iflag, FC = fc, GC = as.matrix(gc))
}
e04kd_monit = function(n, xc, fc, gc, istate, gpjnrm, 
    cond, posdef, niter, nf) {
    
    sprintf("\n Itn Fn evals Fn value Norm of proj gradient\n", 
        "\n")
    
    sprintf(" %3d %5d %20.4f %20.4f\n", niter, nf, fc, gpjnrm, 
        "\n")
    
    sprintf("\n J XJ GJ Status\n", "\n")
    
    for (j in c(1:n)) {
        isj <- istate[j]
        if (isj > 0) {
            
            sprintf("%2d %16.4f%20.4f %s\n", j, xc, j, gc, j, 
                " Free", "\n")
            
        }
        else if (isj == -1) {
            
        }
        else if (isj == -2) {
            
        }
        else if (isj == -3) {
            
        }
    }
    
    if (cond != 0) {
        
        if (cond > 1e+06) {
            
            sprintf("\nEstimated condition number of projected Hessian is more than 1.0e+6\n", 
                "\n")
            
        }
        else {
            
            sprintf("\nEstimated condition number of projected Hessian = %10.2f\n", 
                cond, "\n")
            
        }
        if (!posdef) {
            
            sprintf("\nProjected Hessian matrix is not positive definite\n", 
                "\n")
            
        }
    }
    list()
}

eta <- 0.5
ibound <- 0
bl <- matrix(c(1, -2, -1e+06, 1), nrow = 4, ncol = 1, 
    byrow = TRUE)


bu <- matrix(c(3, 0, 1e+06, 3), nrow = 4, ncol = 1, 
    byrow = TRUE)


x <- matrix(c(3, -1, 0, 1), nrow = 4, ncol = 1, byrow = TRUE)


lh <- 6
iw <- matrix(c(0, 0), nrow = 2, ncol = 1, byrow = TRUE)


w <- as.matrix(mat.or.vec(34, 1))
e04kd(e04kd_funct, e04kd_monit, eta, ibound, bl, bu, 
    x, lh, iw, w) 


[Package NAGFWrappers version 24.0 Index]