hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_opt_qpconvex2_sparse_solve (e04nq)

Purpose

nag_opt_qpconvex2_sparse_solve (e04nq) solves sparse linear programming or convex quadratic programming problems. The initialization function nag_opt_qpconvex2_sparse_init (e04np) must have been called before calling nag_opt_qpconvex2_sparse_solve (e04nq).

Syntax

[hs, x, pi, rc, ns, ninf, sinf, obj, user, cw, iw, rw, ifail] = e04nq(start, qphx, m, n, lenc, ncolh, iobj, objadd, prob, acol, inda, loca, bl, bu, c, names, helast, hs, x, ns, cw, iw, rw, 'ne', ne, 'nname', nname, 'user', user)
[hs, x, pi, rc, ns, ninf, sinf, obj, user, cw, iw, rw, ifail] = nag_opt_qpconvex2_sparse_solve(start, qphx, m, n, lenc, ncolh, iobj, objadd, prob, acol, inda, loca, bl, bu, c, names, helast, hs, x, ns, cw, iw, rw, 'ne', ne, 'nname', nname, 'user', user)
Before calling nag_opt_qpconvex2_sparse_solve (e04nq) or one of the option setting functions nag_opt_qpconvex2_sparse_option_string (e04ns), nag_opt_qpconvex2_sparse_option_integer_set (e04nt) or nag_opt_qpconvex2_sparse_option_double_set (e04nu), nag_opt_qpconvex2_sparse_init (e04np) must be called.
Note: the interface to this routine has changed since earlier releases of the toolbox:
Mark 22: lencw, leniw, lenrw have been removed from the interface
.

Description

nag_opt_qpconvex2_sparse_solve (e04nq) is designed to solve large-scale linear or quadratic programming problems of the form:
minimize f(x)  subject to ​l
(x)
Ax
u,
xRn
minimize xRn f(x)   subject to ​ l x Ax u,
(1)
where xx is an nn-vector of variables, ll and uu are constant lower and upper bounds, AA is an mm by nn sparse matrix and f(x)f(x) is a linear or quadratic objective function that may be specified in a variety of ways, depending upon the particular problem being solved. The optional parameter Maximize may be used to specify a problem in which f(x)f(x) is maximized instead of minimized.
Upper and lower bounds are specified for all variables and constraints. This form allows full generality in specifying various types of constraint. In particular, the jjth constraint may be defined as an equality by setting lj = ujlj=uj. If certain bounds are not present, the associated elements of ll or uu may be set to special values that are treated as - or + +.
The possible forms for the function f(x)f(x) are summarised in Table 1. The most general form for f(x)f(x) is
n nn
f(x) = q + cTx + (1/2)xTHx = q + cjxj + (1/2)xiHijxj
j = 1 i = 1j = 1
f(x) = q + cTx + 12 xTHx = q + j=1 n cj xj + 12 i=1 n j=1 n xi Hij xj
where qq is a constant, cc is a constant nn-vector and HH is a constant symmetric nn by nn matrix with elements {Hij} { Hij } . In this form, ff is a quadratic function of xx and (1) is known as a quadratic program (QP). nag_opt_qpconvex2_sparse_solve (e04nq) is suitable for all convex quadratic programs. The defining feature of a convex QP is that the matrix HH must be positive semidefinite, i.e., it must satisfy xT H x 0 xT H x 0  for all xx. If not, f(x) f(x)  is nonconvex and nag_opt_qpconvex2_sparse_solve (e04nq) will terminate with the error indicator ifail = 11ifail=11. If f(x) f(x)  is nonconvex it may be more appropriate to call nag_opt_nlp2_sparse_solve (e04vh) instead.
Problem type Objective function f(x)f(x) Hessian matrix HH
FP Not applicable q = c = H = 0q=c=H=0
LP q + cTxq+cTx H = 0H=0
QP q + cTx + (1/2)xTHxq+cTx+12xTHx Symmetric positive semidefinite
Table 1
Choices for the objective function f(x)f(x)
If H = 0 H=0 , then f(x) = q + cTx f(x) = q+cTx  and the problem is known as a linear program (LP). In this case, rather than defining an HH with zero elements, you can define HH to have no columns by setting ncolh = 0ncolh=0 (see Section [Parameters]).
If H = 0 H=0 , q = 0 q=0 , and c = 0 c=0 , there is no objective function and the problem is a feasible point problem (FP), which is equivalent to finding a point that satisfies the constraints on xx. In the situation where no feasible point exists, several options are available for finding a point that minimizes the constraint violations (see the description of the optional parameter Elastic Mode).
nag_opt_qpconvex2_sparse_solve (e04nq) is suitable for large LPs and QPs in which the matrix AA is sparse, i.e., when the number of zero elements is sufficiently large that it is worthwhile using algorithms which avoid computations and storage involving zero elements. The matrix AA is input to nag_opt_qpconvex2_sparse_solve (e04nq) by means of the three array parameters acol, inda and loca. This allows you to specify the pattern of nonzero elements in AA.
nag_opt_qpconvex2_sparse_solve (e04nq) exploits structure in HH by requiring HH to be defined implicitly in a function that computes the product HxHx for any given vector xx. In many cases, the product HxHx can be computed very efficiently for any given xx, e.g., HH may be a sparse matrix, or a sum of matrices of rank-one.
For problems in which AA can be treated as a dense matrix, it is usually more efficient to use nag_opt_lp_solve (e04mf), nag_opt_lsq_lincon_solve (e04nc) or nag_opt_qp_dense_solve (e04nf).
There is considerable flexibility allowed in the definition of f(x) f(x)  in Table 1. The vector cc defining the linear term cTxcTx can be input in three ways: as a sparse row of AA; as an explicit dense vector cc; or as both a sparse row and an explicit vector (in which case, cTxcTx will be the sum of two linear terms). When stored in AA, cc is the iobjth row of AA, which is known as the objective row. The objective row must always be a free row of AA in the sense that its lower and upper bounds must be -  and + + . Storing cc as part of AA is recommended if cc is a sparse vector. Storing cc as an explicit vector is recommended for a sequence of problems, each with a different objective (see parameters c and lenc).
The upper and lower bounds on the mm elements of AxAx are said to define the general constraints of the problem. Internally, nag_opt_qpconvex2_sparse_solve (e04nq) converts the general constraints to equalities by introducing a set of slack variables ss, where s = (s1,s2,,sm)Ts=(s1,s2,,sm)T. For example, the linear constraint 52x1 + 3x2 + 52x1+3x2+ is replaced by 2x1 + 3x2s1 = 02x1+3x2-s1=0, together with the bounded slack 5s1 + 5s1+. The problem defined by (1) can therefore be re-written in the following equivalent form:
minimize f(x) subject to ​Axs = 0,  l
(x)
s
u.
x Rn , s Rm
minimize x Rn , s Rm f(x)  subject to ​ Ax-s=0 ,   l x s u.
Since the slack variables ss are subject to the same upper and lower bounds as the elements of AxAx, the bounds on xx and AxAx can simply be thought of as bounds on the combined vector (x,s)(x,s). (In order to indicate their special role in QP problems, the original variables xx are sometimes known as ‘column variables’, and the slack variables ss are known as ‘row variables’.)
Each LP or QP problem is solved using a two-phase iterative procedure (in which the general constraints are satisfied throughout): a feasibility phase (Phase 1), in which the sum of infeasibilities with respect to the bounds on xx and ss is minimized to find a feasible point that satisfies all constraints within a specified feasibility tolerance; and an optimality phase (Phase 2), in which f(x)f(x) is minimized (or maximized) by constructing a sequence of iterates that lies within the feasible region.
Phase 1 involves solving a linear program of the form
Phase 1  
  minimize x, s, v, w   j = 1n + m ( vj + wj ) minimize x, s, v, w j=1 n+m ( vj + wj )  
  subject to ​ Axs = 0 ,  l
(x)
s
v + w u ,  v 0 ,  w 0
subject to ​ Ax-s=0 ,  l x s -v + w u ,  v 0 ,  w 0  
which is equivalent to minimizing the sum of the constraint violations. If the constraints are feasible (i.e., at least one feasible point exists), eventually a point will be found at which both vv and ww are zero. Then the associated value of (x,s) (x,s)  satisfies the original constraints and is used as the starting point for the Phase 2 iterations for minimizing f(x) f(x) .
If the constraints are infeasible (i.e., v0 v0  or w0 w0  at the end of Phase 1), no solution exists for (1) and you have the option of either terminating or continuing in so-called elastic mode (see the discussion of the optional parameter Elastic Mode). In elastic mode, a ‘relaxed’ or ‘perturbed’ problem is solved in which f(x) f(x)  is minimized while allowing some of the bounds to become ‘elastic’, i.e., to change from their specified values. Variables subject to elastic bounds are known as elastic variables. An elastic variable is free to violate one or both of its original upper or lower bounds. You are able to assign which bounds will become elastic if elastic mode is ever started (see the parameter helast in Section [Parameters]).
To make the relaxed problem meaningful, nag_opt_qpconvex2_sparse_solve (e04nq) minimizes f(x) f(x)  while (in some sense) finding the ‘smallest’ violation of the elastic variables. In the situation where all the variables are elastic, the relaxed problem has the form
Phase 2 (γγ)  
  minimize x, s, v, w   f(x) + γ j = 1n + m ( vj + wj ) minimize x, s, v, w f(x) + γ j=1 n+m ( vj + wj )  
  subject to ​ Axs = 0 ,  l
(x)
s
v + w u ,  v 0 ,  w 0
subject to ​ Ax-s=0 ,  l x s -v + w u ,  v 0 ,  w 0 ,
where γγ is a non-negative parameter known as the elastic weight (see the description of the optional parameter Elastic Weight), and f(x) + γ j  (vj + wj) f(x) + γ j ( vj + wj )  is called the composite objective. In the more general situation where only a subset of the bounds are elastic, the vv's and ww's for the non-elastic bounds are fixed at zero.
The elastic weight can be chosen to make the composite objective behave like the original objective f(x) f(x) , the sum of infeasibilities, or anything in-between. If γ = 0 γ=0 , nag_opt_qpconvex2_sparse_solve (e04nq) will attempt to minimize ff subject to the (true) upper and lower bounds on the non-elastic variables (and declare the problem infeasible if the non-elastic variables cannot be made feasible).
At the other extreme, choosing γγ sufficiently large will have the effect of minimizing the sum of the violations of the elastic variables subject to the original constraints on the non-elastic variables. Choosing a large value of the elastic weight is useful for defining a ‘least-infeasible’ point for an infeasible problem.
In Phase 1 and elastic mode, all calculations involving vv and ww are done implicitly in the sense that an elastic variable xj xj  is allowed to violate its lower bound (say) and an explicit value of vv can be recovered as vj = lj xj vj = lj - xj .
A constraint is said to be active or binding at xx if the associated element of either xx or AxAx is equal to one of its upper or lower bounds. Since an active constraint in AxAx has its associated slack variable at a bound, the status of both simple and general upper and lower bounds can be conveniently described in terms of the status of the variables (x,s)(x,s). A variable is said to be nonbasic if it is temporarily fixed at its upper or lower bound. It follows that regarding a general constraint as being active is equivalent to thinking of its associated slack as being nonbasic.
At each iteration of an active-set method, the constraints Axs = 0Ax-s=0 are (conceptually) partitioned into the form
BxB + SxS + NxN = 0,
BxB+SxS+NxN=0,
where xNxN consists of the nonbasic elements of (x,s)(x,s) and the basis matrix BB is square and nonsingular. The elements of xBxB and xSxS are called the basic and superbasic variables respectively; with xNxN they are a permutation of the elements of xx and ss. At a QP solution, the basic and superbasic variables will lie somewhere between their upper or lower bounds, while the nonbasic variables will be equal to one of their bounds. At each iteration, xSxS is regarded as a set of independent variables that are free to move in any desired direction, namely one that will improve the value of the objective function (or sum of infeasibilities). The basic variables are then adjusted in order to ensure that (x,s)(x,s) continues to satisfy Axs = 0Ax-s=0. The number of superbasic variables (nSnS say) therefore indicates the number of degrees of freedom remaining after the constraints have been satisfied. In broad terms, nSnS is a measure of how nonlinear the problem is. In particular, nSnS will always be zero for FP and LP problems.
If it appears that no improvement can be made with the current definition of BB, SS and NN, a nonbasic variable is selected to be added to SS, and the process is repeated with the value of nSnS increased by one. At all stages, if a basic or superbasic variable encounters one of its bounds, the variable is made nonbasic and the value of nSnS is decreased by one.
Associated with each of the mm equality constraints Axs = 0Ax-s=0 is a dual variable πiπi. Similarly, each variable in (x,s)(x,s) has an associated reduced gradient djdj (also known as a reduced cost). The reduced gradients for the variables xx are the quantities gATπg-ATπ, where gg is the gradient of the QP objective function, and the reduced gradients for the slack variables ss are the dual variables ππ. The QP subproblem is optimal if dj0dj0 for all nonbasic variables at their lower bounds, dj0dj0 for all nonbasic variables at their upper bounds and dj = 0dj=0 for all superbasic variables. In practice, an approximate QP solution is found by slightly relaxing these conditions on djdj (see the description of the optional parameter Optimality Tolerance).
The process of computing and comparing reduced gradients is known as pricing (a term first introduced in the context of the simplex method for linear programming). To ‘price’ a nonbasic variable xjxj means that the reduced gradient djdj associated with the relevant active upper or lower bound on xjxj is computed via the formula dj = gjajTπdj=gj-ajTπ, where ajaj is the jjth column of
(AI)
A -I . (The variable selected by such a process and the corresponding value of djdj (i.e., its reduced gradient) are the quantities +SBS and dj in the monitoring file output; see Section [Printed output].) If AA has significantly more columns than rows (i.e., nmnm), pricing can be computationally expensive. In this case, a strategy known as partial pricing can be used to compute and compare only a subset of the djdjs.
nag_opt_qpconvex2_sparse_solve (e04nq) is based on SQOPT, which is part of the SNOPT package described in Gill et al. (2005a). It uses stable numerical methods throughout and includes a reliable basis package (for maintaining sparse LULU factors of the basis matrix BB), a practical anti-degeneracy procedure, efficient handling of linear constraints and bounds on the variables (by an active-set strategy), as well as automatic scaling of the constraints. Further details can be found in Section [Algorithmic Details].

