DQP#

purpose#

The dqp package uses a dual gradient-projection method to solve a given stricly-convex quadratic program. The aim is to minimize the quadratic objective function

\[q(x) = f + g^T x + \frac{1}{2} x^T H x,\]
or the shifted-least-distance objective function
\[s(x) = f + g^T x + \frac{1}{2} \sum_{j=1}^n w_j^2 (x_j - x_j^0)^2,\]
subject to the general linear constraints and simple bounds
\[c_l \leq A x \leq c_u \;\;\mbox{and} \;\; x_l \leq x \leq x_u,\]
where \(H\) and \(A\) are, respectively, given \(n\) by \(n\) symmetric postive-definite and \(m\) by \(n\) matrices, \(g\), \(w\) and \(x^0\) are vectors, \(f\) is a scalar, and any of the components of the vectors \(c_l\), \(c_u\), \(x_l\) or \(x_u\) may be infinite. The method offers the choice of direct and iterative solution of the key regularization subproblems, and is most suitable for problems involving a large number of unknowns \(x\).

See Section 4 of $GALAHAD/doc/dqp.pdf for additional details.

terminology#

Any required solution \(x\) necessarily satisfies the primal optimality conditions

\[A x = c\]
and
\[c_l \leq c \leq c_u, \;\; x_l \leq x \leq x_u,\]
the dual optimality conditions
\[H x + g = A^{T} y + z,\;\; y = y_l + y_u \;\;\mbox{and}\;\; z = z_l + z_u,\]
and
\[y_l \geq 0, \;\; y_u \leq 0, \;\; z_l \geq 0 \;\;\mbox{and}\;\; z_u \leq 0,\]
and the complementary slackness conditions
\[( A x - c_l )^{T} y_l = 0,\;\; ( A x - c_u )^{T} y_u = 0,\;\; (x -x_l )^{T} z_l = 0 \;\;\mbox{and}\;\;(x -x_u )^{T} z_u = 0,\]
where the vectors \(y\) and \(z\) are known as the Lagrange multipliers for the general linear constraints, and the dual variables for the bounds, respectively, and where the vector inequalities hold component-wise.

In the shifted-least-distance case, \(g\) is shifted by \(-W^2 x^0\), and \(H = W^2\), where \(W\) is the diagonal matrix whose entries are the \(w_j\).

method#

Dual gradient-projection methods solve the quadratic programmimg problem by instead solving the dual quadratic program

\[\begin{split}\begin{array}{ll}\mbox{minimize}\;\; q^D(y^{l}, y^{u}, z^{l}, z^{u}) = & \!\!\! \frac{1}{2} [ ( y^{l0} + y^{u} )^T A + ( z^{l} + z^{u} ]^T ) H^{-1} [ A^T ( y^{l} + y^{u} ) + z^{l} + z^{u} ] \\ & - [ ( y^{l} + y^{u} )^T A + ( z^{l} + z^{u} ]^T ) H^{-1} g - ( c^{l T} y^{l} + c^{u T} y^{u} + x^{l T} z^{l} + x^{u T} z^{u}) \\ \mbox{subject to} & ( y^{l}, z^{l} ) \geq 0 \;\;\mbox{and} \;\; (y^{u}, z^{u}) \leq 0,\end{array}\end{split}\]
and then recovering the required solution from the linear system
\[H x = - g + A^T ( y^{l} + y^{u} ) + z^{l} + z^{u}.\]
The dual problem is solved by an accelerated gradient-projection method comprising of alternating phases in which (i) the current projected dual gradient is traced downhill (the ‘arc search’) as far as possible and (ii) the dual variables that are currently on their bounds are temporarily fixed and the unconstrained minimizer of \(q^D(y^{l}, y^{u}, z^{l}, z^{u})\) with respect to the remaining variables is sought; the minimizer in the second phase may itself need to be projected back into the dual feasible region (either using a brute-force backtrack or a second arc search).

