e04ug {NAGFWrappers}R Documentation

e04ug: NLP problem (sparse)

Description

e04ug solves sparse nonlinear programming problems.

Usage

e04ug(confun, objfun, n, m, ncnln, nonln, njnln, iobj, a, ha, ka, bl, bu, start, names, ns, xs, istate, clamda, optlist,
      nnz = nrow(a),
      nname = nrow(names),
      leniz = (1000),
      lenz = (1000))

Arguments

confun

function

confun must calculate the vector F(x) of nonlinear constraint functions and (optionally) its Jacobian ( = ( \partial F)/( \partial x)) for a specified n_1'' ( <= n) element vector x. If there are no nonlinear constraints (i.e., ncnln = 0), confun will never be called by e04ug and confun may be the dummy function e04ugm. (e04ugm is included in the NAG Library.) If there are nonlinear constraints, the first call to confun will occur before the first call to objfun.

(MODE,F,FJAC) = confun(mode,ncnln,njnln,nnzjac,x,fjac,nstate)

objfun

function

objfun must calculate the nonlinear part of the objective function f(x) and (optionally) its gradient ( = ( \partial f)/( \partial x)) for a specified n_1' ( <= n) element vector x. If there are no nonlinear objective variables (i.e., nonln = 0), objfun will never be called by e04ug and objfun may be the dummy function e04ugn. (e04ugn is included in the NAG Library.)

(MODE,OBJF,OBJGRD) = objfun(mode,nonln,x,objgrd,nstate)

n

integer

n

, the number of variables (excluding slacks). This is the number of columns in the full Jacobian matrix A.

m

integer

m

, the number of general constraints (or slacks). This is the number of rows in A, including the free row (if any; see iobj). Note that A must contain at least one row. If your problem has no constraints, or only upper and lower bounds on the variables, then you must include a dummy ‘free’ row consisting of a single (zero) element subject to ‘infinite’ upper and lower bounds. Further details can be found under the descriptions for iobj, nnz, a, ha, ka, bl and bu.

ncnln

integer

n_N

, the number of nonlinear constraints.

nonln

integer

n_1'

, the number of nonlinear objective variables. If the objective function is nonlinear, the leading n_1' columns of A belong to the nonlinear objective variables. (See also the description for njnln.)

njnln

integer

n_1''

, the number of nonlinear Jacobian variables. If there are any nonlinear constraints, the leading n_1'' columns of A belong to the nonlinear Jacobian variables. If n_1' > 0 and n_1'' > 0, the nonlinear objective and Jacobian variables overlap. The total number of nonlinear variables is given by barn = max(n_1', n_1'').

iobj

integer

If iobj > ncnln, row iobj of A is a free row containing the nonzero elements of the linear part of the objective function.

iobj = 0

: There is no free row.

iobj = - 1

: There is a dummy ‘free’ row.

a

double array

The nonzero elements of the Jacobian matrix A, ordered by increasing column index. Since the constraint Jacobian matrix J(x'') must always appear in the top left-hand corner of A, those elements in a column associated with any nonlinear constraints must come before any elements belonging to the linear constraint matrix G and the free row (if any; see iobj).

ha

integer array

ha[i]

must contain the row index of the nonzero element stored in a[i] for i=1 . . . nnz. The row indices for a column may be supplied in any order subject to the condition that those elements in a column associated with any nonlinear constraints must appear before those elements associated with any linear constraints (including the free row, if any). Note that confun must define the Jacobian elements in the same order. If iobj = - 1, set ha[1] = 1.

ka

integer array

ka[j]

must contain the index in a of the start of the jth column for j=1 . . . n. To specify the jth column as empty, set ka[j] = ka[j+1]. Note that the first and last elements of ka must be such that ka[1] = 1 and ka[n+1] = nnz + 1. If iobj = - 1, set ka[j] = 2 for j=2 . . . n.

bl

double array

l

, the lower bounds for all the variables and general constraints, in the following order. The first n elements of bl must contain the bounds on the variables x, the next ncnln elements the bounds for the nonlinear constraints F(x) (if any) and the next (m - ncnln) elements the bounds for the linear constraints Gx and the free row (if any). To specify a nonexistent lower bound (i.e., l_j = - infinity), set bl[j] <= - bigbnd. To specify the jth constraint as an equality, set bl[j] = bu[j] = β, say, where abs(β) < bigbnd. If iobj = - 1, set bl[n+abs(iobj)] <= - bigbnd.

bu

double array

u

, the upper bounds for all the variables and general constraints, in the following order. The first n elements of bu must contain the bounds on the variables x, the next ncnln elements the bounds for the nonlinear constraints F(x) (if any) and the next (m - ncnln) elements the bounds for the linear constraints Gx and the free row (if any). To specify a nonexistent upper bound (i.e., u_j = + infinity), set bu[j] >= bigbnd. To specify the jth constraint as an equality, set bu[j] = bl[j] = β, say, where abs(β) < bigbnd. If iobj = - 1, set bu[n+abs(iobj)] >= bigbnd.

start

string