References

Fourer R (1982) Solving staircase linear programs by the simplex method Math. Programming 23 274–313
Gill P E and Murray W (1978) Numerically stable methods for quadratic programming Math. Programming 14 349–372
Gill P E, Murray W and Saunders M A (1995) User's guide for QPOPT 1.0: a Fortran package for quadratic programming Report SOL 95-4 Department of Operations Research, Stanford University
Gill P E, Murray W and Saunders M A (2005a) Users' guide for SQOPT 7: a Fortran package for large-scale linear and quadratic programming Report NA 05-1 Department of Mathematics, University of California, San Diego ftp://www.cam.ucsd.edu/pub/peg/reports/sqdoc7.pdf
Gill P E, Murray W and Saunders M A (2005b) Users' guide for SNOPT 7.1: a Fortran package for large-scale linear nonlinear programming Report NA 05-2 Department of Mathematics, University of California, San Diego ftp://www.cam.ucsd.edu/pub/peg/reports/sndoc7.pdf
Gill P E, Murray W, Saunders M A and Wright M H (1987) Maintaining LU factors of a general sparse matrix Linear Algebra and its Applics. 88/89 239–270
Gill P E, Murray W, Saunders M A and Wright M H (1989) A practical anti-cycling procedure for linearly constrained optimization Math. Programming 45 437–474
Gill P E, Murray W, Saunders M A and Wright M H (1991) Inertia-controlling methods for general quadratic programming SIAM Rev. 33 1–36
Hall J A J and McKinnon K I M (1996) The simplest examples where the simplex method cycles and conditions where EXPAND fails to prevent cycling Report MS 96–100 Department of Mathematics and Statistics, University of Edinburgh

Parameters

The first nn entries of the parameters bl, bu, hs and x refer to the variables xx. The last mm entries refer to the slacks ss.

Compulsory Input Parameters

1:     start – string (length ≥ 1)
Indicates how a starting basis (and certain other items) will be obtained.
start = 'C'start='C'
Requests that an internal Crash procedure be used to choose an initial basis, unless a Basis file is provided via optional parameters Old Basis File, Insert File or Load File.
start = 'B'start='B'
Is the same as start = 'C'start='C' but is more meaningful when a Basis file is given.
start = 'W'start='W'
Means that a basis is already defined in hs and a start point is already defined in x (probably from an earlier call).
Constraint: start = 'B'start='B', 'C''C' or 'W''W'.
2:     qphx – function handle or string containing name of m-file
For QP problems, you must supply a version of qphx to compute the matrix product HxHx for a given vector xx. If HH has rows and columns of zeros, it is most efficient to order xx so that the nonlinear variables appear first. For example, if x = (y,z)Tx=(y,z)T and only yy enters the objective quadratically then
H x =
(H10)
0 0
(y)
z
=
(H1y)
0
.
H x= H1 0 0 0 y z = H1y 0 .
(2)
In this case, ncolh should be the dimension of yy, and qphx should compute H1yH1y. For FP and LP problems, qphx will never be called by nag_opt_qpconvex2_sparse_solve (e04nq) and hence qphx may be the string 'e04nsh'.
[hx, user] = qphx(ncolh, x, nstate, user)

Input Parameters

1:     ncolh – int64int32nag_int scalar
This is the same parameter ncolh as supplied to nag_opt_qpconvex2_sparse_solve (e04nq).
2:     x(ncolh) – double array
The first ncolh elements of the vector xx.
3:     nstate – int64int32nag_int scalar
Allows you to save computation time if certain data must be read or calculated only once. To preserve this data for a subsequent calculation place it in one of user, user or user .
nstate = 1nstate=1
nag_opt_qpconvex2_sparse_solve (e04nq) is calling qphx for the first time.
nstate = 0nstate=0
There is nothing special about the current call of qphx.
nstate2nstate2
nag_opt_qpconvex2_sparse_solve (e04nq) is calling qphx for the last time. This parameter setting allows you to perform some additional computation on the final solution.
nstate = 2nstate=2
The current xx is optimal.
nstate = 3nstate=3
The problem appears to be infeasible.
nstate = 4nstate=4
The problem appears to be unbounded.
nstate = 5nstate=5
The iterations limit was reached.
4:     user – Any MATLAB object
qphx is called from nag_opt_qpconvex2_sparse_solve (e04nq) with the object supplied to nag_opt_qpconvex2_sparse_solve (e04nq).

Output Parameters

1:     hx(ncolh) – double array
The product HxHx. If ncolh is less than the input parameter n, HxHx is really the product H1yH1y in (2).
2:     user – Any MATLAB object
3:     m – int64int32nag_int scalar
mm, the number of general linear constraints (or slacks). This is the number of rows in the linear constraint matrix AA, including the free row (if any; see iobj). Note that AA must have at least one row. If your problem has no constraints, or only upper or lower bounds on the variables, then you must include a dummy row with sufficiently wide upper and lower bounds (see also acol, inda and loca).
Constraint: m1m1.
4:     n – int64int32nag_int scalar
nn, the number of variables (excluding slacks). This is the number of columns in the linear constraint matrix AA.
Constraint: n1n1.
5:     lenc – int64int32nag_int scalar
The number of elements in the constant objective vector cc.
If lenc > 0lenc>0, the first lenc elements of xx belong to variables corresponding to the constant objective term cc.
Constraint: 0lencn0lencn.
6:     ncolh – int64int32nag_int scalar
nHnH, the number of leading nonzero columns of the Hessian matrix HH. For FP and LP problems, ncolh must be set to zero.
The first ncolh elements of xx belong to variables corresponding to the nonzero block of the QP Hessian.
Constraint: 0ncolhn0ncolhn.
7:     iobj – int64int32nag_int scalar
If iobj > 0iobj>0, row iobj of AA is a free row containing the nonzero elements of the vector cc appearing in the linear objective term cTxcTx.
If iobj = 0iobj=0, there is no free row, and the linear objective vector should be supplied in array c.
Constraint: 0iobjm0iobjm.
8:     objadd – double scalar
The constant qq, to be added to the objective for printing purposes. Typically objadd = 0.0e0objadd=0.0e0.
9:     prob – string (length at least 8) (length ≥ 8)
The name for the problem. It is used in the printed solution and in some functions that output Basis files. A blank name may be used.
10:   acol(ne) – double array
ne, the dimension of the array, must satisfy the constraint 1nen × m1nen×m.
The nonzero elements of AA, ordered by increasing column index. Note that all elements must be assigned a value in the calling program.
11:   inda(ne) – int64int32nag_int array
ne, the dimension of the array, must satisfy the constraint 1nen × m1nen×m.
inda(i)indai must contain the row index of the nonzero element stored in acol(i)acoli, for i = 1,2,,nei=1,2,,ne. Thus a pair of values (acol(i),inda(i))(acoli,indai) contains a matrix element and its corresponding row index.
Note that the row indices for a column may be supplied in any order.
Constraint: 1inda(i)m1indaim, for i = 1,2,,nei=1,2,,ne.
12:   loca(n + 1n+1) – int64int32nag_int array
loca(j)locaj must contain the index in acol and inda of the start of the jjth column, for j = 1,2,,nj=1,2,,n. Thus for j = 1 : nj=1:n, the entries of column jj are held in acol(k : l)acolk:l and their corresponding row indices are in inda(k : l)indak:l, where k = loca(j)k=locaj and l = loca(j + 1)1 l=locaj+1-1 . To specify the jjth column as empty, set loca(j) = loca(j + 1)locaj=locaj+1. Note that the first and last elements of loca must be loca(1) = 1loca1=1 and loca(n + 1) = ne + 1locan+1=ne+1. If your problem has no constraints, or just bounds on the variables, you may include a dummy ‘free’ row with a single (zero) element by setting ne = 1ne=1, acol(1) = 0.0acol1=0.0, inda(1) = 1inda1=1, loca(1) = 1loca1=1, and loca(j) = 2locaj=2, for j = 2 : n + 1j=2:n+1. This row is made ‘free’ by setting its bounds to be bl(n + 1) = bigbndbln+1=-bigbnd and bu(n + 1) = bigbndbun+1=bigbnd, where bigbndbigbnd is the value of the optional parameter Infinite Bound Size.
Constraints:
  • loca(1) = 1loca1=1;
  • loca(j)1locaj1, for j = 2,3,,nj=2,3,,n;
  • loca(n + 1) = ne + 1locan+1=ne+1;
  • 0loca(j + 1)loca(j)m0locaj+1-locajm, for j = 1,2,,nj=1,2,,n.
13:   bl(n + mn+m) – double array
ll, 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 xx, and the next m elements the bounds for the general linear constraints AxAx (which, equivalently, are the bounds for the slacks, ss) and the free row (if any). To fix the jjth variable, set bl(j) = bu(j) = βblj=buj=β, say, where |β| < bigbnd |β|<bigbnd . To specify a nonexistent lower bound (i.e., lj = lj=-), set bl(j)bigbndblj-bigbnd. Here, bigbndbigbnd is the value of the optional parameter Infinite Bound Size. To specify the jjth constraint as an equality, set bl(n + j) = bu(n + j) = βbln+j=bun+j=β, say, where |β| < bigbnd|β|<bigbnd. Note that the lower bound corresponding to the free row must be set to - and stored in bl(n + iobj)bln+iobj.
Constraint: if iobj > 0iobj>0, bl(n + iobj)bigbndbln+iobj-bigbnd
(See also the description for bu.)
14:   bu(n + mn+m) – double array
uu, 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 xx, and the next m elements the bounds for the general linear constraints AxAx (which, equivalently, are the bounds for the slacks, ss) and the free row (if any). To specify a nonexistent upper bound (i.e., uj = + uj=+), set bu(j)bigbndbujbigbnd. Note that the upper bound corresponding to the free row must be set to + + and stored in bu(n + iobj)bun+iobj.
Constraints:
  • if iobj > 0iobj>0, bu(n + iobj)bigbndbun+iobjbigbnd;
  • otherwise bl(i)bu(i)blibui.
15:   c(max (1,lenc)max(1,lenc)) – double array
Note: the dimension of the array c must be at least max (1,lenc)max(1,lenc) if iobj0iobj0, and at least 11 otherwise.
Contains the explicit objective vector cc (if any). If the problem is of type FP, or if lenc = 0lenc=0, then c is not referenced. (In that case, c may be dimensioned (1), or it could be any convenient array.)
16:   names(nname) – cell array of strings
nname, the dimension of the array, must satisfy the constraint nname = 1nname=1 or n + mn+m.
The optional column and row names, respectively.
If nname = 1nname=1, names is not referenced and the printed output will use default names for the columns and rows.
If nname = n + mnname=n+m, the first n elements must contain the names for the columns and the next m elements must contain the names for the rows. Note that the name for the free row (if any) must be stored in names(n + iobj)namesn+iobj.
17:   helast(n + mn+m) – int64int32nag_int array
Defines which variables are to be treated as being elastic in elastic mode. The allowed values of helast are:
helast(j)helastj Status in elastic mode
00 Variable jj is non-elastic and cannot be infeasible
11 Variable jj can violate its lower bound
22 Variable jj can violate its upper bound
33 Variable jj can violate either its lower or upper bound
helast need not be assigned if optional parameter Elastic Mode = 0Elastic Mode=0.
Constraint: if Elastic Mode0Elastic Mode0, helast(j) = 0,1,2,3helastj=0,1,2,3, for j = 1,2,,n + mj=1,2,,n+m.
18:   hs(n + mn+m) – int64int32nag_int array
If start = 'C'start='C' or 'B''B', and a Basis file of some sort is to be input (see the description of the optional parameters Old Basis File, Insert File or Load File), then hs and x need not be set at all.
If start = 'C'start='C' or 'B''B' and there is no Basis file, the first n elements of hs and x must specify the initial states and values, respectively, of the variables xx. (The slacks ss need not be initialized.) An internal Crash procedure is then used to select an initial basis matrix BB. The initial basis matrix will be triangular (neglecting certain small elements in each column). It is chosen from various rows and columns of
(AI)
A -I . Possible values for hs(j)hsj are as follows:
hs(j)hsj State of x(j)xj during Crash procedure
00 or 11 Eligible for the basis
22 Ignored
33 Eligible for the basis (given preference over 00 or 11)
44 or 55 Ignored
If nothing special is known about the problem, or there is no wish to provide special information, you may set hs(j) = 0hsj=0 and x(j) = 0.0xj=0.0, for j = 1,2,,nj=1,2,,n. All variables will then be eligible for the initial basis. Less trivially, to say that the jjth variable will probably be equal to one of its bounds, set hs(j) = 4hsj=4 and x(j) = bl(j)xj=blj or hs(j) = 5hsj=5 and x(j) = bu(j)xj=buj as appropriate.
Following the Crash procedure, variables for which hs(j) = 2hsj=2 are made superbasic. Other variables not selected for the basis are then made nonbasic at the value x(j)xj if bl(j)x(j)bu(j)bljxjbuj, or at the value bl(j)blj or bu(j)buj closest to x(j)xj.
If start = 'W'start='W', hs and x must specify the initial states and values, respectively, of the variables and slacks (x,s)(x,s). If nag_opt_qpconvex2_sparse_solve (e04nq) has been called previously with the same values of n and m, hs already contains satisfactory information.
Constraints:
  • if start = 'C'start='C' or 'B''B', 0hs(j)50hsj5, for j = 1,2,,nj=1,2,,n;
  • if start = 'W'start='W', 0hs(j)30hsj3, for j = 1,2,,n + mj=1,2,,n+m.
