e04jy {NAGFWrappers} | R Documentation |
e04jy is an easy-to-use quasi-Newton algorithm for finding a minimum of a function F(x_1x_2 . . . x_n), subject to fixed upper and lower bounds of the independent variables x_1 , x_2 , . . . , x_n, using function values only.
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).
e04jy(ibound, funct1, bl, bu, x, n=nrow(bl), liw=n+2, lw=max(n*(n-1)/2+12*n,13) )
ibound |
integer Indicates whether the facility for dealing with bounds of special forms is to be used. |
funct1 |
void function You must supply funct1 to calculate the value of the function F(x) at any point x. It should be tested separately before being used with e04jy (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 an estimate of the jth component of the position of the minimum for j=1 . . . n. |
n |
integer: default = nrow(bl) The number n of independent variables. |
liw |
integer: default = n+2 |
lw |
integer: default = max(n*(n-1)/2+12*n,13) |
R interface to the NAG Fortran routine E04JYF.
bl |
double array The lower bounds actually used by e04jy. |
bu |
double array The upper bounds actually used by e04jy. |
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. |
iw |
integer array If ifail =0, ifail =3, ifail =5, the first n elements of iw contain information about which variables are currently on their bounds and which are free. Specifically, if x_i is: -: fixed on its upper bound, iw(i) is - 1; -: fixed on its lower bound, iw(i) is - 2; -: effectively a constant (i.e., l_j = u_j), iw(i) is - 3; -: free, iw(i) gives its position in the sequence of free variables. |
w |
double array If ifail =0, ifail =3, ifail =5, w(i) contains a finite difference approximation to the ith element of the projected gradient vector g_z for i=1 . . . n. In addition, w(n+1) contains an estimate of the condition number of the projected Hessian matrix (i.e., k). The rest of the array is used as workspace. |
NAG
http://www.nag.co.uk/numeric/FL/nagdoc_fl23/pdf/E04/e04jyf.pdf
e04jy_funct1 = function(n, xc, fc) { 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 list(FC = fc) } 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) e04jy(ibound, e04jy_funct1, bl, bu, x)