GALAHAD CQP package#
purpose#
The cqp
package uses a primal-dual interior-point method to solve a
given convex quadratic program.
The aim is to minimize the quadratic objective function
See Section 4 of $GALAHAD/doc/cqp.pdf for additional details.
terminology#
Any required solution \(x\) necessarily satisfies the primal optimality conditions
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#
Primal-dual interior point methods iterate towards a point that satisfies these optimality conditions by ultimately aiming to satisfy (1a), (2a) and (3), while ensuring that (1b) and (2b) are satisfied as strict inequalities at each stage. Appropriate norms of the amounts by which (1a), (2a) and (3) fail to be satisfied are known as the primal and dual infeasibility, and the violation of complementary slackness, respectively. The fact that (1b) and (2b) are satisfied as strict inequalities gives such methods their other title, namely interior-point methods.
The method aims at each stage to reduce the overall violation of (1a), (2a) and (3), rather than reducing each of the terms individually. Given an estimate \(v = (x, \; c, \; y, \; y^{l}, \; y^{u}, \; z, \; z^{l}, \; z^{u})\) of the primal-dual variables, a correction \(\Delta v = \Delta (x, \; c, \; y, \; y^{l}, \; y^{u} ,\;z,\;z^{l} ,\;z^{u} )\) is obtained by solving a suitable linear system of Newton equations for the nonlinear systems (1a), (2a) and a parameterized ``residual trajectory’’ perturbation of (3); residual trajectories proposed by Zhang (1994) and Zhao and Sun (1999) are possibilities. An improved estimate \(v + \alpha \Delta v\) is then used, where the step-size \(\alpha\) is chosen as close to 1.0 as possible while ensuring both that (1b) and (2b) continue to hold and that the individual components which make up the complementary slackness (3) do not deviate too significantly from their average value. The parameter that controls the perturbation of (3) is ultimately driven to zero.
If the algorithm believes that it is close to the solution, it may take a speculative ``pounce’’ extrapolation, based on an estimate of the ultimate active set, to avoid further costly iterations. If the pounce is unsuccessful, the iteration continues, but further pounces may be attempted later.
The Newton equations are solved by applying the matrix factorization
package SBLS
, but there are options
to factorize the matrix as a whole (the so-called “augmented system”
approach), to perform a block elimination first (the “Schur-complement”
approach), or to let the method itself decide which of the two
previous options is more appropriate.
The “Schur-complement” approach is usually to be preferred when all the
weights are nonzero or when every variable is bounded (at least one side),
but may be inefficient if any of the columns of \(A\) is too dense.
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.
references#
The basic algorithm is a generalisation of those of
Y. Zhang, ``On the convergence of a class of infeasible interior-point methods for the horizontal linear complementarity problem’’. SIAM J. Optimization 4(1) (1994) 208-227,
and
with many enhancements described by
N. I. M. Gould, D. Orban and D. P. Robinson, ``Trajectory-following methods for large-scale degenerate convex quadratic programming’’, Mathematical Programming Computation 5(2) (2013) 113-142.
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.
introduction to function calls#
To solve a given problem, functions from the cqp package must be called in the following order:
cqp_initialize - provide default control parameters and set up initial data structures
cqp_read_specfile (optional) - override control values by reading replacement values from a file
cqp_import - set up problem data structures and fixed values
cqp_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved
solve the problem by calling one of
cqp_solve_qp - solve the quadratic program
cqp_solve_sldqp - solve the shifted least-distance problem
cqp_information (optional) - recover information about the solution and solution process
cqp_terminate - deallocate data structures
See the examples section for illustrations of use.
callable functions#
overview of functions provided#
// typedefs typedef float spc_; typedef double rpc_; typedef int ipc_; // structs struct cqp_control_type; struct cqp_inform_type; struct cqp_time_type; // function calls void cqp_initialize(void **data, struct cqp_control_type* control, ipc_ *status); void cqp_read_specfile(struct cqp_control_type* control, const char specfile[]); void cqp_import( struct cqp_control_type* control, void **data, ipc_ *status, ipc_ n, ipc_ m, const char H_type[], ipc_ H_ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[], const char A_type[], ipc_ A_ne, const ipc_ A_row[], const ipc_ A_col[], const ipc_ A_ptr[] ); void cqp_reset_control( struct cqp_control_type* control, void **data, ipc_ *status ); void cqp_solve_qp( void **data, ipc_ *status, ipc_ n, ipc_ m, ipc_ h_ne, const rpc_ H_val[], const rpc_ g[], const rpc_ f, ipc_ a_ne, const rpc_ A_val[], const rpc_ c_l[], const rpc_ c_u[], const rpc_ x_l[], const rpc_ x_u[], rpc_ x[], rpc_ c[], rpc_ y[], rpc_ z[], ipc_ x_stat[], ipc_ c_stat[] ); void cqp_solve_sldqp( void **data, ipc_ *status, ipc_ n, ipc_ m, const rpc_ w[], const rpc_ x0[], const rpc_ g[], const rpc_ f, ipc_ a_ne, const rpc_ A_val[], const rpc_ c_l[], const rpc_ c_u[], const rpc_ x_l[], const rpc_ x_u[], rpc_ x[], rpc_ c[], rpc_ y[], rpc_ z[], ipc_ x_stat[], ipc_ c_stat[] ); void cqp_information(void **data, struct cqp_inform_type* inform, ipc_ *status); void cqp_terminate( void **data, struct cqp_control_type* control, struct cqp_inform_type* inform );
typedefs#
typedef float spc_
spc_
is real single precision
typedef double rpc_
rpc_
is the real working precision used, but may be changed to float
by
defining the preprocessor variable REAL_32
or (if supported) to
__real128
using the variable REAL_128
.
typedef int ipc_
ipc_
is the default integer word length used, but may be changed to
int64_t
by defining the preprocessor variable INTEGER_64
.
function calls#
void cqp_initialize(void **data, struct cqp_control_type* control, ipc_ *status)
Set default control values and initialize private data
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see cqp_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void cqp_read_specfile(struct cqp_control_type* control, const char specfile[])
Read the content of a specification file, and assign values associated with given keywords to the corresponding control parameters. An in-depth discussion of specification files is available, and a detailed list of keywords with associated default values is provided in $GALAHAD/src/cqp/CQP.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/cqp.pdf for a list of how these keywords relate to the components of the control structure.
Parameters:
control |
is a struct containing control information (see cqp_control_type) |
specfile |
is a character string containing the name of the specification file |
void cqp_import( struct cqp_control_type* control, void **data, ipc_ *status, ipc_ n, ipc_ m, const char H_type[], ipc_ H_ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[], const char A_type[], ipc_ A_ne, const ipc_ A_row[], const ipc_ A_col[], const ipc_ A_ptr[] )
Import problem data into internal storage prior to solution.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see cqp_control_type) |
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are:
|
n |
is a scalar variable of type ipc_, that holds the number of variables. |
m |
is a scalar variable of type ipc_, that holds the number of general linear constraints. |
H_type |
is a one-dimensional array of type char that 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’, the latter pair if \(H=0\); lower or upper case variants are allowed. |
H_ne |
is a scalar variable of type ipc_, that 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_row |
is a one-dimensional array of size H_ne and type ipc_, that 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 three schemes, and in this case can be NULL. |
H_col |
is a one-dimensional array of size H_ne and type ipc_, that 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 dense, diagonal or (scaled) identity storage schemes are used, and in this case can be NULL. |
H_ptr |
is a one-dimensional array of size n+1 and type ipc_, that 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 NULL. |
A_type |
is a one-dimensional array of type char that specifies the unsymmetric storage scheme used for the constraint Jacobian, \(A\). It should be one of ‘coordinate’, ‘sparse_by_rows’ or ‘dense; lower or upper case variants are allowed. |
A_ne |
is a scalar variable of type ipc_, that holds the number of entries in \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. |
A_row |
is a one-dimensional array of size A_ne and type ipc_, that holds the row indices of \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes, and in this case can be NULL. |
A_col |
is a one-dimensional array of size A_ne and type ipc_, that 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 or diagonal storage schemes are used, and in this case can be NULL. |
A_ptr |
is a one-dimensional array of size n+1 and type ipc_, that 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 NULL. |
void cqp_reset_control( struct cqp_control_type* control, void **data, ipc_ *status )
Reset control parameters after import if required.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see cqp_control_type) |
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are:
|
void cqp_solve_qp( void **data, ipc_ *status, ipc_ n, ipc_ m, ipc_ h_ne, const rpc_ H_val[], const rpc_ g[], const rpc_ f, ipc_ a_ne, const rpc_ A_val[], const rpc_ c_l[], const rpc_ c_u[], const rpc_ x_l[], const rpc_ x_u[], rpc_ x[], rpc_ c[], rpc_ y[], rpc_ z[], ipc_ x_stat[], ipc_ c_stat[] )
Solve the quadratic program when the Hessian \(H\) is available.
Parameters:
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the entry and exit status from the package. Possible exit values are:
|
n |
is a scalar variable of type ipc_, that holds the number of variables |
m |
is a scalar variable of type ipc_, that holds the number of general linear constraints. |
h_ne |
is a scalar variable of type ipc_, that holds the number of entries in the lower triangular part of the Hessian matrix \(H\). |
H_val |
is a one-dimensional array of size h_ne and type rpc_, that holds the values of the entries of the lower triangular part of the Hessian matrix \(H\) in any of the available storage schemes. |
g |
is a one-dimensional array of size n and type rpc_, that holds the linear term \(g\) of the objective function. The j-th component of g, j = 0, … , n-1, contains \(g_j\). |
f |
is a scalar of type rpc_, that holds the constant term \(f\) of the objective function. |
a_ne |
is a scalar variable of type ipc_, that holds the number of entries in the constraint Jacobian matrix \(A\). |
A_val |
is a one-dimensional array of size a_ne and type rpc_, that holds the values of the entries of the constraint Jacobian matrix \(A\) in any of the available storage schemes. |
c_l |
is a one-dimensional array of size m and type rpc_, that holds the lower bounds \(c^l\) on the constraints \(A x\). The i-th component of c_l, i = 0, … , m-1, contains \(c^l_i\). |
c_u |
is a one-dimensional array of size m and type rpc_, that holds the upper bounds \(c^l\) on the constraints \(A x\). The i-th component of c_u, i = 0, … , m-1, contains \(c^u_i\). |
x_l |
is a one-dimensional array of size n and type rpc_, that holds the lower bounds \(x^l\) on the variables \(x\). The j-th component of x_l, j = 0, … , n-1, contains \(x^l_j\). |
x_u |
is a one-dimensional array of size n and type rpc_, that holds the upper bounds \(x^l\) on the variables \(x\). The j-th component of x_u, j = 0, … , n-1, contains \(x^l_j\). |
x |
is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\). |
c |
is a one-dimensional array of size m and type rpc_, that holds the residual \(c(x)\). The i-th component of c, j = 0, … , n-1, contains \(c_j(x)\). |
y |
is a one-dimensional array of size n and type rpc_, that holds the values \(y\) of the Lagrange multipliers for the general linear constraints. The j-th component of y, j = 0, … , n-1, contains \(y_j\). |
z |
is a one-dimensional array of size n and type rpc_, that holds the values \(z\) of the dual variables. The j-th component of z, j = 0, … , n-1, contains \(z_j\). |
x_stat |
is a one-dimensional array of size n and type ipc_, that gives the optimal status of the problem variables. If x_stat(j) is negative, the variable \(x_j\) most likely lies on its lower bound, if it is positive, it lies on its upper bound, and if it is zero, it lies between its bounds. |
c_stat |
is a one-dimensional array of size m and type ipc_, that gives the optimal status of the general linear constraints. If c_stat(i) is negative, the constraint value \(a_i^Tx\) most likely lies on its lower bound, if it is positive, it lies on its upper bound, and if it is zero, it lies between its bounds. |
void cqp_solve_sldqp( void **data, ipc_ *status, ipc_ n, ipc_ m, const rpc_ w[], const rpc_ x0[], const rpc_ g[], const rpc_ f, ipc_ a_ne, const rpc_ A_val[], const rpc_ c_l[], const rpc_ c_u[], const rpc_ x_l[], const rpc_ x_u[], rpc_ x[], rpc_ c[], rpc_ y[], rpc_ z[], ipc_ x_stat[], ipc_ c_stat[] )
Solve the shifted least-distance quadratic program
Parameters:
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the entry and exit status from the package. Possible exit values are:
|
n |
is a scalar variable of type ipc_, that holds the number of variables |
m |
is a scalar variable of type ipc_, that holds the number of general linear constraints. |
w |
is a one-dimensional array of size n and type rpc_, that holds the values of the weights \(w\). |
x0 |
is a one-dimensional array of size n and type rpc_, that holds the values of the shifts \(x^0\). |
g |
is a one-dimensional array of size n and type rpc_, that holds the linear term \(g\) of the objective function. The j-th component of g, j = 0, … , n-1, contains \(g_j\). |
f |
is a scalar of type rpc_, that holds the constant term \(f\) of the objective function. |
a_ne |
is a scalar variable of type ipc_, that holds the number of entries in the constraint Jacobian matrix \(A\). |
A_val |
is a one-dimensional array of size a_ne and type rpc_, that holds the values of the entries of the constraint Jacobian matrix \(A\) in any of the available storage schemes. |
c_l |
is a one-dimensional array of size m and type rpc_, that holds the lower bounds \(c^l\) on the constraints \(A x\). The i-th component of c_l, i = 0, … , m-1, contains \(c^l_i\). |
c_u |
is a one-dimensional array of size m and type rpc_, that holds the upper bounds \(c^l\) on the constraints \(A x\). The i-th component of c_u, i = 0, … , m-1, contains \(c^u_i\). |
x_l |
is a one-dimensional array of size n and type rpc_, that holds the lower bounds \(x^l\) on the variables \(x\). The j-th component of x_l, j = 0, … , n-1, contains \(x^l_j\). |
x_u |
is a one-dimensional array of size n and type rpc_, that holds the upper bounds \(x^l\) on the variables \(x\). The j-th component of x_u, j = 0, … , n-1, contains \(x^l_j\). |
x |
is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\). |
c |
is a one-dimensional array of size m and type rpc_, that holds the residual \(c(x)\). The i-th component of c, i = 0, … , m-1, contains \(c_i(x)\). |
y |
is a one-dimensional array of size n and type rpc_, that holds the values \(y\) of the Lagrange multipliers for the general linear constraints. The j-th component of y, i = 0, … , m-1, contains \(y_i\). |
z |
is a one-dimensional array of size n and type rpc_, that holds the values \(z\) of the dual variables. The j-th component of z, j = 0, … , n-1, contains \(z_j\). |
x_stat |
is a one-dimensional array of size n and type ipc_, that gives the optimal status of the problem variables. If x_stat(j) is negative, the variable \(x_j\) most likely lies on its lower bound, if it is positive, it lies on its upper bound, and if it is zero, it lies between its bounds. |
c_stat |
is a one-dimensional array of size m and type ipc_, that gives the optimal status of the general linear constraints. If c_stat(i) is negative, the constraint value \(a_i^T x\) most likely lies on its lower bound, if it is positive, it lies on its upper bound, and if it is zero, it lies between its bounds. |
void cqp_information(void **data, struct cqp_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see cqp_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void cqp_terminate( void **data, struct cqp_control_type* control, struct cqp_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see cqp_control_type) |
inform |
is a struct containing output information (see cqp_inform_type) |
available structures#
cqp_control_type structure#
#include <galahad_cqp.h> struct cqp_control_type { // components bool f_indexing; ipc_ error; ipc_ out; ipc_ print_level; ipc_ start_print; ipc_ stop_print; ipc_ maxit; ipc_ infeas_max; ipc_ muzero_fixed; ipc_ restore_problem; ipc_ indicator_type; ipc_ arc; ipc_ series_order; ipc_ sif_file_device; ipc_ qplib_file_device; rpc_ infinity; rpc_ stop_abs_p; rpc_ stop_rel_p; rpc_ stop_abs_d; rpc_ stop_rel_d; rpc_ stop_abs_c; rpc_ stop_rel_c; rpc_ perturb_h; rpc_ prfeas; rpc_ dufeas; rpc_ muzero; rpc_ tau; rpc_ gamma_c; rpc_ gamma_f; rpc_ reduce_infeas; rpc_ obj_unbounded; rpc_ potential_unbounded; rpc_ identical_bounds_tol; rpc_ mu_pounce; rpc_ indicator_tol_p; rpc_ indicator_tol_pd; rpc_ indicator_tol_tapia; rpc_ cpu_time_limit; rpc_ clock_time_limit; bool remove_dependencies; bool treat_zero_bounds_as_general; bool treat_separable_as_general; bool just_feasible; bool getdua; bool puiseux; bool every_order; bool feasol; bool balance_initial_complentarity; bool crossover; bool space_critical; bool deallocate_error_fatal; bool generate_sif_file; bool generate_qplib_file; char sif_file_name[31]; char qplib_file_name[31]; char prefix[31]; struct fdc_control_type fdc_control; struct sbls_control_type sbls_control; struct fit_control_type fit_control; struct roots_control_type roots_control; struct cro_control_type cro_control; };
detailed documentation#
control derived type as a C struct
components#
bool f_indexing
use C or Fortran sparse matrix indexing
ipc_ error
error and warning diagnostics occur on stream error
ipc_ out
general output occurs on stream out
ipc_ print_level
the level of output required is specified by print_level
\(\leq\) 0 gives no output,
= 1 gives a one-line summary for every iteration,
= 2 gives a summary of the inner iteration for each iteration,
\(\geq\) 3 gives increasingly verbose (debugging) output
ipc_ start_print
any printing will start on this iteration
ipc_ stop_print
any printing will stop on this iteration
ipc_ maxit
at most maxit inner iterations are allowed
ipc_ infeas_max
the number of iterations for which the overall infeasibility of the problem is not reduced by at least a factor .reduce_infeas before the problem is flagged as infeasible (see reduce_infeas)
ipc_ muzero_fixed
the initial value of the barrier parameter will not be changed for the first muzero_fixed iterations
ipc_ restore_problem
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
ipc_ indicator_type
specifies the type of indicator function used. Possible values are
1 primal indicator: a constraint is active if and only if the distance to its nearest bound \(\leq\).indicator_p_tol
2 primal-dual indicator: a constraint is active if and only if the distance to its nearest bound \(\leq\).indicator_tol_pd \* size of corresponding multiplier
3 primal-dual indicator: a constraint is active if and only if the distance to its nearest bound \(\leq\).indicator_tol_tapia \* distance to same bound at previous iteration
ipc_ arc
which residual trajectory should be used to aim from the current iterate to the solution. Possible values are
1 the Zhang linear residual trajectory
2 the Zhao-Sun quadratic residual trajectory
3 the Zhang arc ultimately switching to the Zhao-Sun residual trajectory
4 the mixed linear-quadratic residual trajectory
5 the Zhang arc ultimately switching to the mixed linear-quadratic residual trajectory
ipc_ series_order
the order of (Taylor/Puiseux) series to fit to the path data
ipc_ sif_file_device
specifies the unit number to write generated SIF file describing the current problem
ipc_ qplib_file_device
specifies the unit number to write generated QPLIB file describing the current problem
rpc_ infinity
any bound larger than infinity in modulus will be regarded as infinite
rpc_ stop_abs_p
the required absolute and relative accuracies for the primal infeasibility
rpc_ stop_rel_p
see stop_abs_p
rpc_ stop_abs_d
the required absolute and relative accuracies for the dual infeasibility
rpc_ stop_rel_d
see stop_abs_d
rpc_ stop_abs_c
the required absolute and relative accuracies for the complementarity
rpc_ stop_rel_c
see stop_abs_c
rpc_ perturb_h
.perturb_h will be added to the Hessian
rpc_ prfeas
initial primal variables will not be closer than .prfeas from their bounds
rpc_ dufeas
initial dual variables will not be closer than .dufeas from their bounds
rpc_ muzero
the initial value of the barrier parameter. If muzero is not positive, it will be reset to an appropriate value
rpc_ tau
the weight attached to primal-dual infeasibility compared to complementa when assessing step acceptance
rpc_ gamma_c
individual complementarities will not be allowed to be smaller than gamma_c times the average value
rpc_ gamma_f
the average complementarity will not be allowed to be smaller than gamma_f times the primal/dual infeasibility
rpc_ reduce_infeas
if the overall infeasibility of the problem is not reduced by at least a factor .reduce_infeas over .infeas_max iterations, the problem is flagged as infeasible (see infeas_max)
rpc_ obj_unbounded
if the objective function value is smaller than obj_unbounded, it will be flagged as unbounded from below.
rpc_ potential_unbounded
if W=0 and the potential function value is smaller than .potential_unbounded \(\ast\) number of one-sided bounds, the analytic center will be flagged as unbounded
rpc_ identical_bounds_tol
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
rpc_ mu_pounce
start terminal extrapolation when mu reaches mu_pounce
rpc_ indicator_tol_p
if .indicator_type = 1, a constraint/bound will be deemed to be active if and only if the distance to its nearest bound \(\leq\).indicator_p_tol
rpc_ indicator_tol_pd
if .indicator_type = 2, a constraint/bound will be deemed to be active if and only if the distance to its nearest bound \(\leq\).indicator_tol_pd \* size of corresponding multiplier
rpc_ indicator_tol_tapia
if .indicator_type = 3, a constraint/bound will be deemed to be active if and only if the distance to its nearest bound \(\leq\).indicator_tol_tapia \* distance to same bound at previous iteration
rpc_ cpu_time_limit
the maximum CPU time allowed (-ve means infinite)
rpc_ clock_time_limit
the maximum elapsed clock time allowed (-ve means infinite)
bool remove_dependencies
the equality constraints will be preprocessed to remove any linear dependencies if true
bool treat_zero_bounds_as_general
any problem bound with the value zero will be treated as if it were a general value if true
bool treat_separable_as_general
if .just_feasible is true, the algorithm will stop as soon as a feasible point is found. Otherwise, the optimal solution to the problem will be found
bool just_feasible
if .treat_separable_as_general, is true, any separability in the problem structure will be ignored
bool getdua
if .getdua, is true, advanced initial values are obtained for the dual variables
bool puiseux
decide between Puiseux and Taylor series approximations to the arc
bool every_order
try every order of series up to series_order?
bool feasol
if .feasol is true, the final solution obtained will be perturbed so that variables close to their bounds are moved onto these bounds
bool balance_initial_complentarity
if .balance_initial_complentarity is true, the initial complemetarity is required to be balanced
bool crossover
if .crossover is true, cross over the solution to one defined by linearly-independent constraints if possible
bool space_critical
if .space_critical true, every effort will be made to use as little space as possible. This may result in longer computation time
bool deallocate_error_fatal
if .deallocate_error_fatal is true, any array/pointer deallocation error will terminate execution. Otherwise, computation will continue
bool generate_sif_file
if .generate_sif_file is .true. if a SIF file describing the current problem is to be generated
bool generate_qplib_file
if .generate_qplib_file is .true. if a QPLIB file describing the current problem is to be generated
char sif_file_name[31]
name of generated SIF file containing input problem
char qplib_file_name[31]
name of generated QPLIB file containing input problem
char prefix[31]
all output lines will be prefixed by .prefix(2:LEN(TRIM(.prefix))-1) where .prefix contains the required string enclosed in quotes, e.g. “string” or ‘string’
struct fdc_control_type fdc_control
control parameters for FDC
struct sbls_control_type sbls_control
control parameters for SBLS
struct fit_control_type fit_control
control parameters for FIT
struct roots_control_type roots_control
control parameters for ROOTS
struct cro_control_type cro_control
control parameters for CRO
cqp_time_type structure#
#include <galahad_cqp.h> struct cqp_time_type { // components rpc_ total; rpc_ preprocess; rpc_ find_dependent; rpc_ analyse; rpc_ factorize; rpc_ solve; rpc_ clock_total; rpc_ clock_preprocess; rpc_ clock_find_dependent; rpc_ clock_analyse; rpc_ clock_factorize; rpc_ clock_solve; };
detailed documentation#
time derived type as a C struct
components#
rpc_ total
the total CPU time spent in the package
rpc_ preprocess
the CPU time spent preprocessing the problem
rpc_ find_dependent
the CPU time spent detecting linear dependencies
rpc_ analyse
the CPU time spent analysing the required matrices prior to factorization
rpc_ factorize
the CPU time spent factorizing the required matrices
rpc_ solve
the CPU time spent computing the search direction
rpc_ clock_total
the total clock time spent in the package
rpc_ clock_preprocess
the clock time spent preprocessing the problem
rpc_ clock_find_dependent
the clock time spent detecting linear dependencies
rpc_ clock_analyse
the clock time spent analysing the required matrices prior to factorization
rpc_ clock_factorize
the clock time spent factorizing the required matrices
rpc_ clock_solve
the clock time spent computing the search direction
cqp_inform_type structure#
#include <galahad_cqp.h> struct cqp_inform_type { // components ipc_ status; ipc_ alloc_status; char bad_alloc[81]; ipc_ iter; ipc_ factorization_status; int64_t factorization_integer; int64_t factorization_real; ipc_ nfacts; ipc_ nbacts; ipc_ threads; rpc_ obj; rpc_ primal_infeasibility; rpc_ dual_infeasibility; rpc_ complementary_slackness; rpc_ init_primal_infeasibility; rpc_ init_dual_infeasibility; rpc_ init_complementary_slackness; rpc_ potential; rpc_ non_negligible_pivot; bool feasible; ipc_ checkpointsIter[16]; rpc_ checkpointsTime[16]; struct cqp_time_type time; struct fdc_inform_type fdc_inform; struct sbls_inform_type sbls_inform; struct fit_inform_type fit_inform; struct roots_inform_type roots_inform; struct cro_inform_type cro_inform; struct rpd_inform_type rpd_inform; };
detailed documentation#
inform derived type as a C struct
components#
ipc_ status
return status. See CQP_solve for details.
ipc_ alloc_status
the status of the last attempted allocation/deallocation.
char bad_alloc[81]
the name of the array for which an allocation/deallocation error occurred.
ipc_ iter
the total number of iterations required
ipc_ factorization_status
the return status from the factorization
int64_t factorization_integer
the total integer workspace required for the factorization
int64_t factorization_real
the total real workspace required for the factorization
ipc_ nfacts
the total number of factorizations performed
ipc_ nbacts
the total number of “wasted” function evaluations during the linesearch
ipc_ threads
the number of threads used
rpc_ obj
the value of the objective function at the best estimate of the solution determined by CQP_solve
rpc_ primal_infeasibility
the value of the primal infeasibility
rpc_ dual_infeasibility
the value of the dual infeasibility
rpc_ complementary_slackness
the value of the complementary slackness
rpc_ init_primal_infeasibility
these values at the initial point (needed by GALAHAD_CQP)
rpc_ init_dual_infeasibility
see init_primal_infeasibility
rpc_ init_complementary_slackness
see init_primal_infeasibility
rpc_ potential
the value of the logarithmic potential function sum -log(distance to constraint boundary)
rpc_ non_negligible_pivot
the smallest pivot which was not judged to be zero when detecting linear dependent constraints
bool feasible
is the returned “solution” feasible?
ipc_ checkpointsIter[16]
checkpoints(i) records the iteration at which the criticality measures first fall below \(10^{-i-1}\), i = 0, …, 15 (-1 means not achieved)
rpc_ checkpointsTime[16]
see checkpointsIter
struct cqp_time_type time
timings (see above)
struct fdc_inform_type fdc_inform
inform parameters for FDC
struct sbls_inform_type sbls_inform
inform parameters for SBLS
struct fit_inform_type fit_inform
return information from FIT
struct roots_inform_type roots_inform
return information from ROOTS
struct cro_inform_type cro_inform
inform parameters for CRO
struct rpd_inform_type rpd_inform
inform parameters for RPD
example calls#
This is an example of how to use the package to solve a given convex quadratic program; the code is available in $GALAHAD/src/cqp/C/cqpt.c . A variety of supported Hessian and constraint matrix storage formats are shown.
Notice that C-style indexing is used, and that this is flagged by setting
control.f_indexing
to false
. The floating-point type rpc_
is set in galahad_precision.h
to double
by default, but to float
if the preprocessor variable SINGLE
is defined. Similarly, the integer
type ipc_
from galahad_precision.h
is set to int
by default,
but to int64_t
if the preprocessor variable INTEGER_64
is defined.
/* cqpt.c */
/* Full test for the CQP C interface using C sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_cqp.h"
#ifdef REAL_128
#include <quadmath.h>
#endif
int main(void) {
// Derived types
void *data;
struct cqp_control_type control;
struct cqp_inform_type inform;
// Set problem data
ipc_ n = 3; // dimension
ipc_ m = 2; // number of general constraints
ipc_ H_ne = 3; // Hesssian elements
ipc_ H_row[] = {0, 1, 2 }; // row indices, NB lower triangle
ipc_ H_col[] = {0, 1, 2}; // column indices, NB lower triangle
ipc_ H_ptr[] = {0, 1, 2, 3}; // row pointers
rpc_ H_val[] = {1.0, 1.0, 1.0 }; // values
rpc_ g[] = {0.0, 2.0, 0.0}; // linear term in the objective
rpc_ f = 1.0; // constant term in the objective
ipc_ A_ne = 4; // Jacobian elements
ipc_ A_row[] = {0, 0, 1, 1}; // row indices
ipc_ A_col[] = {0, 1, 1, 2}; // column indices
ipc_ A_ptr[] = {0, 2, 4}; // row pointers
rpc_ A_val[] = {2.0, 1.0, 1.0, 1.0 }; // values
rpc_ c_l[] = {1.0, 2.0}; // constraint lower bound
rpc_ c_u[] = {2.0, 2.0}; // constraint upper bound
rpc_ x_l[] = {-1.0, - INFINITY, - INFINITY}; // variable lower bound
rpc_ x_u[] = {1.0, INFINITY, 2.0}; // variable upper bound
// Set output storage
rpc_ c[m]; // constraint values
ipc_ x_stat[n]; // variable status
ipc_ c_stat[m]; // constraint status
char st = ' ';
ipc_ status;
printf(" C sparse matrix indexing\n\n");
printf(" basic tests of qp storage formats\n\n");
for( ipc_ d=1; d <= 7; d++){
// Initialize CQP
cqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
// Start from 0
rpc_ x[] = {0.0,0.0,0.0};
rpc_ y[] = {0.0,0.0};
rpc_ z[] = {0.0,0.0,0.0};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
cqp_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
printf(" case %1" i_ipc_ " break\n",d);
case 2: // sparse by rows
st = 'R';
cqp_import( &control, &data, &status, n, m,
"sparse_by_rows", H_ne, NULL, H_col, H_ptr,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 3: // dense
st = 'D';
ipc_ H_dense_ne = 6; // number of elements of H
ipc_ A_dense_ne = 6; // number of elements of A
rpc_ H_dense[] = {1.0, 0.0, 1.0, 0.0, 0.0, 1.0};
rpc_ A_dense[] = {2.0, 1.0, 0.0, 0.0, 1.0, 1.0};
cqp_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
cqp_solve_qp( &data, &status, n, m, H_dense_ne, H_dense, g, f,
A_dense_ne, A_dense, c_l, c_u, x_l, x_u,
x, c, y, z, x_stat, c_stat );
break;
case 4: // diagonal
st = 'L';
cqp_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 5: // scaled identity
st = 'S';
cqp_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 6: // identity
st = 'I';
cqp_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 7: // zero
st = 'Z';
cqp_import( &control, &data, &status, n, m,
"zero", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
}
cqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: CQP_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
cqp_terminate( &data, &control, &inform );
}
// test shifted least-distance interface
for( ipc_ d=1; d <= 1; d++){
// Initialize CQP
cqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
// Start from 0
rpc_ x[] = {0.0,0.0,0.0};
rpc_ y[] = {0.0,0.0};
rpc_ z[] = {0.0,0.0,0.0};
// Set shifted least-distance data
rpc_ w[] = {1.0,1.0,1.0};
rpc_ x_0[] = {0.0,0.0,0.0};
switch(d){
case 1: // sparse co-ordinate storage
st = 'W';
cqp_import( &control, &data, &status, n, m,
"shifted_least_distance", H_ne, NULL, NULL, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
cqp_solve_sldqp( &data, &status, n, m, w, x_0, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
}
cqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: CQP_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
cqp_terminate( &data, &control, &inform );
}
}
This is the same example, but now fortran-style indexing is used; the code is available in $GALAHAD/src/cqp/C/cqptf.c .
/* cqptf.c */
/* Full test for the CQP C interface using Fortran sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_cqp.h"
#ifdef REAL_128
#include <quadmath.h>
#endif
int main(void) {
// Derived types
void *data;
struct cqp_control_type control;
struct cqp_inform_type inform;
// Set problem data
ipc_ n = 3; // dimension
ipc_ m = 2; // number of general constraints
ipc_ H_ne = 3; // Hesssian elements
ipc_ H_row[] = {1, 2, 3 }; // row indices, NB lower triangle
ipc_ H_col[] = {1, 2, 3}; // column indices, NB lower triangle
ipc_ H_ptr[] = {1, 2, 3, 4}; // row pointers
rpc_ H_val[] = {1.0, 1.0, 1.0 }; // values
rpc_ g[] = {0.0, 2.0, 0.0}; // linear term in the objective
rpc_ f = 1.0; // constant term in the objective
ipc_ A_ne = 4; // Jacobian elements
ipc_ A_row[] = {1, 1, 2, 2}; // row indices
ipc_ A_col[] = {1, 2, 2, 3}; // column indices
ipc_ A_ptr[] = {1, 3, 5}; // row pointers
rpc_ A_val[] = {2.0, 1.0, 1.0, 1.0 }; // values
rpc_ c_l[] = {1.0, 2.0}; // constraint lower bound
rpc_ c_u[] = {2.0, 2.0}; // constraint upper bound
rpc_ x_l[] = {-1.0, - INFINITY, - INFINITY}; // variable lower bound
rpc_ x_u[] = {1.0, INFINITY, 2.0}; // variable upper bound
// Set output storage
rpc_ c[m]; // constraint values
ipc_ x_stat[n]; // variable status
ipc_ c_stat[m]; // constraint status
char st = ' ';
ipc_ status;
printf(" Fortran sparse matrix indexing\n\n");
printf(" basic tests of qp storage formats\n\n");
for( ipc_ d=1; d <= 7; d++){
// Initialize CQP
cqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = true; // Fortran sparse matrix indexing
// Start from 0
rpc_ x[] = {0.0,0.0,0.0};
rpc_ y[] = {0.0,0.0};
rpc_ z[] = {0.0,0.0,0.0};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
cqp_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
printf(" case %1" i_ipc_ " break\n",d);
case 2: // sparse by rows
st = 'R';
cqp_import( &control, &data, &status, n, m,
"sparse_by_rows", H_ne, NULL, H_col, H_ptr,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 3: // dense
st = 'D';
ipc_ H_dense_ne = 6; // number of elements of H
ipc_ A_dense_ne = 6; // number of elements of A
rpc_ H_dense[] = {1.0, 0.0, 1.0, 0.0, 0.0, 1.0};
rpc_ A_dense[] = {2.0, 1.0, 0.0, 0.0, 1.0, 1.0};
cqp_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
cqp_solve_qp( &data, &status, n, m, H_dense_ne, H_dense, g, f,
A_dense_ne, A_dense, c_l, c_u, x_l, x_u,
x, c, y, z, x_stat, c_stat );
break;
case 4: // diagonal
st = 'L';
cqp_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 5: // scaled identity
st = 'S';
cqp_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 6: // identity
st = 'I';
cqp_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
case 7: // zero
st = 'Z';
cqp_import( &control, &data, &status, n, m,
"zero", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
cqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
}
cqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: CQP_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
cqp_terminate( &data, &control, &inform );
}
// test shifted least-distance interface
for( ipc_ d=1; d <= 1; d++){
// Initialize CQP
cqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = true; // Fortran sparse matrix indexing
// Start from 0
rpc_ x[] = {0.0,0.0,0.0};
rpc_ y[] = {0.0,0.0};
rpc_ z[] = {0.0,0.0,0.0};
// Set shifted least-distance data
rpc_ w[] = {1.0,1.0,1.0};
rpc_ x_0[] = {0.0,0.0,0.0};
switch(d){
case 1: // sparse co-ordinate storage
st = 'W';
cqp_import( &control, &data, &status, n, m,
"shifted_least_distance", H_ne, NULL, NULL, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
cqp_solve_sldqp( &data, &status, n, m, w, x_0, g, f,
A_ne, A_val, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat );
break;
}
cqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: CQP_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
cqp_terminate( &data, &control, &inform );
}
}