19:   x(n + mn+m) – double array
The initial values of the variables x x , and, if start = 'W'start='W', the slacks s s , i.e., (x,s)(x,s). (See the description for parameter hs.)
20:   ns – int64int32nag_int scalar
nSnS, the number of superbasics. For QP problems, ns need not be specified if start = 'C'start='C', but must retain its value from a previous call when start = 'W'start='W'. For FP and LP problems, ns need not be initialized.
21:   cw(lencw) – cell array of strings
lencw, the dimension of the array, must satisfy the constraint lencw600lencw600.
Constraint: lencw600lencw600.
22:   iw(leniw) – int64int32nag_int array
leniw, the dimension of the array, must satisfy the constraint leniw600leniw600.
Constraint: leniw600leniw600.
23:   rw(lenrw) – double array
lenrw, the dimension of the array, must satisfy the constraint lenrw600lenrw600.
Constraint: lenrw600lenrw600.

Optional Input Parameters

1:     ne – int64int32nag_int scalar
Default: The dimension of the arrays acol, inda. (An error is raised if these dimensions are not equal.)
The number of nonzero elements in AA.
Constraint: 1nen × m1nen×m.
2:     nname – int64int32nag_int scalar
Default: The dimension of the array names.
The number of column (i.e., variable) and row names supplied in the array names.
nname = 1nname=1
There are no names. Default names will be used in the printed output.
nname = n + mnname=n+m
All names must be supplied.
Constraint: nname = 1nname=1 or n + mn+m.
3:     user – Any MATLAB object
user is not used by nag_opt_qpconvex2_sparse_solve (e04nq), but is passed to qphx. Note that for large objects it may be more efficient to use a global variable which is accessible from the m-files than to use user.

Input Parameters Omitted from the MATLAB Interface

cuser iuser ruser lencw leniw lenrw

Output Parameters

1:     hs(n + mn+m) – int64int32nag_int array
The final states of the variables and slacks (x,s)(x,s). The significance of each possible value of hs(j)hsj is as follows:
hs(j)hsj State of variable jj Normal value of x(j)xj
00 Nonbasic bl(j)blj
11 Nonbasic bu(j)buj
22 Superbasic Between bl(j)blj and bu(j)buj
33 Basic Between bl(j)blj and bu(j)buj
If ninf = 0ninf=0, basic and superbasic variables may be outside their bounds by as much as the value of the optional parameter Feasibility Tolerance. Note that unless the optional parameter Scale Option = 0Scale Option=0 is specified, the optional parameter Feasibility Tolerance applies to the variables of the scaled problem. In this case, the variables of the original problem may be as much as 0.10.1 outside their bounds, but this is unlikely unless the problem is very badly scaled.
Very occasionally some nonbasic variables may be outside their bounds by as much as the optional parameter Feasibility Tolerance, and there may be some nonbasic variables for which x(j)xj lies strictly between its bounds.
If ninf > 0ninf>0, some basic and superbasic variables may be outside their bounds by an arbitrary amount (bounded by sinf if Scale Option = 0Scale Option=0).
2:     x(n + mn+m) – double array
The final values of the variables and slacks (x,s)(x,s).
3:     pi(m) – double array
Contains the dual variables ππ (a set of Lagrange multipliers (shadow prices) for the general constraints).
4:     rc(n + mn+m) – double array
Contains the reduced costs, g − (A − I)
T
π
g-A-ITπ. The vector gg is the gradient of the objective if x is feasible, otherwise it is the gradient of the Phase 1 objective. In the former case, g(i) = 0g(i)=0, for i = n + 1 : mi=n+1:m, hence rc(n + 1 : m) = πrc(n+1:m)=π.
5:     ns – int64int32nag_int scalar
The final number of superbasics. This will be zero for FP and LP problems.
6:     ninf – int64int32nag_int scalar
The number of infeasibilities.
7:     sinf – double scalar
The sum of the scaled infeasibilities. This will be zero if ninf = 0ninf=0, and is most meaningful when Scale Option = 0Scale Option=0.
8:     obj – double scalar
The value of the objective function.
If ninf = 0ninf=0, obj includes the quadratic objective term (1/2)xTHx12xTHx (if any).
If ninf > 0ninf>0, obj is just the linear objective term cTxcTx (if any).
For FP problems, obj is set to zero.
Note that obj does not include contributions from the constant term objadd or the objective row, if any.
9:     user – Any MATLAB object
10:   cw(lencw) – cell array of strings
cw = state . cwcw=state.cwlencw = 600lencw=600.
Communication array, used to store information between calls to nag_opt_qpconvex2_sparse_solve (e04nq).
11:   iw(leniw) – int64int32nag_int array
iw = state . iwiw=state.iwleniw = 600leniw=600.
Communication array, used to store information between calls to nag_opt_qpconvex2_sparse_solve (e04nq).
12:   rw(lenrw) – double array
rw = state . rwrw=state.rwlenrw = 600lenrw=600.
Communication array, used to store information between calls to nag_opt_qpconvex2_sparse_solve (e04nq).
13:   ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).
nag_opt_qpconvex2_sparse_solve (e04nq) returns with ifail = 0ifail=0 if the reduced gradient (rgNorm; see Section [Printed output]) is negligible, the Lagrange multipliers (Lagr Mult; see Section [Printed output]) are optimal, xx satisfies the constraints to the accuracy requested by the value of the optional parameter Feasibility Tolerance and the reduced Hessian factor RR (see Section [Definition of the Working Set and Search Direction]) is nonsingular.

Error Indicators and Warnings

Note: nag_opt_qpconvex2_sparse_solve (e04nq) may return useful information for one or more of the following detected errors or warnings.
Errors or warnings detected by the function:

Cases prefixed with W are classified as warnings and do not generate an error of type NAG:error_n. See nag_issue_warnings.

  ifail = 1ifail=1
The initialization function nag_opt_qpconvex2_sparse_init (e04np) has not been called or at least one of lencw, leniw and lenrw is less than 600600.
  ifail = 2ifail=2
An input parameter is invalid.
  ifail = 3ifail=3
The requested accuracy could not be achieved.
W ifail = 4ifail=4
Weak QP solution found. The final xx is not unique.
This exit will occur when the following are true:
(i) the problem is feasible;
(ii) the reduced gradient is negligible;
(iii) the Lagrange multipliers are optimal; and
(iv) the reduced Hessian is singular or there are some very small multipliers.
This exit cannot occur if HH is positive definite (i.e., f(x) f(x)  is strictly convex).
W ifail = 5ifail=5
The problem is infeasible. The general constraints cannot all be satisfied simultaneously to within the value of the optional parameter Feasibility Tolerance.
Feasibility is measured with respect to the upper and lower bounds on the variables and slacks. The message tells us that among all the points satisfying the general constraints Axs = 0 Ax-s=0 , there is apparently no point that satisfies the bounds on x x  and s s . Violations as small as the Feasibility Tolerance are ignored, but at least one component of x x  or s s  violates a bound by more than the tolerance.
Note:  although the objective function is the sum of infeasibilities (when ninf > 0ninf>0), this sum will not necessarily have been minimized when Elastic Mode = 1 Elastic Mode = 1 .
If Elastic Mode 0 Elastic Mode 0 , nag_opt_qpconvex2_sparse_solve (e04nq) will optimize the QP objective and the sum of infeasibilities, suitably weighted using the optional parameter Elastic Mode. The function will tend to determine a ‘good’ infeasible point if the elastic weight is sufficiently large.
W ifail = 6ifail=6
The problem is unbounded (or badly scaled). For a minimization problem, the objective function is not bounded below in the feasible region.
For linear problems, unboundedness is detected by the simplex method when a nonbasic variable can be increased or decreased by an arbitrary amount without causing a basic variable to violate a bound. Consider adding an upper or lower bound to the variable. Also, examine the constraints that have nonzeros in the associated column, to see if they have been formulated as intended.
Very rarely, the scaling of the problem could be so poor that numerical error will give an erroneous indication of unboundedness. Consider using the optional parameter Scale Option.
  ifail = 7ifail=7
Too many iterations. The value of the optional parameter Iterations Limit is too small.
The Iterations limit was exceeded before the required solution could be found. Check the iteration log to be sure that progress was being made. If so, restart the run using a Basis file that was saved at the end of the run.
  ifail = 8ifail=8
The value of the optional parameter Superbasics Limit is too small. The current set of basic and superbasic variables have been optimized as much as possible and a pricing operation is necessary to continue, but there are already Superbasics Limit superbasics (and no room for any more).
In general, raise the Superbasics Limit ss by a reasonable amount, bearing in mind the storage needed for reduced Hessian (see Section [Definition of the Working Set and Search Direction]). (The Reduced Hessian Dimension hh will also increase to ss unless specified otherwise, and the associated storage will be about (1/2)s212s2 words.) In some cases you may have to set h < sh<s to conserve storage, but beware that the rate of convergence will probably fall off severely.
  ifail = 9ifail=9
The basis is singular after several attempts to factorize it (adding slacks where necessary). Either the problem is badly scaled or the value of the optional parameter LU Factor Tolerance is too large.
  ifail = 10ifail=10
Numerical error in trying to satisfy the general constraints. The basis is very ill-conditioned.
An LU LU  factorization of the basis has just been obtained and used to recompute the basic variables xB xB , given the present values of the superbasic and nonbasic variables. However, a row check has revealed that the resulting solution does not satisfy the current constraints Axs = 0 Ax-s=0  sufficiently well.
This probably means that the current basis is very ill-conditioned. Request the Scale Option if there are any linear constraints and variables.
For certain highly structured basis matrices (notably those with band structure), a systematic growth may occur in the factor U U . Consult the description of Umax, Umin and Growth in Section [Description of Monitoring Information], and set the optional parameter LU Factor Tolerance to 2.02.0 (or possibly even smaller, but not less than 1.01.0).
  ifail = 11ifail=11
An indefinite matrix was detected during the computation of the reduced Hessian factor RR (see Section [Definition of the Working Set and Search Direction]). This may be caused by HH being indefinite. Check also that qphx has been coded correctly and that all relevant elements of HxHx have been assigned their correct values. If qphx is coded correctly and HH is positive semidefinite, the failure may be caused by ill conditioning. Try reducing the values of the optional parameters LU Factor Tolerance and LU Update Tolerance. If there are very large values in HH, check the scaling of the variables and constraints.
  ifail = 12ifail=12
Internal memory allocation failed when attempting to obtain the required workspace. Please contact NAG.
  ifail = 13ifail=13
Internal memory allocation was insufficient. Please contact NAG.
  ifail = 14ifail=14
An error has occurred in the basis package, perhaps indicating incorrect setup of arrays inda and loca. Set the optional parameter Print File and examine the output carefully for further information.
  ifail = 15ifail=15
An unexpected error has occurred. Set the optional parameter Print File and examine the output carefully for further information.

Accuracy

nag_opt_qpconvex2_sparse_solve (e04nq) implements a numerically stable active-set strategy and returns solutions that are as accurate as the condition of the problem warrants on the machine.

Further Comments

This section contains a description of the printed output.

Description of the Printed Output