Indicates how a starting basis is to be obtained.

start='C'

: An internal Crash procedure will be used to choose an initial basis.

start='W'

: A basis is already defined in istate and ns (probably from a previous call).

names

string array

Specifies the column and row names to be used in the printed output.

ns

integer

n_S

, the number of superbasics. It need not be specified if start='C', but must retain its value from a previous call when start='W'.

xs

double array

The initial values of the variables and slacks (xs). (See the description for istate.)

istate

integer array

If start='C', the first n elements of istate and xs must specify the initial states and values, respectively, of the variables x. (The slacks s need not be initialized.) An internal Crash procedure is then used to select an initial basis matrix B. The initial basis matrix will be triangular (neglecting certain small elements in each column). It is chosen from various rows and columns of ( A -I ) . Possible values for istate[j] are as follows:

clamda

double array

If ncnln > 0, clamda[j] must contain a Lagrange multiplier estimate for the jth nonlinear constraint F_j(x) for j=n+1 . . . n+ncnln. If nothing special is known about the problem, or there is no wish to provide special information, you may set clamda[j] = 0.0. The remaining elements need not be set.

optlist

options list

Optional parameters may be listed, as shown in the following table:

Name Type Default
Central Difference Interval double Default = sqrt[3](functionprecision)
Check Frequency integer Default = 60
Crash Option integer Default = 0 or 3
Crash Tolerance double Default = 0.1
Defaults
Derivative Level integer Default = 3
Derivative Linesearch Default
Nonderivative Linesearch
Elastic Weight double Default = 1.0 or 100.0
Expand Frequency integer Default = 10000
Factorization Frequency integer Default = 50 or 100
Infeasible Exit Default
Feasible Exit
Minimize Default
Maximize
Feasible Point
Forward Difference Interval double Default = sqrt(functionprecision)
Function Precision double Default = ε^0.8
Hessian Frequency integer Default = 99999999
Hessian Full Memory Default when barn < 75
Hessian Limited Memory Default when barn >= 75
Hessian Updates integer Default = 20 or 99999999
Infinite Bound Size double Default = 10^20
Iteration Limit integer Default = 10000
Linesearch Tolerance double Default = 0.9
List Default for e04ug = list
Nolist Default for e04ug = nolist
LU Density Tolerance double Default = 0.6
LU Singularity Tolerance double Default = ε^0.67
LU Factor Tolerance double Default = 5.0 or 100.0
LU Update Tolerance double Default = 5.0 or 10.0
Major Feasibility Tolerance double Default = sqrt(ε)
Major Iteration Limit integer Default = 1000
Major Optimality Tolerance double Default = sqrt(ε)
Optimality Tolerance double
Major Print Level integer = 0
Print Level
Major Step Limit double Default = 2.0
Minor Feasibility Tolerance double Default = sqrt(ε)
Feasibility Tolerance double
Minor Iteration Limit integer Default = 500
Minor Optimality Tolerance double Default = sqrt(ε)
Minor Print Level integer Default = 0
Monitoring File integer Default = - 1
Partial Price integer Default = 1 or 10
Pivot Tolerance double Default = ε^0.67
Scale Option integer Default = 1 or 2
Scale Tolerance double Default = 0.9
Start Objective Check At Column integer Default = 1
Stop Objective Check At Column integer Default = n_1'
Start Constraint Check At Column integer Default = 1
Stop Constraint Check At Column integer Default = n_1''
Superbasics Limit integer Default = min(500, barn + 1)
Unbounded Objective double Default = 10^15
Unbounded Step Size double Default = max(bigbnd , 10^20)
Verify Level integer Default = 0
Violation Limit double Default = 10.0
nnz

integer: default = nrow(a)

The number of nonzero elements in A (including the Jacobian for any nonlinear constraints). If iobj = - 1, set nnz = 1.

nname

integer: default = nrow(names)

The number of column (i.e., variable) and row (i.e., constraint) names supplied in names.

nname = 1

: There are no names. Default names will be used in the printed output.

nname = n + m

: All names must be supplied.

leniz

integer: default = (max(500,(n+m)))

integer: default = (max(500,(n+m)))

lenz

integer: default = (500)

integer: default = (500)

Details

R interface to the NAG Fortran routine E04UGF.

Value

A

double array

Elements in the nonlinear part corresponding to nonlinear Jacobian variables are overwritten.

NS

integer

The final number of superbasics.

XS

double array

The final values of the variables and slacks (xs).

ISTATE

integer array

The final states of the variables and slacks (xs). The significance of each possible value of istate[j] is as follows:

CLAMDA

double array

A set of Lagrange multipliers for the bounds on the variables (reduced costs) and the general constraints (shadow costs). More precisely, the first n elements contain the multipliers for the bounds on the variables, the next ncnln elements contain the multipliers for the nonlinear constraints F(x) (if any) and the next (m - ncnln) elements contain the multipliers for the linear constraints Gx and the free row (if any).

MINIZ

integer