Both phases require the solution of sparse systems of symmetric linear equations, and these are handled by the matrix factorization package SBLS or the conjugate-gradient package GLTR. The systems are commonly singular, and this leads to a requirement to find the Fredholm Alternative for the given matrix and its right-hand side. In the non-singular case, there is an option to update existing factorizations using the “Schur-complement” approach given by the package SCU.

Optionally, the problem may be pre-processed temporarily to eliminate dependent constraints using the package FDC. This may improve the performance of the subsequent iteration.

reference#

The basic algorithm is described in

N. I. M. Gould and D. P. Robinson, ``A dual gradient-projection method for large-scale strictly-convex quadratic problems’’, Computational Optimization and Applications 67(1) (2017) 1-38.

matrix storage#

The unsymmetric \(m\) by \(n\) matrix \(A\) may be presented and stored in a variety of convenient input formats.

Dense storage format: The matrix \(A\) is stored as a compact dense matrix by rows, that is, the values of the entries of each row in turn are stored in order within an appropriate real one-dimensional array. In this case, component \(n \ast i + j\) of the storage array A_val will hold the value \(A_{ij}\) for \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\). The string A_type = ‘dense’ should be specified.

Dense by columns storage format: The matrix \(A\) is stored as a compact dense matrix by columns, that is, the values of the entries of each column in turn are stored in order within an appropriate real one-dimensional array. In this case, component \(m \ast j + i\) of the storage array A_val will hold the value \(A_{ij}\) for \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\). The string A_type = ‘dense_by_columns’ should be specified.

Sparse co-ordinate storage format: Only the nonzero entries of the matrices are stored. For the \(l\)-th entry, \(0 \leq l \leq ne-1\), of \(A\), its row index i, column index j and value \(A_{ij}\), \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\), are stored as the \(l\)-th components of the integer arrays A_row and A_col and real array A_val, respectively, while the number of nonzeros is recorded as A_ne = \(ne\). The string A_type = ‘coordinate’should be specified.

Sparse row-wise storage format: Again only the nonzero entries are stored, but this time they are ordered so that those in row i appear directly before those in row i+1. For the i-th row of \(A\) the i-th component of the integer array A_ptr holds the position of the first entry in this row, while A_ptr(m) holds the total number of entries. The column indices j, \(0 \leq j \leq n-1\), and values \(A_{ij}\) of the nonzero entries in the i-th row are stored in components l = A_ptr(i), \(\ldots\), A_ptr(i+1)-1, \(0 \leq i \leq m-1,\) of the integer array A_col, and real array A_val, respectively. For sparse matrices, this scheme almost always requires less storage than its predecessor. The string A_type = ‘sparse_by_rows’ should be specified.

Sparse column-wise storage format: Once again only the nonzero entries are stored, but this time they are ordered so that those in column j appear directly before those in column j+1. For the j-th column of \(A\) the j-th component of the integer array A_ptr holds the position of the first entry in this column, while A_ptr(n) holds the total number of entries. The row indices i, \(0 \leq i \leq m-1\), and values \(A_{ij}\) of the nonzero entries in the j-th columns are stored in components l = A_ptr(j), \(\ldots\), A_ptr(j+1)-1, \(0 \leq j \leq n-1\), of the integer array A_row, and real array A_val, respectively. As before, for sparse matrices, this scheme almost always requires less storage than the co-ordinate format. The string A_type = ‘sparse_by_columns’ should be specified.

The symmetric \(n\) by \(n\) matrix \(H\) may also be presented and stored in a variety of formats. But crucially symmetry is exploited by only storing values from the lower triangular part (i.e, those entries that lie on or below the leading diagonal).

Dense storage format: The matrix \(H\) is stored as a compact dense matrix by rows, that is, the values of the entries of each row in turn are stored in order within an appropriate real one-dimensional array. Since \(H\) is symmetric, only the lower triangular part (that is the part \(H_{ij}\) for \(0 \leq j \leq i \leq n-1\)) need be held. In this case the lower triangle should be stored by rows, that is component \(i * i / 2 + j\) of the storage array H_val will hold the value \(H_{ij}\) (and, by symmetry, \(H_{ji}\)) for \(0 \leq j \leq i \leq n-1\). The string H_type = ‘dense’ should be specified.