If Print Level > 0 Print Level > 0 , one line of information is output to the Print File every k k th iteration, where k k  is the specified Print Frequency. A heading is printed before the first such line following a basis factorization. The heading contains the items described below. In this description, a pricing operation is defined to be the process by which one or more nonbasic variables are selected to become superbasic (in addition to those already in the superbasic set). The variable selected will be denoted by jq. If the problem is purely linear, variable jq will usually become basic immediately (unless it should happen to reach its opposite bound and return to the nonbasic set).
If optional parameter Partial Price is in effect, variable jq is selected from App App  or Ipp Ipp , the ppth segments of the constraint matrix
(AI)
A -I .
Label Description
Itn is the iteration count.
pp is the partial-price indicator. The variable selected by the last pricing operation came from the ppth partition of AA and I-I. Note that pp is reset to zero whenever the basis is refactorized.
dj is the value of the reduced gradient (or reduced cost) for the variable selected by the pricing operation at the start of the current iteration.
Algebraically, dj is dj = gj πT aj dj = gj - πT aj , for j = jq j = jq , where gj gj  is the gradient of the current objective function, π π  is the vector of dual variables, and aj aj  is the j j th column of the constraint matrix
(AI)
A -I .
Note that dj is the norm of the reduced-gradient vector at the start of the iteration, just after the pricing operation.
+SBS is the variable jq selected by the pricing operation to be added to the superbasic set.
-SBS is the variable chosen to leave the superbasic set. It has become basic if the entry under -B is nonzero, otherwise it becomes nonbasic.
-BS is the variable removed from the basis to become nonbasic.
Step is the value of the step length αα taken along the current search direction pp. The variables xx have just been changed to x + αpx+αp. If a variable is made superbasic during the current iteration (i.e., +SBS is positive), Step will be the step to the nearest bound. During the optimality phase, the step can be greater than unity only if the reduced Hessian is not positive definite.
Pivot is the rrth element of a vector yy satisfying By = aqBy=aq whenever aqaq (the qqth column of the constraint matrix
(AI)
A -I  replaces the rrth column of the basis matrix BB. Wherever possible, Step is chosen so as to avoid extremely small values of Pivot (since they may cause the basis to be nearly singular). In extreme cases, it may be necessary to increase the value of the optional parameter Pivot Tolerance to exclude very small elements of yy from consideration during the computation of Step.
nInf is the number of violated constraints (infeasibilities) before the present iteration. This number will not increase unless iterations are in elastic mode.
sInf is the sum of infeasibilities before the present iteration. It will usually decrease at each nonzero step, but if nInf decreases by 22 or more, sInf may occasionally increase. However, in elastic mode it will decrease monotonically.
Objective is the value of the current objective function after the present iteration. Note, if Elastic Mode is 22, the heading is Composite Obj.
L+U L is the number of nonzeros in the basis factor LL. Immediately after a basis factorization B = LUB=LU, L contains lenL (see Section [Description of Monitoring Information]). Further nonzeros are added to L when various columns of BB are later replaced. (Thus, L increases monotonically.) U is the number of nonzeros in the basis factor UU. Immediately after a basis factorization B = LUB=LU, U contains lenU (see Section [Description of Monitoring Information]). As columns of BB are replaced, the matrix UU is maintained explicitly (in sparse form). The value of U may fluctuate up or down; in general, it will tend to increase.
ncp is the number of compressions required to recover workspace in the data structure for UU. This includes the number of compressions needed during the previous basis factorization. Normally, ncp should increase very slowly.
The following will be output if the problem is QP or if the superbasic set is non-empty.
Label Description
rgNorm is the largest reduced-gradient among the superbasic variables after the current iteration. During the optimality phase, this will be approximately zero after a unit step.
nS is the current number of superbasic variables.
condHz is a lower bound on the condition number of the reduced Hessian (see Section [Definition of the Working Set and Search Direction]). The larger this number, the more difficult the problem. Attention should be given to the scaling of the variables and the constraints to guard against high values of condHz.

Example

function nag_opt_qpconvex2_sparse_solve_example
start = 'C';
m = int64(8);
n = int64(7);
lenc = int64(0);
ncolh = int64(7);
iobj = int64(8);
objadd = 0;
prob = '        ';
acol = [0.02;
     0.02;
     0.03;
     1;
     0.7;
     0.02;
     0.15;
     -200;
     0.06;
     0.75;
     0.03;
     0.04;
     0.05;
     0.04;
     1;
     -2000;
     0.02;
     1;
     0.01;
     0.08;
     0.08;
     0.8;
     -2000;
     1;
     0.12;
     0.02;
     0.02;
     0.75;
     0.04;
     -2000;
     0.01;
     0.8;
     0.02;
     1;
     0.02;
     0.06;
     0.02;
     -2000;
     1;
     0.01;
     0.01;
     0.97;
     0.01;
     400;
     0.97;
     0.03;
     1;
     400];
inda = [int64(7);5;3;1;6;4;2;8;7;6;5;4;3;2;1;8;2;1;4;3;7;6;8;1;7;3;4;6;2;8;5;6;7;1;2;3;4;8;1;2;3;6;7;8;7;2;1;8];
loca = [int64(1);9;17;24;31;39;45;49];
bl = [0;0;400;100;0;0;0;2000;-1e25;-1e25;-1e25;-1e25;1500;250;-1e25];
bu = [200;2500;800;700;1500;1e25;1e25;2000;60;100;40;30;1e25;300;1e25];
c = [0];
names = {'...X1...'; '...X2...'; '...X3...'; '...X4...'; '...X5...'; ...
         '...X6...'; '...X7...'; '..ROW1..'; '..ROW2..'; '..ROW3..'; ...
         '..ROW4..'; '..ROW5..'; '..ROW6..'; '..ROW7..'; '..COST..'};
helast = [int64(0);0;0;0;0;0;0;0;0;0;0;0;0;0;0];
hs = [int64(0);0;0;0;0;0;0;0;0;0;0;0;0;0;0];
x = [0;0;0;0;0;0;0;0;0;0;0;0;0;0;0];
ns = int64(0);
[cw, iw, rw, ifail] = nag_opt_qpconvex2_sparse_init();
[hsOut, xOut, pi, rc, nsOut, ninf, sinf, obj, cuser, cwOut, iwOut, rwOut, ifail] = ...
     nag_opt_qpconvex2_sparse_solve(start, @qphx, m, n, lenc, ncolh, iobj, objadd, prob, acol, ...
     inda, loca, bl, bu, c, names, helast, hs, x, ns, cw, iw, rw);
 hsOut, xOut, pi, rc, nsOut, ninf, sinf, obj, ifail

function [hx, user] = qphx(ncolh, x, nstate, user)
  hx = zeros(ncolh, 1);

  hx(1) = 2*x(1);
  hx(2) = 2*x(2);
  hx(3) = 2*(x(3)+x(4));
  hx(4) = hx(3);
  hx(5) = 2*x(5);
  hx(6) = 2*(x(6)+x(7));
  hx(7) = hx(6);
 

hsOut =

                    0
                    3
                    2
                    2
                    3
                    3
                    3
                    0
                    3
                    1
                    3
                    3
                    0
                    0
                    3


xOut =

   1.0e+06 *

         0
    0.0003
    0.0006
    0.0002
    0.0004
    0.0003
    0.0002
    0.0020
    0.0000
    0.0001
    0.0000
    0.0000
    0.0015
    0.0002
   -2.9887


pi =

   1.0e+04 *

   -1.2901
         0
   -0.2325
         0
         0
    1.4455
    1.4581
   -0.0001


rc =

   1.0e+04 *

    0.2361
   -0.0000
   -0.0000
    0.0000
   -0.0000
   -0.0000
    0.0000
   -1.2901
         0
   -0.2325
         0
         0
    1.4455
    1.4581
   -0.0001


nsOut =

                    2


ninf =

                    0


sinf =

     0


obj =

  -1.8478e+06


ifail =

                    0


function e04nq_example
start = 'C';
m = int64(8);
n = int64(7);
lenc = int64(0);
ncolh = int64(7);
iobj = int64(8);
objadd = 0;
prob = '        ';
acol = [0.02;
     0.02;
     0.03;
     1;
     0.7;
     0.02;
     0.15;
     -200;
     0.06;
     0.75;
     0.03;
     0.04;
     0.05;
     0.04;
     1;
     -2000;
     0.02;
     1;
     0.01;
     0.08;
     0.08;
     0.8;
     -2000;
     1;
     0.12;
     0.02;
     0.02;
     0.75;
     0.04;
     -2000;
     0.01;
     0.8;
     0.02;
     1;
     0.02;
     0.06;
     0.02;
     -2000;
     1;
     0.01;
     0.01;
     0.97;
     0.01;
     400;
     0.97;
     0.03;
     1;
     400];
inda = [int64(7);5;3;1;6;4;2;8;7;6;5;4;3;2;1;8;2;1;4;3;7;6; ...
        8;1;7;3;4;6;2;8;5;6;7;1;2;3;4;8;1;2;3;6;7;8;7;2;1;8];
loca = [int64(1);9;17;24;31;39;45;49];
bl = [0;0;400;100;0;0;0;2000;-1e25;-1e25;-1e25;-1e25;1500;250;-1e25];
bu = [200;2500;800;700;1500;1e25;1e25;2000;60;100;40;30;1e25;300;1e25];
c = [0];
names = {'...X1...'; '...X2...'; '...X3...'; '...X4...'; '...X5...'; ...
         '...X6...'; '...X7...'; '..ROW1..'; '..ROW2..'; '..ROW3..'; ...
         '..ROW4..'; '..ROW5..'; '..ROW6..'; '..ROW7..'; '..COST..'};
helast = [int64(0);0;0;0;0;0;0;0;0;0;0;0;0;0;0];
hs = [int64(0);0;0;0;0;0;0;0;0;0;0;0;0;0;0];
x = [0;0;0;0;0;0;0;0;0;0;0;0;0;0;0];
ns = int64(0);
[cw, iw, rw, ifail] = e04np;
[hsOut, xOut, pi, rc, nsOut, ninf, sinf, obj, cuser, cwOut, iwOut, rwOut, ifail] = ...
     e04nq(start, @qphx, m, n, lenc, ncolh, iobj, objadd, prob, acol, ...
     inda, loca, bl, bu, c, names, helast, hs, x, ns, cw, iw, rw);
 hsOut, xOut, pi, rc, nsOut, ninf, sinf, obj, ifail

function [hx, user] = qphx(ncolh, x, nstate, user)
  hx = zeros(ncolh, 1);

  hx(1) = 2*x(1);
  hx(2) = 2*x(2);
  hx(3) = 2*(x(3)+x(4));
  hx(4) = hx(3);
  hx(5) = 2*x(5);
  hx(6) = 2*(x(6)+x(7));
  hx(7) = hx(6);
 

hsOut =

                    0
                    3
                    2
                    2
                    3
                    3
                    3
                    0
                    3
                    1
                    3
                    3
                    0
                    0
                    3


xOut =

   1.0e+06 *

         0
    0.0003
    0.0006
    0.0002
    0.0004
    0.0003
    0.0002
    0.0020
    0.0000
    0.0001
    0.0000
    0.0000
    0.0015
    0.0002
   -2.9887


pi =

   1.0e+04 *

   -1.2901
         0
   -0.2325
         0
         0
    1.4455
    1.4581
   -0.0001


rc =

   1.0e+04 *

    0.2361
   -0.0000
   -0.0000
    0.0000
   -0.0000
   -0.0000
    0.0000
   -1.2901
         0
   -0.2325
         0
         0
    1.4455
    1.4581
   -0.0001


nsOut =

                    2


ninf =

                    0


sinf =

     0


obj =

  -1.8478e+06


ifail =

                    0


the remainder of this document is intended for more advanced users. Section [Algorithmic Details] contains a detailed description of the algorithm which may be needed in order to understand Sections [Optional Parameters] and [Description of Monitoring Information]. Section [Optional Parameters] describes the optional parameters which may be set by calls to nag_opt_qpconvex2_sparse_option_string (e04ns), nag_opt_qpconvex2_sparse_option_integer_set (e04nt) and/or nag_opt_qpconvex2_sparse_option_double_set (e04nu). Section [Description of Monitoring Information] describes the quantities which can be requested to monitor the course of the computation.

Algorithmic Details

This section contains a detailed description of the method used by nag_opt_qpconvex2_sparse_solve (e04nq).

Overview

nag_opt_qpconvex2_sparse_solve (e04nq) is based on an inertia-controlling method that maintains a Cholesky factorization of the reduced Hessian (see below). The method is similar to that of Gill and Murray (1978), and is described in detail by Gill et al. (1991). Here we briefly summarise the main features of the method. Where possible, explicit reference is made to the names of variables that are parameters of the function or appear in the printed output.
The method used has two distinct phases: finding an initial feasible point by minimizing the sum of infeasibilities (the feasibility phase), and minimizing the quadratic objective function within the feasible region (the optimality phase). The computations in both phases are performed by the same functions. The two-phase nature of the algorithm is reflected by changing the function being minimized from the sum of infeasibilities (the printed quantity sInf; see Section [Printed output]) to the quadratic objective function (the printed quantity Objective; see Section [Printed output]).
In general, an iterative process is required to solve a quadratic program. Given an iterate (x,s)(x,s) in both the original variables xx and the slack variables ss, a new iterate (x,s)(x-,s-) is defined by
(x)
s
=
(x)
s
+ αp,
x- s- = x s +αp,
(3)
where the step length αα is a non-negative scalar (the printed quantity Step; see Section [Description of Monitoring Information]), and pp is called the search direction. (For simplicity, we shall consider a typical iteration and avoid reference to the index of the iteration.) Once an iterate is feasible (i.e., satisfies the constraints), all subsequent iterates remain feasible.

Definition of the Working Set and Search Direction

At each iterate (x,s)(x,s), a working set of constraints is defined to be a linearly independent subset of the constraints that are satisfied ‘exactly’ (to within the value of the optional parameter Feasibility Tolerance). The working set is the current prediction of the constraints that hold with equality at a solution of the LP or QP problem. Let mWmW denote the number of constraints in the working set (including bounds), and let WW denote the associated mWmW by (n + m)(n+m) working set matrix consisting of the mWmW gradients of the working set constraints.
The search direction is defined so that constraints in the working set remain unaltered for any value of the step length. It follows that pp must satisfy the identity
Wp = 0.
Wp=0.
(4)
This characterisation allows pp to be computed using any nn by nZnZ full-rank matrix ZZ that spans the null space of WW. (Thus, nZ = nmWnZ=n-mW and WZ = 0WZ=0.) The null space matrix ZZ is defined from a sparse LULU factorization of part of WW (see (7) and (8)). The direction pp will satisfy (4) if
p = ZpZ,
p=ZpZ,
(5)
where pZpZ is any nZnZ-vector.
The working set contains the constraints Axs = 0Ax-s=0 and a subset of the upper and lower bounds on the variables (x,s)(x,s). Since the gradient of a bound constraint xjljxjlj or xjujxjuj is a vector of all zeros except for ± 1±1 in position jj, it follows that the working set matrix contains the rows of
(AI)
A -I  and the unit rows associated with the upper and lower bounds in the working set.
The working set matrix WW can be represented in terms of a certain column partition of the matrix
(AI)
A -I  by (conceptually) partitioning the constraints Axs = 0Ax-s=0 so that
BxB + SxS + NxN = 0,
BxB+SxS+NxN=0,
(6)
where BB is a square nonsingular basis and xBxB, xSxS and xNxN are the basic, superbasic and nonbasic variables respectively. The nonbasic variables are equal to their upper or lower bounds at (x,s)(x,s), and the superbasic variables are independent variables that are chosen to improve the value of the current objective function. The number of superbasic variables is nSnS (the printed quantity nS; see Section [Printed output]). Given values of xNxN and xSxS, the basic variables xBxB are adjusted so that (x,s)(x,s) satisfies (6).
If PP is a permutation matrix such that
(AI)
P =
(BSN)
A -I P= B S N , then WW satisfies
WP =
(BSN)
0 0 IN
,
WP= B S N 0 0 IN ,
(7)
where ININ is the identity matrix with the same number of columns as NN.
The null space matrix ZZ is defined from a sparse LULU factorization of part of WW. In particular, ZZ is maintained in ‘reduced gradient’ form, using the LUSOL package (see Gill et al. (1991)) to maintain sparse LULU factors of the basis matrix BB as the BSNBSN partition changes. Given the permutation PP, the null space basis is given by
Z = P
  − B − 1S I 0  
.
Z=P -B-1S I 0 .
(8)
This matrix is used only as an operator, i.e., it is never computed explicitly. Products of the form ZvZv and ZTgZTg are obtained by solving with BB or BTBT. This choice of ZZ implies that nZnZ, the number of ‘degrees of freedom’ at (x,s)(x,s), is the same as nSnS, the number of superbasic variables.
Let gZgZ and HZHZ denote the reduced gradient and reduced Hessian of the objective function:
gZ = ZTg  and  HZ = ZTHZ,
gZ=ZTg  and  HZ=ZTHZ,
(9)
where gg is the objective gradient at (x,s)(x,s). Roughly speaking, gZgZ and HZHZ describe the first and second derivatives of an nSnS-dimensional unconstrained problem for the calculation of pZpZ. (The condition estimator of HZHZ is the quantity condHz in the monitoring file output; see Section [Printed output].)
At each iteration, an upper triangular factor RR is available such that HZ = RTRHZ=RTR. Normally, RR is computed from RTR = ZTHZRTR=ZTHZ at the start of the optimality phase and then updated as the QP working set changes. For efficiency, the dimension of RR should not be excessive (say, nS1000nS1000). This is guaranteed if the number of nonlinear variables is ‘moderate’.
If the QP problem contains linear variables, HH is positive semidefinite and RR may be singular with at least one zero diagonal element. However, an inertia-controlling strategy is used to ensure that only the last diagonal element of RR can be zero. (See Gill et al. (1991) for a discussion of a similar strategy for indefinite quadratic programming.)
If the initial RR is singular, enough variables are fixed at their current value to give a nonsingular RR. This is equivalent to including temporary bound constraints in the working set. Thereafter, RR can become singular only when a constraint is deleted from the working set (in which case no further constraints are deleted until RR becomes nonsingular).

Main Iteration

If the reduced gradient is zero, (x,s)(x,s) is a constrained stationary point on the working set. During the feasibility phase, the reduced gradient will usually be zero only at a vertex (although it may be zero elsewhere in the presence of constraint dependencies). During the optimality phase, a zero reduced gradient implies that xx minimizes the quadratic objective function when the constraints in the working set are treated as equalities. At a constrained stationary point, Lagrange multipliers λλ are defined from the equations
WTλ = g(x).
WTλ=g(x).
(10)
A Lagrange multiplier, λjλj, corresponding to an inequality constraint in the working set is said to be optimal if λjσλjσ when the associated constraint is at its upper bound, or if λjσλj-σ when the associated constraint is at its lower bound, where σσ depends on the value of the optional parameter Optimality Tolerance. If a multiplier is nonoptimal, the objective function (either the true objective or the sum of infeasibilities) can be reduced by continuing the minimization with the corresponding constraint excluded from the working set. (This step is sometimes referred to as ‘deleting’ a constraint from the working set.) If optimal multipliers occur during the feasibility phase but the sum of infeasibilities is nonzero, there is no feasible point and the function terminates immediately with ifail = 3ifail=3.
The special form (7) of the working set allows the multiplier vector λλ, the solution of (10), to be written in terms of the vector
d =
(g)
0
(AT)
I
π =
(gATπ)
π
,
d= g 0 - AT -I π= g-ATπ π ,
(11)
where ππ satisfies the equations BTπ = gBBTπ=gB, and gBgB denotes the basic elements of gg. The elements of ππ are the Lagrange multipliers λjλj associated with the equality constraints Axs = 0Ax-s=0. The vector dNdN of nonbasic elements of dd consists of the Lagrange multipliers λjλj associated with the upper and lower bound constraints in the working set. The vector dSdS of superbasic elements of dd is the reduced gradient gZgZ in (9). The vector dBdB of basic elements of dd is zero, by construction. (The Euclidean norm of dSdS and the final values of dSdS, gg and ππ are the quantities rgNorm, Reduced Gradnt, Obj Gradient and Dual Activity in the monitoring file output; see Section [Description of Monitoring Information].)
If the reduced gradient is not zero, Lagrange multipliers need not be computed and the search direction is given by p = ZpZp=ZpZ (see (8) and (12)). The step length is chosen to maintain feasibility with respect to the satisfied constraints.
There are two possible choices for pZpZ, depending on whether or not HZHZ is singular. If HZHZ is nonsingular, RR is nonsingular and pZpZ in (5) is computed from the equations
RTRpZ = gZ,
RTRpZ=-gZ,
(12)
where gZgZ is the reduced gradient at xx. In this case, (x,s) + p(x,s)+p is the minimizer of the objective function subject to the working set constraints being treated as equalities. If (x,s) + p(x,s)+p is feasible, αα is defined to be unity. In this case, the reduced gradient at (x,s)(x-,s-) will be zero, and Lagrange multipliers are computed at the next iteration. Otherwise, αα is set to αNαN, the step to the ‘nearest’ constraint along pp. This constraint is then added to the working set at the next iteration.
If HZHZ is singular, then RR must also be singular, and an inertia-controlling strategy is used to ensure that only the last diagonal element of RR is zero. (See Gill et al. (1991) for a discussion of a similar strategy for indefinite quadratic programming.) In this case, pZpZ satisfies
pZT HZ pZ = 0   and   gZT pZ0,
pZT HZ pZ=0   and   gZT pZ0,
(13)
which allows the objective function to be reduced by any step of the form (x,s) + α p (x,s) + α p , where α > 0 α > 0 . The vector p = Z pZ p = Z pZ  is a direction of unbounded descent for the QP problem in the sense that the QP objective is linear and decreases without bound along p p . If no finite step of the form (x,s) + α p (x,s) + α p  (where α > 0 α > 0 ) reaches a constraint not in the working set, the QP problem is unbounded and the function terminates immediately with ifail = 6ifail=6. Otherwise, α α  is defined as the maximum feasible step along p p  and a constraint active at (x,s) + α p (x,s) + α p  is added to the working set for the next iteration.
nag_opt_qpconvex2_sparse_solve (e04nq) makes explicit allowance for infeasible constraints. Infeasible linear constraints are detected first by solving a problem of the form
minimize x,v,w eT (v+w)   subject to ​ l x Gx-v+w u,   v0,   w0,
minimizex,v,w eT(v + w)  subject to ​l ≤ (x) Gx − v + w ≤ u,  v ≥ 0,  w ≥ 0,
(14)
where eT = (1,1,,1) eT = (1,1,,1) . This is equivalent to minimizing the sum of the general linear constraint violations subject to the simple bounds. (In the linear programming literature, the approach is often called elastic programming.)

Miscellaneous

If the basis matrix is not chosen carefully, the condition of the null space matrix ZZ in (8) could be arbitrarily high. To guard against this, the function implements a ‘basis repair’ feature in which the LUSOL package (see Gill et al. (1991)) is used to compute the rectangular factorization
(BS)
T
= LU,
B S T=LU,
(15)
returning just the permutation PP that makes PLPTPLPT unit lower triangular. The pivot tolerance is set to require |PLPT|ij2|PLPT|ij2, and the permutation is used to define PP in (7). It can be shown that ZZ is likely to be little more than unity. Hence, ZZ should be well-conditioned regardless of the condition of WW. This feature is applied at the beginning of the optimality phase if a potential BSB-S ordering is known.
The EXPAND procedure (see Gill et al. (1989)) is used to reduce the possibility of cycling at a point where the active constraints are nearly linearly dependent. Although there is no absolute guarantee that cycling will not occur, the probability of cycling is extremely small (see Hall and McKinnon (1996)). The main feature of EXPAND is that the feasibility tolerance is increased at the start of every iteration. This allows a positive step to be taken at every iteration, perhaps at the expense of violating the bounds on (x,s)(x,s) by a small amount.
Suppose that the value of the optional parameter Feasibility Tolerance is δδ. Over a period of KK iterations (where KK is the value of the optional parameter Expand Frequency), the feasibility tolerance actually used by the function (i.e., the working feasibility tolerance) increases from 0.5δ0.5δ to δδ (in steps of 0.5δ / K0.5δ/K).
At certain stages the following ‘resetting procedure’ is used to remove small constraint infeasibilities. First, all nonbasic variables are moved exactly onto their bounds. A count is kept of the number of nontrivial adjustments made. If the count is nonzero, the basic variables are recomputed. Finally, the working feasibility tolerance is reinitialized to 0.5δ0.5δ.
If a problem requires more than KK iterations, the resetting procedure is invoked and a new cycle of iterations is started. (The decision to resume the feasibility phase or optimality phase is based on comparing any constraint infeasibilities with δδ.)
The resetting procedure is also invoked when the function reaches an apparently optimal, infeasible or unbounded solution, unless this situation has already occurred twice. If any nontrivial adjustments are made, iterations are continued.
The EXPAND procedure not only allows a positive step to be taken at every iteration, but also provides a potential choice of constraints to be added to the working set. All constraints at a distance αα (where ααNααN) along pp from the current point are then viewed as acceptable candidates for inclusion in the working set. The constraint whose normal makes the largest angle with the search direction is added to the working set. This strategy helps keep the basis matrix BB well-conditioned.

Optional Parameters

Several optional parameters in nag_opt_qpconvex2_sparse_solve (e04nq) define choices in the problem specification or the algorithm logic. In order to reduce the number of formal parameters of nag_opt_qpconvex2_sparse_solve (e04nq) these optional parameters have associated default values that are appropriate for most problems. Therefore, you need only specify those optional parameters whose values are to be different from their default values.
The remainder of this section can be skipped if you wish to use the default values for all optional parameters.
The following is a list of the optional parameters available. A full description of each optional parameter is provided in Section [Description of the Optional s].
Optional parameters may be specified by calling one, or any, of the functions nag_opt_qpconvex2_sparse_option_string (e04ns), nag_opt_qpconvex2_sparse_option_integer_set (e04nt) and nag_opt_qpconvex2_sparse_option_double_set (e04nu) before a call to nag_opt_qpconvex2_sparse_solve (e04nq), but after a call to nag_opt_qpconvex2_sparse_init (e04np).
nag_opt_qpconvex2_sparse_option_string (e04ns), nag_opt_qpconvex2_sparse_option_integer_set (e04nt) or nag_opt_qpconvex2_sparse_option_double_set (e04nu) can be called to supply options directly, one call being necessary for each optional parameter. nag_opt_qpconvex2_sparse_option_string (e04ns), nag_opt_qpconvex2_sparse_option_integer_set (e04nt) or nag_opt_qpconvex2_sparse_option_double_set (e04nu) should be consulted for a full description of this method of supplying optional parameters.
All optional parameters not specified by you are set to their default values. Optional parameters specified by you are unaltered by nag_opt_qpconvex2_sparse_solve (e04nq) (unless they define invalid values) and so remain in effect for subsequent calls unless altered by you.

Description of the Optional Parameters

For each option, we give a summary line, a description of the optional parameter and details of constraints.
The summary line contains:
Keywords and character values are case and white space insensitive.
Check Frequency  ii
Default = 60=60
Every iith iteration after the most recent basis factorization, a numerical test is made to see if the current solution (x,s)(x,s) satisfies the linear constraints Axs = 0Ax-s=0. If the largest element of the residual vector r = Axsr=Ax-s is judged to be too large, the current basis is refactorized and the basic variables recomputed to satisfy the constraints more accurately. If i0i0, the value i = 99999999i=99999999 is used and effectively no checks are made.
Check Frequency = 1Check Frequency=1 is useful for debugging purposes, but otherwise this option should not be needed.
Crash Option  ii
Default = 3=3
Crash Tolerance  rr
Default = 0.1=0.1
Note that these options do not apply when start = 'W'start='W' (see Section [Parameters]).
If start = 'C'start='C', an internal Crash procedure is used to select an initial basis from various rows and columns of the constraint matrix
(AI)
A -I . The value of ii determines which rows and columns of AA are initially eligible for the basis, and how many times the Crash procedure is called. Columns of I-I are used to pad the basis where necessary.
ii Meaning
00 The initial basis contains only slack variables: B = IB=I.
11 The Crash procedure is called once, looking for a triangular basis in all rows and columns of the matrix AA.
22 The Crash procedure is called once, looking for a triangular basis in rows.
33 The Crash procedure is called twice, treating linear equalities and linear inequalities separately.
If i1i1, certain slacks on inequality rows are selected for the basis first. (If i2i2, numerical values are used to exclude slacks that are close to a bound.) The Crash procedure then makes several passes through the columns of AA, searching for a basis matrix that is essentially triangular. A column is assigned to ‘pivot’ on a particular row if the column contains a suitably large element in a row that has not yet been assigned. (The pivot elements ultimately form the diagonals of the triangular basis.) For remaining unassigned rows, slack variables are inserted to complete the basis.
The Crash Tolerance allows the Crash procedure to ignore certain ‘small’ nonzero elements in each column of AA. If amax amax  is the largest element in column jj, other nonzeros aij aij  in the column are ignored if |aij| amax × r |aij| amax × r . (To be meaningful, rr should be in the range 0 r < 1 0 r < 1 .)
When r > 0.0r>0.0, the basis obtained by the Crash procedure may not be strictly triangular, but it is likely to be nonsingular and almost triangular. The intention is to obtain a starting basis containing more columns of AA and fewer (arbitrary) slacks. A feasible solution may be reached sooner on some problems.
For example, suppose the first mm columns of AA form the matrix shown under LU Factor Tolerance; i.e., a tridiagonal matrix with entries 1-1, 44, 1-1. To help the Crash procedure choose all mm columns for the initial basis, we would specify a Crash Tolerance of rr for some value of r > 0.5r>0.5.
Defaults  
This special keyword may be used to reset all optional parameters to their default values.
Dump File  i1i1
Default = 0=0
Load File  i2i2
Default = 0=0
Optional parameters Dump File and Load File are similar to optional parameters Punch File and Insert File, but they record solution information in a manner that is more direct and more easily modified. A full description of information recorded in optional parameters Dump File and Load File is given in Gill et al. (2005a).
If i1 > 0i1>0, the last solution obtained will be output to the file with unit number ii.
If i2 > 0 i2>0 , the Load File containing basis information will be read. The file will usually have been output previously as a Dump File. The file will not be accessed if optional parameters Old Basis File or Insert File are specified.
Elastic Mode  ii
Default = 1=1
This parameter determines if (and when) elastic mode is to be started. Three elastic modes are available as follows:
ii Meaning
00 Elastic mode is never invoked. nag_opt_qpconvex2_sparse_solve (e04nq) will terminate as soon as infeasibility is detected. There may be other points with significantly smaller sums of infeasibilities.
11 Elastic mode is invoked only if the constraints are found to be infeasible (the default). If the constraints are infeasible, continue in elastic mode with the composite objective determined by the values of the optional parameters Elastic Objective and Elastic Weight.
22 The iterations start and remain in elastic mode. This option allows you to minimize the composite objective function directly without first performing Phase 1 iterations.
The success of this option will depend critically on your choice of Elastic Weight. If Elastic Weight is sufficiently large and the constraints are feasible, the minimizer of the composite objective and the solution of the original problem are identical. However, if the Elastic Weight is not sufficiently large, the minimizer of the composite function may be infeasible, even if a feasible point exists.
Elastic Objective  ii
Default = 1=1
This determines the form of the composite objective f(x) + γ j  (vj + wj) f(x) + γ j ( vj + wj )  in Phase 2 (γγ). Three types of composite objectives are available.
ii Meaning
00 Include only the true objective f(x) f(x)  in the composite objective. This option sets γ = 0γ=0 in the composite objective and allows nag_opt_qpconvex2_sparse_solve (e04nq) to ignore the elastic bounds and find a solution that minimizes f(x)f(x) subject to the non-elastic constraints. This option is useful if there are some ‘soft’ constraints that you would like to ignore if the constraints are infeasible.
11 Use a composite objective defined with γγ determined by the value of Elastic Weight. This value is intended to be used in conjunction with Elastic Mode = 2Elastic Mode=2.
22 Include only the elastic variables in the composite objective. The elastics are weighted by γ = 1γ=1. This choice minimizes the violations of the elastic variables at the expense of possibly increasing the true objective. This option can be used to find a point that minimizes the sum of the violations of a subset of constraints specified by the input array helast.
Elastic Weight  rr
Default = 1.0=1.0
This defines the value of γγ in the composite objective in Phase 2 (γγ).
At each iteration of elastic mode, the composite objective is defined to be
minimize​ ​σ ​ ​ f(x) + γ ​ (sum of infeasibilities);
minimize​ ​σ ​ ​ f(x) + γ ​ (sum of infeasibilities);
where σ = 1σ=1 for Minimize, σ = 1σ=-1 for Maximize, and f(x)f(x) is the quadratic objective.
Note that the effect of γγ is not disabled once a feasible point is obtained.
Expand Frequency  ii
Default = 10000=10000
This option is part of an anti-cycling procedure (see Section [Miscellaneous]) designed to allow progress even on highly degenerate problems.
The strategy is to force a positive step at every iteration, at the expense of violating the constraints by a small amount. Suppose that the value of the optional parameter Feasibility Tolerance is δδ. Over a period of ii iterations, the feasibility tolerance actually used by nag_opt_qpconvex2_sparse_solve (e04nq) (i.e., the working feasibility tolerance) increases from 0.5δ0.5δ to δδ (in steps of 0.5δ / i0.5δ/i).
Increasing the value of ii helps reduce the number of slightly infeasible nonbasic variables (most of which are eliminated during the resetting procedure). However, it also diminishes the freedom to choose a large pivot element (see the description of the optional parameter Pivot Tolerance).
If i0i0, the value i = 99999999i=99999999 is used and effectively no anti-cycling procedure is invoked.
Factorization Frequency  ii
Default = 100(LP)=100(LP) or 50(QP)50(QP)
If i > 0i>0, at most ii basis changes will occur between factorizations of the basis matrix.
For LP problems, the basis factors are usually updated at every iteration. Higher values of i i  may be more efficient on problems that are extremely sparse and well scaled.
For QP problems, fewer basis updates will occur as the solution is approached. The number of iterations between basis factorizations will therefore increase. During these iterations a test is made regularly according to the value of optional parameter Check Frequency to ensure that the linear constraints Axs = 0Ax-s=0 are satisfied. Occasionally, the basis will be refactorized before the limit of ii updates is reached. If i0i0, the default value is used.
Feasibility Tolerance  rr
Default = max {106,sqrt(ε)} = max{10-6,ε}
A feasible problem is one in which all variables satisfy their upper and lower bounds to within the absolute tolerance rr. (This includes slack variables. Hence, the general constraints are also satisfied to within rr.)
nag_opt_qpconvex2_sparse_solve (e04nq) attempts to find a feasible solution before optimizing the objective function. If the sum of infeasibilities cannot be reduced to zero, the problem is assumed to be infeasible. Let sInf be the corresponding sum of infeasibilities. If sInf is quite small, it may be appropriate to raise rr by a factor of 1010 or 100100. Otherwise, some error in the data should be suspected.
Note that if sInf is not small and you have not asked nag_opt_qpconvex2_sparse_solve (e04nq) to minimize the violations of the elastic variables (i.e., you have not specified Elastic Objective = 2Elastic Objective=2), there may be other points that have a significantly smaller sum of infeasibilities. nag_opt_qpconvex2_sparse_solve (e04nq) will not attempt to find the solution that minimizes the sum unless Elastic Objective = 2Elastic Objective=2.
If the constraints and variables have been scaled (see the description of the optional parameter Scale Option), then feasibility is defined in terms of the scaled problem (since it is more likely to be meaningful).
Infinite Bound Size  rr
Default = 1020=1020
If r0r0, rr defines the ‘infinite’ bound infbndinfbnd in the definition of the problem constraints. Any upper bound greater than or equal to infbndinfbnd will be regarded as + + (and similarly any lower bound less than or equal to infbnd-infbnd will be regarded as -). If r < 0r<0, the default value is used.
Iterations Limit  ii
Default = max {10000, 10 max {m,n} } = max{10000, 10 max{m,n} }
The value of ii specifies the maximum number of iterations allowed before termination. Setting i = 0i=0 and Print Level > 0Print Level>0 means that: the workspace needed to start solving the problem will be computed and printed; and feasibility and optimality will be checked. No iterations will be performed. If i < 0i<0, the default value is used.
LU Density Tolerance  r1r1
Default = 0.6 = 0.6  
LU Singularity Tolerance  r2r2
Default = ε(2/3) = ε23  
The density tolerance r1r1 is used during LULU factorization of the basis matrix. Columns of LL and rows of UU are formed one at a time, and the remaining rows and columns of the basis are altered appropriately. At any stage, if the density of the remaining matrix exceeds r1r1, the Markowitz strategy for choosing pivots is terminated. The remaining matrix is factored by a dense LULU procedure. Raising the density tolerance towards 1.01.0 may give slightly sparser LULU factors, with a slight increase in factorization time.
If r2 > 0r2>0, r2r2 defines the singularity tolerance used to guard against ill-conditioned basis matrices. After BB is refactorized, the diagonal elements of UU are tested as follows. If |ujj|r2|ujj|r2 or |ujj| < r2maxi |uij||ujj|<r2maxi|uij|, the jjth column of the basis is replaced by the corresponding slack variable. If r20r20, the default value is used.
LU Factor Tolerance  r1r1
Default = 100.0=100.0
LU Update Tolerance  r2r2
Default = 10.0=10.0
The values of r1r1 and r2r2 affect the stability and sparsity of the basis factorization B = LUB=LU, during refactorization and updates respectively. The lower triangular matrix LL is a product of matrices of the form
(10)
μ 1
1 0 μ 1
where the multipliers μμ will satisfy |μ|ri|μ|ri. The default values of r1r1 and r2r2 usually strike a good compromise between stability and sparsity. They must satisfy r1r1, r21.0r21.0.
For large and relatively dense problems, r1 = 10.0​ or ​5.0r1=10.0​ or ​5.0 (say) may give a useful improvement in stability without impairing sparsity to a serious degree.
For certain very regular structures (e.g., band matrices) it may be necessary to reduce r1​ and/or ​r2r1​ and/or ​r2 in order to achieve stability. For example, if the columns of AA include a sub-matrix of the form
  4 − 1 − 1 4 − 1 − 1 4 − 1 … … … − 1 4 − 1 − 1 4  
,
4 -1 -1 4 -1 -1 4 -1 -1 4 -1 -1 4 ,
one should set both r1r1 and r2r2 to values in the range 1.0ri < 4.01.0ri<4.0.
LU Partial Pivoting  
Default
LU Complete Pivoting  
LU Rook Pivoting  
The LULU factorization implements a Markowitz-type search for pivots that locally minimize the fill-in subject to a threshold pivoting stability criterion. The default option is to use threshold partial pivoting. The options LU Complete Pivoting and LU Rook Pivoting are more expensive but more stable and better at revealing rank, as long as the LU Factor Tolerance is not too large (say < 2.0<2.0).
Minimize  
Default
Maximize  
Feasible Point  
This option specifies the required direction of the optimization. It applies to both linear and nonlinear terms (if any) in the objective function. Note that if two problems are the same except that one minimizes f(x)f(x) and the other maximizes f(x)-f(x), their solutions will be the same but the signs of the dual variables πiπi and the reduced gradients djdj (see Section [Main Iteration]) will be reversed.
The option Feasible Point means ‘ignore the objective function, while finding a feasible point for the linear constraints’. It can be used to check that the constraints are feasible without altering the call to nag_opt_qpconvex2_sparse_solve (e04nq).
New Basis File  i1i1
Default = 0=0
Backup Basis File  i2i2
Default = 0=0
Save Frequency  i3i3
Default = 100=100
Optional parameters New Basis File and Backup Basis File are sometimes referred to as basis maps. They contain the most compact representation of the state of each variable. They are intended for restarting the solution of a problem at a point that was reached by an earlier run. For nontrivial problems, it is advisable to save basis maps at the end of a run, in order to restart the run if necessary.
If i1 > 0 i1>0 , a basis map will be saved on file i1 i1  every i3 i3 th iteration, where i3 i3  is the Save Frequency. The first record of the file will contain the word PROCEEDING if the run is still in progress. A basis map will also be saved at the end of a run, with some other word indicating the final solution status.
Use of i2 > 0 i2>0  is intended as a safeguard against losing the results of a long run. Suppose that a New Basis File is being saved every 100100 (Save Frequency) iterations, and that nag_opt_qpconvex2_sparse_solve (e04nq) is about to save such a basis at iteration 20002000. It is conceivable that the run may be interrupted during the next few milliseconds (in the middle of the save). In this case the Basis file will be corrupted and the run will have been essentially wasted.
To eliminate this risk, both a New Basis File and a Backup Basis File may be specified. The following would be suitable for the above example:
 
Backup Basis File 11 
New Basis File 12
The current basis will then be saved every 100100 iterations, first on File 12 and then immediately on File 11. If the run is interrupted at iteration 20002000 during the save on File 12, there will still be a usable basis on File 11 (corresponding to iteration 19001900).
Note that a new basis will be saved in New Basis File at the end of a run if it terminates normally, but it will not be saved in Backup Basis File. In the above example, if an optimum solution is found at iteration 20502050 (or if the iteration limit is 20502050), the final basis on File 12 will correspond to iteration 20502050, but the last basis saved on File 11 will be the one for iteration 20002000.
A full description of information recorded in New Basis File and Backup Basis File is given in Gill et al. (2005a).
Nolist  
Default
List  
Normally each optional parameter specification is printed to unit Print File as it is supplied. Optional parameter Nolist may be used to suppress the printing and optional parameter List may be used to restore printing.
Old Basis File  ii
Default = 0=0
If i > 0 i>0 , the basis maps information will be obtained from this file. The file will usually have been output previously as a New Basis File or Backup Basis File. A full description of information recorded in New Basis File and Backup Basis File is given in Gill et al. (2005a).
The file will not be acceptable if the number of rows or columns in the problem has been altered.
Optimality Tolerance  rr
Default = max {106,sqrt(ε)} = max{10-6,ε}
This is used to judge the size of the reduced gradients dj = gj ajT π dj=gj - ajT π , where gj gj  is the jjth component of the gradient, aj aj  is the associated column of the constraint matrix
(AI)
A -I , and π π  is the set of dual variables.
By construction, the reduced gradients for basic variables are always zero. The problem will be declared optimal if the reduced gradients for nonbasic variables at their lower or upper bounds satisfy
dj / π r  or  dj / π r
dj / π -r  or  dj / π r
respectively, and if |dj| / π r |dj| / π r  for superbasic variables.
In the above tests, π π  is a measure of the size of the dual variables. It is included to make the tests independent of a scale factor on the objective function. The quantity π π  actually used is defined by
m
π = max (σ / sqrt(m),1), where ​σ = |πi|,
i = 1
π=max(σ/m,1) , where ​ σ= i=1 m |πi|,
so that only large scale factors are allowed for.
If the objective is scaled down to be very small, the optimality test reduces to comparing dj dj  against 0.01r 0.01r .
Partial Price  ii
Default = 10(LP)=10(LP) or 1(QP)1(QP)
This option is recommended for large FP or LP problems that have significantly more variables than constraints (i.e., nmnm). It reduces the work required for each pricing operation (i.e., when a nonbasic variable is selected to enter the basis). If i = 1i=1, all columns of the constraint matrix
(AI)
A -I  are searched. If i > 1i>1, AA and II are partitioned to give ii roughly equal segments Aj,IjAj,Ij, for j = 1,2,,ij=1,2,,i (modulo ii). If the previous pricing search was successful on Aj1,Ij1Aj-1,Ij-1, the next search begins on the segments AjAj and IjIj. If a reduced gradient is found that is larger than some dynamic tolerance, the variable with the largest such reduced gradient (of appropriate sign) is selected to enter the basis. If nothing is found, the search continues on the next segments Aj + 1,Ij + 1Aj+1,Ij+1, and so on. If i0i0, the default value is used.
Pivot Tolerance  rr
Default = ε(2/3)=ε23
Broadly speaking, the pivot tolerance is used to prevent columns entering the basis if they would cause the basis to become almost singular.
When x x  changes to x + αp x+αp  for some search direction p p , a ‘ratio test’ determines which component of x x  reaches an upper or lower bound first. The corresponding element of p p  is called the pivot element. Elements of p p  are ignored (and therefore cannot be pivot elements) if they are smaller than the pivot tolerance r r .
It is common for two or more variables to reach a bound at essentially the same time. In such cases, the optional parameter Feasibility Tolerance (say t t ) provides some freedom to maximize the pivot element and thereby improve numerical stability. Excessively small values of t t  should therefore not be specified. To a lesser extent, the optional parameter Expand Frequency (say f f ) also provides some freedom to maximize the pivot element. Excessively large values of f f  should therefore not be specified.
Print File  ii
Default = 0=0
If i > 0 i>0 , the following information is output to i i  during the solution of each problem:
a listing of the optional parameters;
some statistics about the problem;
the amount of storage available for the LU LU  factorization of the basis matrix;
notes about the initial basis resulting from a Crash procedure or a Basis file;
the iteration log;
basis factorization statistics;
the exit ifail condition and some statistics about the solution obtained;
the printed solution, if requested.
The last four items are described in Sections [Further Comments] and [Description of Monitoring Information]. Further brief output may be directed to the Summary File.
Print Frequency  ii
Default = 100=100
If i > 0 i>0 , one line of the iteration log will be printed every i i th iteration. A value such as i = 10 i=10  is suggested for those interested only in the final solution. If i0i0, the value of i = 99999999i=99999999 is used and effectively no checks are made.
Print Level  ii
Default = 1=1
This controls the amount of printing produced by nag_opt_qpconvex2_sparse_solve (e04nq) as follows.
ii Meaning
0 No output except error messages. If you want to suppress all output, set Print File = 0Print File=0.
= 1=1 The set of selected options, problem statistics, summary of the scaling procedure, information about the initial basis resulting from a Crash or a Basis file, a single line of output at each iteration (controlled by the optional parameter Print Frequency), and the exit condition with a summary of the final solution.
1010 Basis factorization statistics.
Punch File  i1i1
Default = 0=0
Insert File  i2i2
Default = 0=0
These files provide compatibility with commercial mathematical programming systems. The Punch File from a previous run may be used as an Insert File for a later run on the same problem. A full description of information recorded in Insert File and Punch File is given in Gill et al. (2005a).
If i1 > 0 i1>0 , the final solution obtained will be output to file i1 i1 . For linear programs, this format is compatible with various commercial systems.
If i2 > 0 i2>0 , the Insert File containing basis information will be read. The file will usually have been output previously as a Punch File. The file will not be accessed if Old Basis File is specified.
QPSolver Cholesky  
Default
QPSolver CG  
QPSolver QN  
Specifies the active-set algorithm used to solve the quadratic program in Phase 2 (γγ). QPSolver Cholesky holds the full Cholesky factor RR of the reduced Hessian ZTHZZTHZ. As the QP iterations proceed, the dimension of RR changes with the number of superbasic variables. If the number of superbasic variables needs to increase beyond the value of Reduced Hessian Dimension, the reduced Hessian cannot be stored and the solver switches to QPSolver CG. The Cholesky solver is reactivated if the number of superbasics stabilizes at a value less than Reduced Hessian Dimension.
QPSolver QN solves the QP using a quasi-Newton method. In this case, RR is the factor of a quasi-Newton approximate Hessian.
QPSolver CG uses an active-set method similar to QPSolver QN, but uses the conjugate-gradient method to solve all systems involving the reduced Hessian.
The Cholesky QP solver is the most robust, but may require a significant amount of computation if there are many superbasics.
The quasi-Newton QP solver does not require computation of the exact RR at the start of Phase 2 (γγ). It may be appropriate when the number of superbasics is large but relatively few iterations are needed to reach a solution (e.g., if nag_opt_qpconvex2_sparse_solve (e04nq) is called with a Warm Start).
The conjugate-gradient QP solver is appropriate for problems with many degrees of freedom (say, more than 20002000 superbasics).
Reduced Hessian Dimension  ii
Default = 1 (LP) ​ or ​ min (2000,nH + 1,n) (QP) = 1 (LP) ​ or ​ min(2000,nH+1,n) (QP)
This specifies that an ii by ii triangular matrix RR (to define the reduced Hessian according to RT R = ZT HZ RT R = ZT HZ ). is to be available for use by the Cholesky QP solver.
Scale Option  ii
Default = 2=2
Scale Tolerance  rr
Default = 0.9=0.9
Scale Print  
Three scale options are available as follows:
ii Meaning
0 No scaling. This is recommended if it is known that x x  and the constraint matrix never have very large elements (say, larger than 100100).
1 The constraints and variables are scaled by an iterative procedure that attempts to make the matrix coefficients as close as possible to 1.01.0 (see Fourer (1982)). This will sometimes improve the performance of the solution procedures.
2 The constraints and variables are scaled by the iterative procedure. Also, a certain additional scaling is performed that may be helpful if the right-hand side b b  or the solution x x  is large. This takes into account columns of
(AI)
A -I  that are fixed or have positive lower bounds or negative upper bounds.
Optional parameter Scale Tolerance affects how many passes might be needed through the constraint matrix. On each pass, the scaling procedure computes the ratio of the largest and smallest nonzero coefficients in each column:
ρj = max |aij| / min |aij|(aij0).
j i
ρj=maxj |aij| / mini |aij| ( aij 0 ) .
If maxj  ρj maxj ρj is less than r r  times its previous value, another scaling pass is performed to adjust the row and column scales. Raising r r  from 0.90.9 to 0.990.99 (say) usually increases the number of scaling passes through A A . At most 1010 passes are made. The value of rr should lie in the range 0 < r < 10<r<1.
Scale Print causes the row scales r(i)r(i) and column scales c(j)c(j) to be printed to Print File, if System Information Yes has been specified. The scaled matrix coefficients are aij = aij c(j) / r(i) a- ij = aij c(j) / r(i), and the scaled bounds on the variables and slacks are lj = lj / c(j) l-j = lj / c(j) , uj = uj / c(j) u-j = uj / c(j) , where c(j) = r(jn) c(j) = r(j-n)  if j > nj>n.
Solution Yes  
Default
Solution No  
This option determines if the final obtained solution is to be output to the Print File. Note that the Solution File option operates independently.
Solution File  ii
Default = 0=0
If i > 0 i>0 , the final solution will be output to file i i  (whether optimal or not).
To see more significant digits in the printed solution, it will sometimes be useful to make i i  refer to the system Print File.
Summary File  i1i1
Default = 0=0
Summary Frequency  i2i2
Default = 100=100
If i1 > 0 i1>0 , the Summary File is output to file i1 i1 , including a line of the iteration log every i2 i2 th iteration. In an interactive environment, it is useful to direct this output to the terminal, to allow a run to be monitored online. (If something looks wrong, the run can be manually terminated.) Further details are given in Section [Description of Monitoring Information]. If i20i20, the value of i2 = 99999999i2=99999999 is used and effectively no checks are made.
Superbasics Limit  ii
Default = 1 (LP) ​ or ​ min { nH + 1 ,n} (QP) = 1 (LP) ​ or ​ min{ nH + 1 ,n} (QP)
This places a limit on the storage allocated for superbasic variables. Ideally, i i  should be set slightly larger than the ‘number of degrees of freedom’ expected at an optimal solution.
For linear programs, an optimum is normally a basic solution with no degrees of freedom. (The number of variables lying strictly between their bounds is no more than m m , the number of general constraints.) The default value of i i  is therefore 11.
For quadratic problems, the number of degrees of freedom is often called the ‘number of independent variables’. Normally, i i  need not be greater than nH + 1 nH+1 , where nHnH is the number of leading nonzero columns of H H . For many problems, i i  may be considerably smaller than nHnH. This will save storage if nHnH is very large.
Suppress Parameters  
Normally nag_opt_qpconvex2_sparse_solve (e04nq) prints the options file as it is being read, and then prints a complete list of the available keywords and their final values. The optional parameter Suppress Parameters tells nag_opt_qpconvex2_sparse_solve (e04nq) not to print the full list.
System Information No  
Default
System Information Yes  
This option prints additional information on the progress of major and minor iterations, and Crash statistics. See Section [Description of Monitoring Information].
Timing Level  ii
Default = 0=0
If i > 0 i>0 , some timing information will be output to the Print file, if Print File > 0Print File>0.
Unbounded Step Size  rr
Default = infbnd=infbnd
If r > 0r>0, rr specifies the magnitude of the change in variables that will be considered a step to an unbounded solution. (Note that an unbounded solution can occur only when the Hessian is not positive definite.) If the change in xx during an iteration would exceed the value of rr, the objective function is considered to be unbounded below in the feasible region. If r0r0, the default value is used. See Infinite Bound Size for the definition of infbndinfbnd.

Description of Monitoring Information

This section describes the intermediate printout and final printout which constitutes the monitoring information produced by nag_opt_qpconvex2_sparse_solve (e04nq). (See also the description of the optional parameters Print File and Print Level.) You can control the level of printed output.

Crash Statistics

When Print Level 10 Print Level 10 , Print File > 0 Print File > 0  and System Information Yes has been specified, the following lines of intermediate printout (less than 120120 characters) are produced on the unit number specified by optional parameter Print File whenever start = 'C'start='C' (see Section [Parameters]). They refer to the number of columns selected by the Crash procedure during each of several passes through AA, whilst searching for a triangular basis matrix.
Label Description
Slacks is the number of slacks selected initially.
Free cols is the number of free columns in the basis, including those whose bounds are rather far apart.
Preferred is the number of ‘preferred’ columns in the basis (i.e., hs(j) = 3hsj=3 for some jnjn). It will be a subset of the columns for which hs(j) = 3hsj=3 was specified.
Unit is the number of unit columns in the basis.
Double is the number of double columns in the basis.
Triangle is the number of triangular columns in the basis.
Pad is the number of slacks used to pad the basis (to make it a nonsingular triangle).

Basis Factorization Statistics

When Print Level 10 Print Level 10  and Print File > 0 Print File > 0 , the first seven items of intermediate printout in the list below are produced on the unit number specified by optional parameter Print File whenever the matrix BB or BS = (BS)
T
BS= B S T is factorized. Gaussian elimination is used to compute an LULU factorization of BB or BSBS, where PLPTPLPT is a lower triangular matrix and PUQPUQ is an upper triangular matrix for some permutation matrices PP and QQ. The factorization is stabilized in the manner described under the optional parameter LU Factor Tolerance. In addition, if System Information Yes has been specified, the entries from Elems onwards are also output.
Label Description
Factor the number of factorizations since the start of the run.
Demand a code giving the reason for the present factorization, as follows:
Code Meaning
0 First LU LU  factorization.
1 The number of updates reached the Factorization Frequency.
2 The nonzeros in the updated factors have increased significantly.
7 Not enough storage to update factors.
10 Row residuals are too large (see the description of the optional parameter Check Frequency).
11 Ill-conditioning has caused inconsistent results.
Itn is the current minor iteration number.
Nonlin is the number of nonlinear variables in the current basis B B .
Linear is the number of linear variables in B B .
Slacks is the number of slack variables in B B .
B, BR, BS or BT factorize is the type of LU LU  factorization.
B periodic factorization of the basis B B .
BR more careful rank-revealing factorization of B B  using threshold rook pivoting. This occurs mainly at the start, if the first basis factors seem singular or ill-conditioned. Followed by a normal B factorize.
BS BS BS  is factorized to choose a well-conditioned B B  from the current (B S) ( B S ) . Followed by a normal B factorize.
BT same as BS except the current B B  is tried first and accepted if it appears to be not much more ill-conditioned than after the previous BS factorize.
m is the number of rows in B B  or BS BS .
n is the number of columns in B B  or BS BS . Preceded by ‘=’ or ‘>’ respectively.
Elems is the number of nonzero elements in B B  or BS BS .
Amax is the largest nonzero in B B  or BS BS .
Density is the percentage nonzero density of B B  or BS BS .
Merit/MerRP/MerCP Merit is the average Markowitz merit count for the elements chosen to be the diagonals of PUQ PUQ . Each merit count is defined to be (c1) (r1) (c-1) (r-1)  where c c  and r r  are the number of nonzeros in the column and row containing the element at the time it is selected to be the next diagonal. Merit is the average of n such quantities. It gives an indication of how much work was required to preserve sparsity during the factorization. If LU Complete Pivoting or LU Rook Pivoting has been selected, this heading is changed to MerCP, respectively MerRP.
lenL is the number of nonzeros in L L .
L+U is the number of nonzeros representing the basis factors L L  and U U . Immediately after a basis factorization B = LU B=LU , this is lenL+lenU, the number of subdiagonal elements in the columns of a lower triangular matrix and the number of diagonal and superdiagonal elements in the rows of an upper-triangular matrix. Further nonzeros are added to L when various columns of B B  are later replaced. As columns of B B  are replaced, the matrix U U  is maintained explicitly (in sparse form). The value of L will steadily increase, whereas the value of U may fluctuate up or down. Thus the value of L+U may fluctuate up or down (in general, it will tend to increase).
Cmpressns is the number of times the data structure holding the partially factored matrix needed to be compressed to recover unused storage. Ideally this number should be zero. If it is more than 33 or 44, the amount of workspace available to nag_opt_qpconvex2_sparse_solve (e04nq) should be increased for efficiency.
Incres is the percentage increase in the number of nonzeros in L L  and U U  relative to the number of nonzeros in B B  or BS BS .
Utri is the number of triangular rows of B B  or BS BS  at the top of U U .
lenU the number of nonzeros in U U , including its diagonals.
Ltol is the largest subdiagonal element allowed in L L . This is the specified LU Factor Tolerance or a smaller value that is currently being used for greater stability.
Umax the maximum nonzero element in U U .
Ugrwth is the ratio Umax / Amax Umax/Amax , which ideally should not be substantially larger than 10.010.0 or 100.0100.0. If it is orders of magnitude larger, it may be advisable to reduce the LU Factor Tolerance to 5.05.0, 4.04.0, 3.03.0 or 2.02.0, say (but bigger than 1.01.0).
As long as Lmax is not large (say 5.05.0 or less), max (Amax,Umax) / DUmin max(Amax,Umax) / DUmin  gives an estimate of the condition number B B . If this is extremely large, the basis is nearly singular. Slacks are used to replace suspect columns of B B  and the modified basis is refactored.
Ltri is the number of triangular columns of B B  or BS BS  at the left of L L .
dense1 is the number of columns remaining when the density of the basis matrix being factorized reached 0.30.3.
Lmax is the actual maximum subdiagonal element in L L  (bounded by Ltol).
Akmax is the largest nonzero generated at any stage of the LU LU  factorization. (Values much larger than Amax indicate instability.) Akmax is not printed if LU Partial Pivoting is selected.
Agrwth is the ratio Akmax / Amax Akmax/Amax . Values much larger than 100100 (say) indicate instability. Agrwth is not printed if LU Partial Pivoting is selected.
bump is the size of the block to be factorized nontrivially after the triangular rows and columns of B B  or BS BS  have been removed.
dense2 is the number of columns remaining when the density of the basis matrix being factorized reached 0.60.6. (The Markowitz pivot strategy searches fewer columns at that stage.)
DUmax is the largest diagonal of PUQ PUQ .
DUmin is the smallest diagonal of PUQ PUQ .
condU the ratio DUmax / DUmin DUmax/DUmin , which estimates the condition number of U U  (and of B B  if Ltol is less than 5.05.0, say).

Basis Map

When Print Level 10 Print Level 10  and Print File > 0 Print File > 0 , the following lines of intermediate printout (less than 8080 characters) are produced on the unit number specified by optional parameter Print File. They refer to the elements of the names array (see Section [Parameters]).
Label Description
Name gives the name for the problem (blank if problem unnamed).
Infeasibilities gives the number of infeasibilities. Printed only if the final point is infeasible.
Objective Value gives the objective value at the final point (or the value of the sum of infeasibilities). Printed only if the final point is feasible.
Status gives the exit status for the problem (i.e., Optimal soln, Weak soln, Unbounded, Infeasible, Excess itns, Error condn or Feasble soln) followed by details of the direction of the optimization (i.e., (Min) or (Max)).
Iteration gives the iteration number when the file was created.
Superbasics gives the number of superbasic variables.
Objective gives the name of the free row for the problem (blank if objective unnamed).
RHS gives the name of the constraint right-hand side for the problem (blank if objective unnamed).
Ranges gives the name of the ranges for the problem (blank if objective unnamed).
Bounds gives the name of the bounds for the problem (blank if objective unnamed).

Solution Output

At the end of a run, the final solution will be output to the Print file. Some header information appears first to identify the problem and the final state of the optimization procedure. A ROWS section and a COLUMNS section then follow, giving one line of information for each row and column.

The ROWS section

General constraints take the form l Ax u l Ax u . The i i th constraint is therefore of the form
α νi x β ,
α νi x β ,
where νi νi  is the i i th row of A A .
Internally, the constraints take the form Ax s = 0 Ax - s = 0 , where s s  is the set of slack variables (which happen to satisfy the bounds l s u l s u ). For the i i th constraint, the slack variable si si  is directly available, and it is sometimes convenient to refer to its state. It should satisfy α si β α si β . A fullstop (.) is printed for any numerical value that is exactly zero.
Label Description
Number is the value of n + in+i. (This is used internally to refer to sisi in the intermediate output.)
Row gives the name of νiνi.
State the state of νi νi  (the state of si si  relative to the bounds α α  and β β ). The various states possible are as follows:
LL si si  is nonbasic at its lower limit, α α .
UL si si  is nonbasic at its upper limit, β β .
EQ si si  is nonbasic and fixed at the value α = β α=β .
FR si si  is nonbasic and currently zero, even though it is free to take any value between its bounds α α  and β β .
BS si si  is basic.
SBS si si  is superbasic.
A key is sometimes printed before State. Note that unless the optional parameter Scale Option = 0Scale Option=0 is specified, the tests for assigning a key are applied to the variables of the scaled problem.
A Alternative optimum possible. The variable is nonbasic, but its reduced gradient is essentially zero. This means that if the variable were allowed to start moving away from its bound, there would be no change in the value of the objective function. The values of the other free variables might change, giving a genuine alternative solution. However, if there are any degenerate variables (labelled D), the actual change might prove to be zero, since one of them could encounter a bound immediately. In either case, the values of the Lagrange multipliers might also change.
D Degenerate. The variable is basic or superbasic, but it is equal (or very close) to one of its bounds.
I Infeasible. The variable is basic or superbasic and is currently violating one of its bounds by more than the value of the Feasibility Tolerance.
N Not precisely optimal. If the slack is superbasic, the dual variable πiπi is not sufficiently small, as measured by the Optimality Tolerance. If the slack is nonbasic, πiπi is not sufficiently positive or negative. If a loose Optimality Tolerance has been used, or if iterations were terminated before optimality, this key might be helpful in deciding whether or not to restart the run.
Activity is the value of νixνix at the final iterate.
Slack Activity is the value by which the row differs from its nearest bound. (For the free row (if any), it is set to Activity.)
Lower Limit is αα, the lower bound specified for the variable sisi. None indicates that bl(j)infbndblj-infbnd.
Upper Limit is ββ, the upper bound specified for the variable sisi. None indicates that bu(j)infbndbujinfbnd.
Dual Activity is the value of the dual variable πiπi (the Lagrange multiplier for νiνi; see Section [Main Iteration]). For FP problems, πiπi is set to zero.
i gives the index ii of the iith row.

The COLUMNS Section

Let the j j th component of x x  be the variable xj xj  and assume that it satisfies the bounds α xj β α xj β . A fullstop (.) is printed for any numerical value that is exactly zero.
Label Description
Number is the column number jj. (This is used internally to refer to xjxj in the intermediate output.)
Column gives the name of xjxj.
State the state of xj xj  relative to the bounds α α  and β β . The various states possible are as follows:
LL xj xj  is nonbasic at its lower limit, α α .
UL xj xj  is nonbasic at its upper limit, β β .
EQ xj xj  is nonbasic and fixed at the value α = β α = β .
FR xj xj  is nonbasic and currently zero, even though it is free to take any value between its bounds α α  and β β .
BS xj xj  is basic.
SBS xj xj  is superbasic.
A key is sometimes printed before State. Note that unless the optional parameter Scale Option = 0Scale Option=0 is specified, the tests for assigning a key are applied to the variables of the scaled problem.
A Alternative optimum possible. The variable is nonbasic, but its reduced gradient is essentially zero. This means that if the variable were allowed to start moving away from its bound, there would be no change in the value of the objective function. The values of the other free variables might change, giving a genuine alternative solution. However, if there are any degenerate variables (labelled D), the actual change might prove to be zero, since one of them could encounter a bound immediately. In either case, the values of the Lagrange multipliers might also change.
D Degenerate. The variable is basic or superbasic, but it is equal (or very close) to one of its bounds.
I Infeasible. The variable is basic or superbasic and is currently violating one of its bounds by more than the value of the Feasibility Tolerance.
N Not precisely optimal. If the slack is superbasic, the dual variable πiπi is not sufficiently small, as measured by the Optimality Tolerance. If the slack is nonbasic, πiπi is not sufficiently positive or negative. If a loose Optimality Tolerance has been used, or if iterations were terminated before optimality, this key might be helpful in deciding whether or not to restart the run.
Activity is the value of xjxj at the final iterate.
Obj Gradient is the value of gjgj at the final iterate. For FP problems, gjgj is set to zero.
Lower Limit is the lower bound specified for the variable. None indicates that bl(j)infbndblj-infbnd.
Upper Limit is the upper bound specified for the variable. None indicates that bu(j)infbndbujinfbnd.
Reduced Gradnt is the value of djdj at the final iterate (see Section [Main Iteration]). For FP problems, djdj is set to zero.
m + j is the value of m + jm+j.
Note:  if two problems are the same except that one minimizes f(x) f(x)  and the other maximizes f(x) -f(x) , their solutions will be the same but the signs of the dual variables πi πi  and the reduced gradients djdj will be reversed.

The Solution File

If a positive Solution File is specified, the information contained in a printed solution may also be output to the relevant file (which may be the Print file if so desired). Infinite Upper and Lower limits appear as ± 1020 ±1020  rather than None. Other real values are output with format 1pe16.6. The maximum line length is 111111 characters, including what would be the carriage-control character if the file were printed.
A Solution file is intended to be read from disk by a self-contained program that extracts and saves certain values as required for possible further computation. Typically the first 1414 lines would be ignored. Each subsequent line may be read using
FORMAT (i8, 2x, 2a4, 1x, a1, 1x, a3, 5e16.6, i7)
adapted to suit the occasion. The end of the ROWS section is marked by a line that starts with a 11 and is otherwise blank. If this and the next 44 lines are skipped, the COLUMNS section (see Section [The COLUMNS Section]) can then be read under the same format. (There should be no need to use any BACKSPACE statements.)

The Summary File

If Summary File ff is specified with f > 0f>0, certain brief information will be output to unit ff. When nag_opt_qpconvex2_sparse_solve (e04nq) is run interactively, unit ff will usually be the terminal. For batch jobs a disk file should be used, to retain a concise log of each run if desired. (A Summary File is more easily perused than the associated Print file).
A Summary file (like the Print file) is not rewound after a problem has been processed. The maximum line length is 7272 characters, including a carriage-control character in column 1.
The following information is included:
  1. The optional parameters supplied via the option setting functions, if any;
  2. The Basis file loaded, if any;
  3. The status of the solution after each basis factorization (whether feasible; the objective value; the number of function calls so far);
  4. The same information every k k th iteration, where k k  is the specified Summary Frequency;
  5. Warnings and error messages;
  6. The exit condition and a summary of the final solution.
Item 4. is preceded by a blank line, but item 5. is not.
The meaning of the printout for linear constraints is the same as that given above for variables, with ‘variable’ replaced by ‘constraint’, nn replaced by mm, names(j)namesj replaced by names(n + j)namesn+j, bl(j)blj and bu(j)buj are replaced by bl(n + j)bln+j and bu(n + j)bun+j respectively, and with the following change in the heading:
Constrnt gives the name of the linear constraint.
Note that movement off a constraint (as opposed to a variable moving away from its bound) can be interpreted as allowing the entry in the Residual column to become positive.
Numerical values are output with a fixed number of digits; they are not guaranteed to be accurate to this precision.

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2013