e04ly {NAGFWrappers} | R Documentation |
e04ly is an easy-to-use modified-Newton algorithm for finding a minimum of a function, F(x_1x_2 . . . x_n) subject to fixed upper and lower bounds on the independent variables, x_1 , x_2 , . . . , x_n when first and second derivatives of F are available. It is intended for functions which are continuous and which have continuous first and second derivatives (although it will usually work even if the derivatives have occasional discontinuities).
e04ly(ibound, funct2, hess2, bl, bu, x, n=nrow(bl) )
ibound |
integer Indicates 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 you are supplying all the l_j and u_j individually. ibound = 1: If there are no bounds on any x_j. ibound = 2: If all the bounds are of the form 0 <= x_j. ibound = 3: If l_1 = l_2 = \cdots = l_n and u_1 = u_2 = \cdots = u_n. |
funct2 |
void function You must supply this function to calculate the values of the function F(x) and its first derivatives ( \partial F)/( \partial x_j) at any point x. It should be tested separately before being used in conjunction with e04ly (see the E04 chapter introduction in the Fortran Library documentation). |
hess2 |
void function You must supply this function to evaluate the elements H_ij = (\partial^2F)/( \partial x_i \partial x_j) of the matrix of second derivatives of F(x) at any point x. It should be tested separately before being used in conjunction with e04ly (see the E04 chapter introduction in the Fortran Library documentation). |
bl |
double array The lower bounds l_j. |
bu |
double array The 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. The function checks the gradient and the Hessian matrix at the starting point, and is more likely to detect any error in your programming if the initial x(j) are nonzero and mutually distinct. |
n |
integer: default = nrow(bl) The number n of independent variables. |
R interface to the NAG Fortran routine E04LYF.
bl |
double array The lower bounds actually used by e04ly. |
bu |
double array The upper bounds actually used by e04ly. |
x |
double array The lowest point found during the calculations. Thus, if ifail =0 on exit, x(j) is the jth component of the position of the minimum. |
f |
double The value of F(x) corresponding to the final point stored in x. |
g |
double array The value of ( \partial F)/( \partial x_j) corresponding to the final point stored in x for j=1 . . . n; the value of g(j) for variables not on a bound should normally be close to zero. |
NAG
http://www.nag.co.uk/numeric/FL/nagdoc_fl23/pdf/E04/e04lyf.pdf
e04ly_funct2 = function(n, xc, fc, gc) { gc <- as.matrix(mat.or.vec(n, 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(FC = fc, GC = as.matrix(gc)) } e04ly_hess2 = function(n, xc, heslc, lh, hesdc) { heslc <- as.matrix(mat.or.vec(lh, 1)) hesdc <- as.matrix(mat.or.vec(n, 1)) hesdc[1] <- 2 + 120 * (xc[1] - xc[4])^2 hesdc[2] <- 200 + 12 * (xc[2] - 2 * xc[3])^2 hesdc[3] <- 10 + 48 * (xc[2] - 2 * xc[3])^2 hesdc[4] <- 10 + 120 * (xc[1] - xc[4])^2 heslc[1] <- 20 heslc[2] <- 0 heslc[3] <- -24 * (xc[2] - 2 * xc[3])^2 heslc[4] <- -120 * (xc[1] - xc[4])^2 heslc[5] <- 0 heslc[6] <- -10 list(HESLC = as.matrix(heslc), HESDC = as.matrix(hesdc)) } 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) e04ly(ibound, e04ly_funct2, e04ly_hess2, bl, bu, x)