Sparse co-ordinate storage format: Only the nonzero entries of the matrices are stored. For the \(l\)-th entry, \(0 \leq l \leq ne-1\), of \(H\), its row index i, column index j and value \(H_{ij}\), \(0 \leq j \leq i \leq n-1\), are stored as the \(l\)-th components of the integer arrays H_row and H_col and real array H_val, respectively, while the number of nonzeros is recorded as H_ne = \(ne\). Note that only the entries in the lower triangle should be stored. The string H_type = ‘coordinate’ should be specified.

Sparse row-wise storage format: Again only the nonzero entries are stored, but this time they are ordered so that those in row i appear directly before those in row i+1. For the i-th row of \(H\) the i-th component of the integer array H_ptr holds the position of the first entry in this row, while H_ptr(n) holds the total number of entries. The column indices j, \(0 \leq j \leq i\), and values \(H_{ij}\) of the entries in the i-th row are stored in components l = H_ptr(i), …, H_ptr(i+1)-1, \(0 \leq i \leq n-1,\) of the integer array H_col, and real array H_val, respectively. Note that as before only the entries in the lower triangle should be stored. For sparse matrices, this scheme almost always requires less storage than its predecessor. The string H_type = ‘sparse_by_rows’ should be specified.

Diagonal storage format: If \(H\) is diagonal (i.e., \(H_{ij} = 0\) for all \(0 \leq i \neq j \leq n-1\)) only the diagonals entries \(H_{ii}\), \(0 \leq i \leq n-1\) need be stored, and the first n components of the array H_val may be used for the purpose. The string H_type = ‘diagonal’ should be specified.

Multiples of the identity storage format: If \(H\) is a multiple of the identity matrix, (i.e., \(H = \alpha I\) where \(I\) is the n by n identity matrix and \(\alpha\) is a scalar), it suffices to store \(\alpha\) as the first component of H_val. The string H_type = ‘scaled_identity’ should be specified.

The identity matrix format: If \(H\) is the identity matrix, no values need be stored. The string H_type = ‘identity’ should be specified.

The zero matrix format: The same is true if \(H\) is the zero matrix, but now the string H_type = ‘zero’ or ‘none’ should be specified.

functions#

dqp.initialize()#

Set default option values and initialize private data

Returns:

optionsdict
dictionary containing default control options:
errorint

error and warning diagnostics occur on stream error.

outint

general output occurs on stream out.

print_levelint

the level of output required is specified by print_level. Possible values are

  • <=0

    gives no output,

  • 1

    gives a one-line summary for every iteration.

  • 2

    gives a summary of the inner iteration for each iteration.

  • >=3

    gives increasingly verbose (debugging) output.

start_printint

any printing will start on this iteration.

stop_printint

any printing will stop on this iteration.

print_gapint

printing will only occur every print_gap iterations.

dual_starting_pointint

which starting point should be used for the dual problem. Possible values are

  • -1

    user supplied comparing primal vs dual variables.

  • 0

    user supplied.

  • 1

    minimize linearized dual.

  • 2

    minimize simplified quadratic dual.

  • 3

    all free (= all active primal costraints).

  • 4

    all fixed on bounds (= no active primal costraints).

maxitint

at most maxit inner iterations are allowed.

max_scint

the maximum permitted size of the Schur complement before a refactorization is performed (used in the case where there is no Fredholm Alternative, 0 = refactor every iteration).

cauchy_onlyint

a subspace step will only be taken when the current Cauchy step has changed no more than than cauchy_only active constraints; the subspace step will always be taken if cauchy_only < 0.

arc_search_maxitint

