GALAHAD DQP package#
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
See Section 4 of $GALAHAD/doc/dqp.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#
Dual gradient-projection methods solve the quadratic programmimg problem by instead solving the dual quadratic program
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.
introduction to function calls#
To solve a given problem, functions from the dqp package must be called in the following order:
dqp_initialize - provide default control parameters and set up initial data structures
dqp_read_specfile (optional) - override control values by reading replacement values from a file
dqp_import - set up problem data structures and fixed values
dqp_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved
solve the problem by calling one of
dqp_solve_qp - solve the quadratic program
dqp_solve_sldqp - solve the shifted least-distance problem
dqp_information (optional) - recover information about the solution and solution process
dqp_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 dqp_control_type; struct dqp_inform_type; struct dqp_time_type; // function calls void dqp_initialize(void **data, struct dqp_control_type* control, ipc_ *status); void dqp_read_specfile(struct dqp_control_type* control, const char specfile[]); void dqp_import( struct dqp_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 dqp_reset_control( struct dqp_control_type* control, void **data, ipc_ *status ); void dqp_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 dqp_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 dqp_information(void **data, struct dqp_inform_type* inform, ipc_ *status); void dqp_terminate( void **data, struct dqp_control_type* control, struct dqp_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 SINGLE
.
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 dqp_initialize(void **data, struct dqp_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 dqp_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void dqp_read_specfile(struct dqp_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/dqp/DQP.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/dqp.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 dqp_control_type) |
specfile |
is a character string containing the name of the specification file |
void dqp_import( struct dqp_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 dqp_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’, or ‘identity’; 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 dqp_reset_control( struct dqp_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 dqp_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 dqp_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 dqp_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, 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 dqp_information(void **data, struct dqp_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see dqp_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void dqp_terminate( void **data, struct dqp_control_type* control, struct dqp_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see dqp_control_type) |
inform |
is a struct containing output information (see dqp_inform_type) |
available structures#
dqp_control_type structure#
#include <galahad_dqp.h> struct dqp_control_type { // components bool f_indexing; ipc_ error; ipc_ out; ipc_ print_level; ipc_ start_print; ipc_ stop_print; ipc_ print_gap; ipc_ dual_starting_point; ipc_ maxit; ipc_ max_sc; ipc_ cauchy_only; ipc_ arc_search_maxit; ipc_ cg_maxit; ipc_ explore_optimal_subspace; ipc_ restore_problem; ipc_ sif_file_device; ipc_ qplib_file_device; rpc_ rho; 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_ stop_cg_relative; rpc_ stop_cg_absolute; rpc_ cg_zero_curvature; rpc_ max_growth; rpc_ identical_bounds_tol; rpc_ cpu_time_limit; rpc_ clock_time_limit; rpc_ initial_perturbation; rpc_ perturbation_reduction; rpc_ final_perturbation; bool factor_optimal_matrix; bool remove_dependencies; bool treat_zero_bounds_as_general; bool exact_arc_search; bool subspace_direct; bool subspace_alternate; bool subspace_arc_search; bool space_critical; bool deallocate_error_fatal; bool generate_sif_file; bool generate_qplib_file; char symmetric_linear_solver[31]; char definite_linear_solver[31]; char unsymmetric_linear_solver[31]; char sif_file_name[31]; char qplib_file_name[31]; char prefix[31]; struct fdc_control_type fdc_control; struct sls_control_type sls_control; struct sbls_control_type sbls_control; struct gltr_control_type gltr_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
ipc_ start_print
any printing will start on this iteration
ipc_ stop_print
any printing will stop on this iteration
ipc_ print_gap
printing will only occur every print_gap iterations
ipc_ dual_starting_point
which starting point should be used for the dual problem
-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)
ipc_ maxit
at most maxit inner iterations are allowed
ipc_ max_sc
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)
ipc_ cauchy_only
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
ipc_ arc_search_maxit
how many iterations are allowed per arc search (-ve = as many as require
ipc_ cg_maxit
how many CG iterations to perform per DQP iteration (-ve reverts to n+1)
ipc_ explore_optimal_subspace
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
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_ 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_ rho
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)
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 infeasibilies
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_ stop_cg_relative
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 )
rpc_ stop_cg_absolute
see stop_cg_relative
rpc_ cg_zero_curvature
threshold below which curvature is regarded as zero if CG is used
rpc_ max_growth
maximum growth factor allowed without a refactorization
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_ cpu_time_limit
the maximum CPU time allowed (-ve means infinite)
rpc_ clock_time_limit
the maximum elapsed clock time allowed (-ve means infinite)
rpc_ initial_perturbation
the initial penalty weight (for DLP only)
rpc_ perturbation_reduction
the penalty weight reduction factor (for DLP only)
rpc_ final_perturbation
the final penalty weight (for DLP only)
bool factor_optimal_matrix
are the factors of the optimal augmented matrix required? (for DLP only)
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 exact_arc_search
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
bool subspace_direct
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.
bool subspace_alternate
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
bool subspace_arc_search
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
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 symmetric_linear_solver[31]
the name of the symmetric-indefinite linear equation solver used. Possible choices are currently: ‘sils’, ‘ma27’, ‘ma57’, ‘ma77’, ‘ma86’, ‘ma97’, ‘ssids’, ‘mumps’, ‘pardiso’, ‘mkl_pardiso’, ‘pastix’, ‘wsmp’, and ‘sytr’, although only ‘sytr’ and, for OMP 4.0-compliant compilers, ‘ssids’ are installed by default; others are easily installed (see README.external). More details of the capabilities of each solver are provided in the documentation for galahad_sls.
char definite_linear_solver[31]
the name of the definite linear equation solver used. Possible choices are currently: ‘sils’, ‘ma27’, ‘ma57’, ‘ma77’, ‘ma86’, ‘ma87’, ‘ma97’, ‘ssids’, ‘mumps’, ‘pardiso’, ‘mkl_pardiso’, ‘pastix’, ‘wsmp’, ‘potr’, ‘sytr’ and ‘pbtr’, although only ‘potr’, ‘sytr’, ‘pbtr’ and, for OMP 4.0-compliant compilers, ‘ssids’ are installed by default; others are easily installed (see README.external). More details of the capabilities of each solver are provided in the documentation for galahad_sls.
char unsymmetric_linear_solver[31]
the name of the unsymmetric linear equation solver used. Possible choices are currently: ‘gls’, ‘ma48’ and ‘getr’, although only ‘getr’ is installed by default; others are easily installed (see README.external). More details of the capabilities of each solver are provided in the documentation for galahad_uls.
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 sls_control_type sls_control
control parameters for SLS
struct sbls_control_type sbls_control
control parameters for SBLS
struct gltr_control_type gltr_control
control parameters for GLTR
dqp_time_type structure#
#include <galahad_dqp.h> struct dqp_time_type { // components rpc_ total; rpc_ preprocess; rpc_ find_dependent; rpc_ analyse; rpc_ factorize; rpc_ solve; rpc_ search; rpc_ clock_total; rpc_ clock_preprocess; rpc_ clock_find_dependent; rpc_ clock_analyse; rpc_ clock_factorize; rpc_ clock_solve; rpc_ clock_search; };
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_ search
the CPU time spent in the linesearch
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
rpc_ clock_search
the clock time spent in the linesearch
dqp_inform_type structure#
#include <galahad_dqp.h> struct dqp_inform_type { // components ipc_ status; ipc_ alloc_status; char bad_alloc[81]; ipc_ iter; ipc_ cg_iter; ipc_ factorization_status; int64_t factorization_integer; int64_t factorization_real; ipc_ nfacts; ipc_ threads; rpc_ obj; rpc_ primal_infeasibility; rpc_ dual_infeasibility; rpc_ complementary_slackness; rpc_ non_negligible_pivot; bool feasible; ipc_ checkpointsIter[16]; rpc_ checkpointsTime[16]; struct dqp_time_type time; struct fdc_inform_type fdc_inform; struct sls_inform_type sls_inform; struct sbls_inform_type sbls_inform; struct gltr_inform_type gltr_inform; ipc_ scu_status; struct scu_inform_type scu_inform; struct rpd_inform_type rpd_inform; };
detailed documentation#
inform derived type as a C struct
components#
ipc_ status
return status. See DQP_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_ cg_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_ threads
the number of threads used
rpc_ obj
the value of the objective function at the best estimate of the solution determined by DQP_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_ non_negligible_pivot
the smallest pivot that was not judged to be zero when detecting linearly 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 dqp_time_type time
timings (see above)
struct fdc_inform_type fdc_inform
inform parameters for FDC
struct sls_inform_type sls_inform
inform parameters for SLS
struct sbls_inform_type sbls_inform
inform parameters for SBLS
struct gltr_inform_type gltr_inform
return information from GLTR
ipc_ scu_status
inform parameters for SCU
struct scu_inform_type scu_inform
see scu_status
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/dqp/C/dqpt.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.
/* dqpt.c */
/* Full test for the DQP C interface using C sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include <string.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_dqp.h"
int main(void) {
// Derived types
void *data;
struct dqp_control_type control;
struct dqp_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
// ipc_ H_ne = 4; // Hesssian elements
// ipc_ H_row[] = {0, 1, 2, 2 }; // row indices, NB lower triangle
// ipc_ H_col[] = {0, 1, 1, 2}; // column indices, NB lower triangle
// ipc_ H_ptr[] = {0, 1, 3, 4}; // row pointers
// rpc_ H_val[] = {1.0, 2.0, 1.0, 3.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 <= 6; d++){
// Initialize DQP
dqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
strcpy(control.symmetric_linear_solver, "sytr ") ;
strcpy(control.definite_linear_solver, "sytr ") ;
strcpy(control.fdc_control.symmetric_linear_solver, "sytr ") ;
control.fdc_control.use_sls = true;
strcpy(control.sbls_control.symmetric_linear_solver, "sytr ") ;
strcpy(control.sbls_control.definite_linear_solver, "sytr ") ;
// 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';
dqp_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
dqp_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';
dqp_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 );
dqp_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};
dqp_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
dqp_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';
dqp_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
dqp_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';
dqp_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
dqp_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';
dqp_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
dqp_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;
}
dqp_information( &data, &inform, &status );
if(inform.status == 0){
printf("%c:%6" i_ipc_ " iterations. Optimal objective value = %5.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
}else{
printf("%c: DQP_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
dqp_terminate( &data, &control, &inform );
}
printf("\n separate test of sldqp\n\n");
// test shifted least-distance interface
// for( ipc_ d=1; d <= 0; d++){
for( ipc_ d=1; d <= 1; d++){
// Initialize DQP
dqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
strcpy(control.symmetric_linear_solver, "sytr ") ;
strcpy(control.definite_linear_solver, "sytr ") ;
strcpy(control.fdc_control.symmetric_linear_solver, "sytr ") ;
control.fdc_control.use_sls = true;
strcpy(control.sbls_control.symmetric_linear_solver, "sytr ") ;
strcpy(control.sbls_control.definite_linear_solver, "sytr ") ;
// 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';
dqp_import( &control, &data, &status, n, m,
"shifted_least_distance", H_ne, NULL, NULL, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
dqp_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;
}
dqp_information( &data, &inform, &status );
if(inform.status == 0){
printf("%c:%6" i_ipc_ " iterations. Optimal objective value = %5.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
}else{
printf("%c: DQP_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
dqp_terminate( &data, &control, &inform );
}
}
This is the same example, but now fortran-style indexing is used; the code is available in $GALAHAD/src/dqp/C/dqptf.c .
/* dqptf.c */
/* Full test for the DQP C interface using Fortran sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_dqp.h"
int main(void) {
// Derived types
void *data;
struct dqp_control_type control;
struct dqp_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 <= 6; d++){
// Initialize DQP
dqp_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';
dqp_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
dqp_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';
dqp_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 );
dqp_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};
dqp_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
dqp_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';
dqp_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
dqp_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';
dqp_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
dqp_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';
dqp_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
dqp_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;
}
dqp_information( &data, &inform, &status );
if(inform.status == 0){
printf("%c:%6" i_ipc_ " iterations. Optimal objective value = %5.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
}else{
printf("%c: DQP_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
dqp_terminate( &data, &control, &inform );
}
// test shifted least-distance interface
for( ipc_ d=1; d <= 1; d++){
// Initialize DQP
dqp_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';
dqp_import( &control, &data, &status, n, m,
"shifted_least_distance", H_ne, NULL, NULL, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
dqp_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;
}
dqp_information( &data, &inform, &status );
if(inform.status == 0){
printf("%c:%6" i_ipc_ " iterations. Optimal objective value = %5.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
}else{
printf("%c: DQP_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
dqp_terminate( &data, &control, &inform );
}
}