The minimum value of leniz required to start solving the problem. If ifail =12, e04ug may be called again with leniz suitably larger than miniz. (The bigger the better, since it is not certain how much workspace the basis factors need.)

MINZ

integer

The minimum value of lenz required to start solving the problem. If ifail =13, e04ug may be called again with lenz suitably larger than minz. (The bigger the better, since it is not certain how much workspace the basis factors need.)

NINF

integer

The number of constraints that lie outside their bounds by more than the value of the optional argument minorfeasibilitytolerance.

SINF

double

The sum of the infeasibilities of constraints that lie outside their bounds by more than the value of the optional argument minorfeasibilitytolerance.

OBJ

double

The value of the objective function.

IFAIL

integer

ifail =0

unless the function detects an error or a warning has been flagged (see the Errors section in Fortran library documentation).

Author(s)

NAG

References

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

Examples


optlist <- list()

ifail <- 0
confun = function(mode, ncnln, njnln, nnzjac, x, fjac, 
    nstate) {
    
    f <- as.matrix(mat.or.vec(ncnln, 1))
    
    if (mode == 0 || mode == 2) {
        
        f[1] <- 1000 %*% sin(-x[1] - 0.25) + 1000 %*% sin(-x[2] - 
            0.25)
        
        f[2] <- 1000 %*% sin(x[1] - 0.25) + 1000 %*% sin(x[1] - 
            x[2] - 0.25)
        
        f[3] <- 1000 %*% sin(x[2] - x[1] - 0.25) + 1000 %*% sin(x[2] - 
            0.25)
        
    }
    
    if (mode == 1 || mode == 2) {
        
        fjac[1] <- -1000 %*% cos(-x[1] - 0.25)
        
        fjac[2] <- 1000 %*% cos(x[1] - 0.25) + 1000 %*% cos(x[1] - 
            x[2] - 0.25)
        
        fjac[3] <- -1000 %*% cos(x[2] - x[1] - 0.25)
        
        fjac[4] <- -1000 %*% cos(-x[2] - 0.25)
        
        fjac[5] <- -1000 %*% cos(x[1] - x[2] - 0.25)
        
        fjac[6] <- 1000 %*% cos(x[2] - x[1] - 0.25) + 1000 %*% 
            cos(x[2] - 0.25)
        
    }
    list(MODE = as.integer(mode), F = as.matrix(f), FJAC = as.matrix(fjac))
}
objfun = function(mode, nonln, x, objgrd, nstate) {
    
    
    if (mode == 0 || mode == 2) {
        
        objf <- 1e-06 %*% x[3]^3 + 2e-06 %*% x[4]^3/3
        
    }
    
    if (mode == 1 || mode == 2) {
        
        objgrd[1] <- 0
        
        objgrd[2] <- 0
        
        objgrd[3] <- 3e-06 %*% x[3]^2
        
        objgrd[4] <- 2e-06 %*% x[4]^2
        
    }
    list(MODE = as.integer(mode), OBJF = objf, OBJGRD = as.matrix(objgrd))
}

n <- 4

m <- 6

ncnln <- 3

nonln <- 4

njnln <- 2

iobj <- 6

a <- matrix(c(1e+25, 1e+25, 1e+25, 1, -1, 1e+25, 1e+25, 
    1e+25, -1, 1, 3, -1, -1, 2), nrow = 14, ncol = 1, byrow = TRUE)



ha <- matrix(c(1, 2, 3, 5, 4, 1, 2, 3, 5, 4, 6, 1, 
    2, 6), nrow = 14, ncol = 1, byrow = TRUE)



ka <- matrix(c(1, 6, 11, 13, 15), nrow = 5, ncol = 1, 
    byrow = TRUE)



bl <- matrix(c(-0.55, -0.55, 0, 0, -894.8, -894.8, 
    -1294.8, -0.55, -0.55, -1e+25), nrow = 10, ncol = 1, byrow = TRUE)



bu <- matrix(c(0.55, 0.55, 1200, 1200, -894.8, -894.8, 
    -1294.8, 1e+25, 1e+25, 1e+25), nrow = 10, ncol = 1, byrow = TRUE)



start <- "C"

names <- matrix(c("Varble 1", "Varble 2", "Varble 3", 
    "Varble 4", "NlnCon 1", "NlnCon 2", "NlnCon 3", "LinCon 1", 
    "LinCon 2", "Free Row"), nrow = 10, byrow = TRUE)



ns <- 0

xs <- matrix(c(0, 0, 0, 0, 0, 0, 0, 0, 0, 0), nrow = 10, 
    ncol = 1, byrow = TRUE)



istate <- as.matrix(mat.or.vec(10, 1))

clamda <- matrix(c(0, 0, 0, 0, 0, 0, 0, 0, 0, 0), 
    nrow = 10, ncol = 1, byrow = TRUE)



leniz <- 1000

lenz <- 1000

e04ug(confun, objfun, n, m, ncnln, nonln, njnln, 
    iobj, a, ha, ka, bl, bu, start, names, ns, xs, istate, clamda, 
    optlist)

[Package NAGFWrappers version 24.0 Index]