how many iterations are allowed per arc search (-ve = as many as require.

cg_maxitint

how many CG iterations to perform per DQP iteration (-ve reverts to n+1).

explore_optimal_subspaceint

once a potentially optimal subspace has been found, investigate it

  • 0

    as per an ordinary subspace.

  • 1

    by increasing the maximum number of allowed CG iterations

  • 2

    by switching to a direct method.

restore_problemint

indicate whether and how much of the input problem should be restored on output. Possible values are

  • 0

    nothing restored.

  • 1

    scalar and vector parameters.

  • 2

    all parameters.

sif_file_deviceint

specifies the unit number to write generated SIF file describing the current problem.

qplib_file_deviceint

specifies the unit number to write generated QPLIB file describing the current problem.

rhofloat

the penalty weight, rho. The general constraints are not enforced explicitly, but instead included in the objective as a penalty term weighted by rho when rho > 0. If rho <= 0, the general constraints are explicit (that is, there is no penalty term in the objective function).

infinityfloat

any bound larger than infinity in modulus will be regarded as infinite.

stop_abs_pfloat

the required absolute and relative accuracies for the primal infeasibilies.

stop_rel_pfloat

see stop_abs_p.

stop_abs_dfloat

the required absolute and relative accuracies for the dual infeasibility.

stop_rel_dfloat

see stop_abs_d.

stop_abs_cfloat

the required absolute and relative accuracies for the complementarity.

stop_rel_cfloat

see stop_abs_c.

stop_cg_relativefloat

the CG iteration will be stopped as soon as the current norm of the preconditioned gradient is smaller than max( stop_cg_relative * initial preconditioned gradient, stop_cg_absolute ).

stop_cg_absolutefloat

see stop_cg_relative.

cg_zero_curvaturefloat

threshold below which curvature is regarded as zero if CG is used.

max_growthfloat

maximum growth factor allowed without a refactorization.

identical_bounds_tolfloat

any pair of constraint bounds (c_l,c_u) or (x_l,x_u) that are closer than identical_bounds_tol will be reset to the average of their values.

cpu_time_limitfloat

the maximum CPU time allowed (-ve means infinite).

clock_time_limitfloat

the maximum elapsed clock time allowed (-ve means infinite).

initial_perturbationfloat

the initial penalty weight (for DLP only).

perturbation_reductionfloat

the penalty weight reduction factor (for DLP only).

final_perturbationfloat

the final penalty weight (for DLP only).

factor_optimal_matrixbool

are the factors of the optimal augmented matrix required? (for DLP only).

remove_dependenciesbool

the equality constraints will be preprocessed to remove any linear dependencies if True.

treat_zero_bounds_as_generalbool

any problem bound with the value zero will be treated as if it were a general value if True.

exact_arc_searchbool

if exact_arc_search is True, an exact piecewise arc search will be performed. Otherwise an ineaxt search using a backtracing Armijo strategy will be employed.

subspace_directbool

if subspace_direct is True, the subspace step will be calculated using a direct (factorization) method, while if it is False, an iterative (conjugate-gradient) method will be used.

subspace_alternatebool

if subspace_alternate is True, the subspace step will alternate between a direct (factorization) method and an iterative (GLTR conjugate-gradient) method. This will override subspace_direct.

subspace_arc_searchbool

if subspace_arc_search is True, a piecewise arc search will be performed along the subspace step. Otherwise the search will stop at the firstconstraint encountered.

space_criticalbool

if space_critical is True, every effort will be made to use as little space as possible. This may result in longer computation time.

deallocate_error_fatalbool

if deallocate_error_fatal is True, any array/pointer deallocation error will terminate execution. Otherwise, computation will continue.

generate_sif_filebool

if generate_sif_file is True, a SIF file describing the current problem is to be generated.

generate_qplib_filebool

if generate_qplib_file is True, a QPLIB file describing the current problem is to be generated.

symmetric_linear_solverstr

indefinite linear equation solver set in symmetric_linear_solver.

definite_linear_solverstr

definite linear equation solver.

unsymmetric_linear_solverstr

unsymmetric linear equation solver.

sif_file_namestr

name of generated SIF file containing input problem.

qplib_file_namestr

name of generated QPLIB file containing input problem.

prefixstr

all output lines will be prefixed by the string contained in quotes within prefix, e.g. ‘word’ (note the qutoes) will result in the prefix word.

fdc_optionsdict

default control options for FDC (see fdc.initialize).

sls_optionsdict

default control options for SLS (see sls.initialize).

sbls_optionsdict

default control options for SBLS (see sbls.initialize).

gltr_optionsdict

default control options for GLTR (see gltr.initialize).

dqp.load(n, m, H_type, H_ne, H_row, H_col, H_ptr, A_type, A_ne, A_row, A_col, A_ptr, options=None)#

Import problem data into internal storage prior to solution.

Parameters:

nint

holds the number of variables.

mint

holds the number of constraints.

H_typestring

specifies the symmetric storage scheme used for the Hessian \(H\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, ‘diagonal’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’; lower or upper case variants are allowed.

H_neint

holds the number of entries in the lower triangular part of \(H\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes.

H_rowndarray(H_ne)

holds the row indices of the lower triangular part of \(H\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes, and in this case can be None.

H_colndarray(H_ne)

holds the column indices of the lower triangular part of \(H\) in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the other storage schemes are used, and in this case can be None.

H_ptrndarray(n+1)

holds the starting position of each row of the lower triangular part of \(H\), as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the other schemes are used, and in this case can be None.

A_typestring

specifies the unsymmetric storage scheme used for the constraints Jacobian \(A\). It should be one of ‘coordinate’, ‘sparse_by_rows’ or ‘dense’; lower or upper case variants are allowed.

A_neint

holds the number of entries in \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other two schemes.

A_rowndarray(A_ne)

holds the row indices of \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other two schemes, and in this case can be None.

A_colndarray(A_ne)

holds the column indices of \(A\) in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the dense storage scheme is used, and in this case can be None.

A_ptrndarray(m+1)

holds the starting position of each row of \(A\), as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the other schemes are used, and in this case can be None.

optionsdict, optional

dictionary of control options (see dqp.initialize).

dqp.solve_qp(n, m, f, g, H_ne, H_val, A_ne, A_val, c_l, c_u, x_l, x_u, x, y, z)#

Find a solution to the strictly convex quadratic program involving the quadratic objective function \(q(x)\).

Parameters:

nint

holds the number of variables.

mint

holds the number of residuals.

ffloat

holds the constant term \(f\) in the objective function.

gndarray(n)

holds the values of the linear term \(g\) in the objective function.

H_neint

holds the number of entries in the lower triangular part of the Hessian \(H\).

H_valndarray(H_ne)

holds the values of the nonzeros in the lower triangle of the Hessian \(H\) in the same order as specified in the sparsity pattern in dqp.load.

A_neint

holds the number of entries in the constraint Jacobian \(A\).

A_valndarray(A_ne)

holds the values of the nonzeros in the constraint Jacobian \(A\) in the same order as specified in the sparsity pattern in dqp.load.

c_lndarray(m)

holds the values of the lower bounds \(c_l\) on the constraints The lower bound on any component of \(A x\) that is unbounded from below should be set no larger than minus options.infinity.

c_undarray(m)

holds the values of the upper bounds \(c_l\) on the constraints The upper bound on any component of \(A x\) that is unbounded from above should be set no smaller than options.infinity.

x_lndarray(n)

holds the values of the lower bounds \(x_l\) on the variables. The lower bound on any component of \(x\) that is unbounded from below should be set no larger than minus options.infinity.

x_undarray(n)

holds the values of the upper bounds \(x_l\) on the variables. The upper bound on any component of \(x\) that is unbounded from above should be set no smaller than options.infinity.

xndarray(n)

holds the initial estimate of the minimizer \(x\), if known. This is not crucial, and if no suitable value is known, then any value, such as \(x=0\), suffices and will be adjusted accordingly.

yndarray(m)

holds the initial estimate of the Lagrange multipliers \(y\) associated with the general constraints, if known. This is not crucial, and if no suitable value is known, then any value, such as \(y=0\), suffices and will be adjusted accordingly.

zndarray(n)

holds the initial estimate of the dual variables \(z\) associated with the simple bound constraints, if known. This is not crucial, and if no suitable value is known, then any value, such as \(z=0\), suffices and will be adjusted accordingly.

Returns:

xndarray(n)

holds the values of the approximate minimizer \(x\) after a successful call.

cndarray(m)

holds the values of the residuals \(c(x) = Ax\).

yndarray(m)

holds the values of the Lagrange multipliers associated with the general linear constraints.

zndarray(n)

holds the values of the dual variables associated with the simple bound constraints.

c_statndarray(m)

holds the return status for each constraint. The i-th component will be negative if the value of the \(i\)-th constraint \((Ax)_i\)) lies on its lower bound, positive if it lies on its upper bound, and zero if it lies between bounds.

x_statndarray(n)

holds the return status for each variable. The i-th component will be negative if the \(i\)-th variable lies on its lower bound, positive if it lies on its upper bound, and zero if it lies between bounds.

dqp.solve_sldqp(n, m, f, g, w, x0, A_ne, A_val, c_l, c_u, x_l, x_u, x, y, z)#

Find a solution to the quadratic program involving the shifted least-distance objective function \(s(x)\).

Parameters:

nint

holds the number of variables.

mint

holds the number of residuals.

ffloat

holds the constant term \(f\) in the objective function.

gndarray(n)

holds the values of the linear term \(g\) in the objective function.

wndarray(n)

holds the values of the weights \(w\) in the objective function.

x0ndarray(n)

holds the values of the shifts \(x^0\) in the objective function.

A_neint

holds the number of entries in the constraint Jacobian \(A\).

A_valndarray(A_ne)

holds the values of the nonzeros in the constraint Jacobian \(A\) in the same order as specified in the sparsity pattern in cqp.load.

c_lndarray(m)

holds the values of the lower bounds \(c_l\) on the constraints The lower bound on any component of \(A x\) that is unbounded from below should be set no larger than minus options.infinity.

c_undarray(m)

holds the values of the upper bounds \(c_l\) on the constraints The upper bound on any component of \(A x\) that is unbounded from above should be set no smaller than options.infinity.

x_lndarray(n)

holds the values of the lower bounds \(x_l\) on the variables. The lower bound on any component of \(x\) that is unbounded from below should be set no larger than minus options.infinity.

x_undarray(n)

holds the values of the upper bounds \(x_l\) on the variables. The upper bound on any component of \(x\) that is unbounded from above should be set no smaller than options.infinity.

xndarray(n)

holds the initial estimate of the minimizer \(x\), if known. This is not crucial, and if no suitable value is known, then any value, such as \(x=0\), suffices and will be adjusted accordingly.

yndarray(m)

holds the initial estimate of the Lagrange multipliers \(y\) associated with the general constraints, if known. This is not crucial, and if no suitable value is known, then any value, such as \(y=0\), suffices and will be adjusted accordingly.

zndarray(n)

holds the initial estimate of the dual variables \(z\) associated with the simple bound constraints, if known. This is not crucial, and if no suitable value is known, then any value, such as \(z=0\), suffices and will be adjusted accordingly.

Returns:

xndarray(n)

holds the values of the approximate minimizer \(x\) after a successful call.

cndarray(m)

holds the values of the residuals \(c(x) = Ax\).

yndarray(m)

holds the values of the Lagrange multipliers associated with the general linear constraints.

zndarray(n)

holds the values of the dual variables associated with the simple bound constraints.

c_statndarray(m)

holds the return status for each constraint. The i-th component will be negative if the value of the \(i\)-th constraint \((Ax)_i\)) lies on its lower bound, positive if it lies on its upper bound, and zero if it lies between bounds.

xstatndarray(n)

holds the return status for each variable. The i-th component will be negative if the \(i\)-th variable lies on its lower bound, positive if it lies on its upper bound, and zero if it lies between bounds.

[optional] dqp.information()

Provide optional output information

Returns:

informdict
dictionary containing output information:
statusint

return status. Possible values are:

  • 0

    The run was successful.

  • -1

    An allocation error occurred. A message indicating the offending array is written on unit options[‘error’], and the returned allocation status and a string containing the name of the offending array are held in inform[‘alloc_status’] and inform[‘bad_alloc’] respectively.

  • -2

    A deallocation error occurred. A message indicating the offending array is written on unit options[‘error’] and the returned allocation status and a string containing the name of the offending array are held in inform[‘alloc_status’] and inform[‘bad_alloc’] respectively.

  • -3

    The restriction n > 0 or m > 0 or requirement that type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, ‘diagonal’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’ has been violated.

  • -4

    The bound constraints are inconsistent.

  • -5

    The constraints appear to have no feasible point.

  • -9

    The analysis phase of the factorization failed; the return status from the factorization package is given by inform[‘factor_status’].

  • -10

    The factorization failed; the return status from the factorization package is given by inform[‘factor_status’].

  • -16

    The problem is so ill-conditioned that further progress is impossible.

  • -18

    Too many iterations have been performed. This may happen if options[‘maxit’] is too small, but may also be symptomatic of a badly scaled problem.

  • -19

    The CPU time limit has been reached. This may happen if options[‘cpu_time_limit’] is too small, but may also be symptomatic of a badly scaled problem.

  • -20

    The Hessian \(H\) appears not to be positive definite.

  • -23

    An entry from the strict upper triangle of \(H\) has been specified.

alloc_statusint

the status of the last attempted allocation/deallocation.

bad_allocstr

the name of the array for which an allocation/deallocation error occurred.

iterint

the total number of iterations required.

cg_iterint

the total number of iterations required.

factorization_statusint

the return status from the factorization.

factorization_integerlong

the total integer workspace required for the factorization.

factorization_reallong

the total real workspace required for the factorization.

nfactsint

the total number of factorizations performed.

threadsint

the number of threads used.

objfloat

the value of the objective function at the best estimate of the solution determined by DQP_solve.

primal_infeasibilityfloat

the value of the primal infeasibility.

dual_infeasibilityfloat

the value of the dual infeasibility.

complementary_slacknessfloat

the value of the complementary slackness.

non_negligible_pivotfloat

the smallest pivot that was not judged to be zero when detecting linearly dependent constraints.

feasiblebool

is the returned “solution” feasible?.

checkpointsIterndarray(16)

checkpointsIter(i) records the iteration at which the criticality measures first fall below \(10^{-i-1}, i = 0, \ldots 15\) (where -1 means not achieved).

checkpointsTimendarray(16)

checkpointsTime(i) records the CPU time at which the criticality measures first fall below \(10^{-i-1}, i = 0, \ldots 15\) (where -1 means not achieved).

timedict
dictionary containing timing information:
totalfloat

the total CPU time spent in the package.

preprocessfloat

the CPU time spent preprocessing the problem.

find_dependentfloat

the CPU time spent detecting linear dependencies.

analysefloat

the CPU time spent analysing the required matrices prior to factorization.

factorizefloat

the CPU time spent factorizing the required matrices.

solvefloat

the CPU time spent computing the search direction.

searchfloat

the CPU time spent in the linesearch.

clock_totalfloat

the total clock time spent in the package.

clock_preprocessfloat

the clock time spent preprocessing the problem.

clock_find_dependentfloat

the clock time spent detecting linear dependencies.

clock_analysefloat

the clock time spent analysing the required matrices prior to factorization.

clock_factorizefloat

the clock time spent factorizing the required matrices.

clock_solvefloat

the clock time spent computing the search direction.

clock_searchfloat

the clock time spent in the linesearch.

fdc_informdict

inform parameters for FDC (see fdc.information).

sls_informdict

inform parameters for SLS (see sls.information).

sbls_informdict

inform parameters for SBLS (see sbls.information).

gltr_informdict

return information from GLTR (see gltr.information).

scu_statusint

status value for SCU (see scu.status).

scu_informdict

inform parameters for SCU (see scu.information).

rpd_informdict

inform parameters for RPD (see rpd.information).

dqp.terminate()#

Deallocate all internal private storage.

example code#

from galahad import dqp
import numpy as np
np.set_printoptions(precision=4,suppress=True,floatmode='fixed')
print("\n** python test: dqp")

# set parameters
n = 3
m = 2
infinity = float("inf")

#  describe objective function

f = 1.0
g = np.array([0.0,2.0,0.0])
H_type = 'coordinate'
H_ne = 4
H_row = np.array([0,1,2,2])
H_col = np.array([0,1,1,2])
H_ptr = None
H_val = np.array([1.0,2.0,1.0,3.0])

#  describe constraints

A_type = 'coordinate'
A_ne = 4
A_row = np.array([0,0,1,1])
A_col = np.array([0,1,1,2])
A_ptr = None
A_val = np.array([2.0,1.0,1.0,1.0])
c_l = np.array([1.0,2.0])
c_u = np.array([2.0,2.0])
x_l = np.array([-1.0,-infinity,-infinity])
x_u = np.array([1.0,infinity,2.0])

# allocate internal data and set default options
options = dqp.initialize()

# set some non-default options
options['print_level'] = 0
options['symmetric_linear_solver'] = 'sytr '
options['definite_linear_solver'] = 'sytr '
options['sbls_options']['symmetric_linear_solver'] = 'sytr '
options['sbls_options']['definite_linear_solver'] = 'sytr '
#print("options:", options)

# load data (and optionally non-default options)
dqp.load(n, m, H_type, H_ne, H_row, H_col, H_ptr,
         A_type, A_ne, A_row, A_col, A_ptr, options)

#  provide starting values (not crucial)

x = np.array([0.0,0.0,0.0])
y = np.array([0.0,0.0])
z = np.array([0.0,0.0,0.0])

# find optimum of qp
print("\n 1st problem: solve qp")
x, c, y, z, x_stat, c_stat \
  = dqp.solve_qp(n, m, f, g, H_ne, H_val,
                 A_ne, A_val, c_l, c_u, x_l, x_u, x, y, z)
print(" x:",x)
print(" c:",c)
print(" y:",y)
print(" z:",z)
print(" x_stat:",x_stat)
print(" c_stat:",c_stat)

# get information
inform = dqp.information()
print(" f: %.4f" % inform['obj'])

# deallocate internal data

dqp.terminate()

#  describe shifted-least-distance qp

w = np.array([1.0,1.0,1.0])
x0 = np.array([1.0,1.0,1.0])
H_type = 'shifted_least_distance'

# allocate internal data
dqp.initialize()

# load data (and optionally non-default options)
dqp.load(n, m, H_type, H_ne, H_row, H_col, H_ptr,
         A_type, A_ne, A_row, A_col, A_ptr, options)

#  provide starting values (not crucial)

x = np.array([0.0,0.0,0.0])
y = np.array([0.0,0.0])
z = np.array([0.0,0.0,0.0])

# find optimum of sldqp
print("\n 2nd problem: solve sldqp")
x, c, y, z, x_stat, c_stat \
  = dqp.solve_sldqp(n, m, f, g, w, x0,
                    A_ne, A_val, c_l, c_u, x_l, x_u, x, y, z)
print(" x:",x)
print(" c:",c)
print(" y:",y)
print(" z:",z)
print(" x_stat:",x_stat)
print(" c_stat:",c_stat)

# get information
inform = dqp.information()
print(" f: %.4f" % inform['obj'])
print('** dqp exit status:', inform['status'])

# deallocate internal data

dqp.terminate()

This example code is available in $GALAHAD/src/dqp/Python/test_dqp.py .