GALAHAD QPB package#
purpose#
The qpb
package uses a
primal-dual interior-point method to solve a given
non-convex quadratic program.
The aim is to minimize the quadratic objective function
If the matrix \(H\) is positive semi-definite, a global solution is found. However, if \(H\) is indefinite, the procedure may find a (weak second-order) critical point that is not the global solution to the given problem.
See Section 4 of $GALAHAD/doc/qpb.pdf for additional details.
terminolgy#
Any required solution \(x\) necessarily satisfies the primal optimality conditions
method#
Primal-dual interior point methods iterate towards a point that satisfies these conditions by ultimately aiming to satisfy (1), (3) and (5), while ensuring that (2) and (4) are satisfied as strict inequalities at each stage. Appropriate norms of the amounts by which (1), (3) and (5) fail to be satisfied are known as the primal and dual infeasibility, and the violation of complementary slackness, respectively. The fact that (2) and (4) are satisfied as strict inequalities gives such methods their other title, namely interior-point methods.
The problem is solved in two phases. The goal of the first
“initial feasible point” phase is
to find a strictly interior point which is primal feasible, that is that
(1) is satisfied. The package LSQP
is used for this
purpose, and offers the options of either accepting the first
strictly feasible point found, or preferably of aiming for the
so-called “analytic center” of the feasible region.
Having found such a suitable initial feasible point, the second “optimality”
phase ensures that (1) remains satisfied while iterating to
satisfy dual feasibility (3) and complementary slackness (5).
The optimality phase proceeds by approximately minimizing a
sequence of barrier functions
Each of the barrier subproblems is solved using a trust-region method. Such a method generates a trial correction step \(\Delta (x, c)\) to the current iterate \((x, c)\) by replacing the nonlinear barrier function locally by a suitable quadratic model, and approximately minimizing this model in the intersection of (1) and a trust region \(\|\Delta (x, c)\| \leq \Delta\) for some appropriate strictly positive trust-region radius \(\Delta\) and norm \(\| \cdot \|\). The step is accepted/rejected and the radius adjusted on the basis of how accurately the model reproduces the value of barrier function at the trial step. If the step proves to be unacceptable, a linesearch is performed along the step to obtain an acceptable new iterate. In practice, the natural primal “Newton” model of the barrier function is frequently less successful than an alternative primal-dual model, and consequently the primal-dual model is usually to be preferred.
Once a barrier subproblem has been solved, extrapolation based on values and derivatives encountered on the central path is optionally used to determine a good starting point for the next subproblem. Traditional Taylor-series extrapolation has been superceded by more accurate Puiseux-series methods as these are particularly suited to deal with degeneracy.
The trust-region subproblem is approximately solved using the combined
conjugate-gradient/Lanczos method implemented in the package GLTR
.
Such a method requires a suitable preconditioner, and in our case, the
only flexibility we have is in approximating the model of the
Hessian. Although using a fixed form of preconditioning is sometimes
effective, we have provided the option of an automatic choice, that aims
to balance the cost of applying the preconditioner against the needs for
an accurate solution of the trust-region subproblem. The preconditioner
is applied using the matrix factorization package SBLS
, but options
at this stage are to factorize the preconditioner as a whole (the
so-called “augmented system” approach), or to perform a block
elimination first (the “Schur-complement” approach). The latter is
usually to be prefered when a (non-singular) diagonal preconditioner is
used, but may be inefficient if any of the columns of \(A\) is too dense.
In order to make the solution as efficient as possible, the
variables and constraints are reordered internally
by the package QPP
prior to solution.
In particular, fixed variables, and
free (unbounded on both sides) constraints are temporarily removed.
reference#
The method is described in detail in
A. R. Conn, N. I. M. Gould, D. Orban and Ph. L. Toint, ``A primal-dual trust-region algorithm for minimizing a non-convex function subject to general inequality and linear equality constraints’’. Mathematical Programming 87 (1999) 215-249.
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 columnsare 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 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 qpb package must be called in the following order:
qpb_initialize - provide default control parameters and set up initial data structures
qpb_read_specfile (optional) - override control values by reading replacement values from a file
qpb_import - set up problem data structures and fixed values
qpb_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved
qpb_solve_qp - solve the quadratic program
qpb_information (optional) - recover information about the solution and solution process
qpb_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 qpb_control_type; struct qpb_inform_type; struct qpb_time_type; // function calls void qpb_initialize(void **data, struct qpb_control_type* control, ipc_ *status); void qpb_read_specfile(struct qpb_control_type* control, const char specfile[]); void qpb_import( struct qpb_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 qpb_reset_control( struct qpb_control_type* control, void **data, ipc_ *status ); void qpb_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 qpb_information(void **data, struct qpb_inform_type* inform, ipc_ *status); void qpb_terminate( void **data, struct qpb_control_type* control, struct qpb_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 qpb_initialize(void **data, struct qpb_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 qpb_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void qpb_read_specfile(struct qpb_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/qpb/QPB.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/qpb.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 qpb_control_type) |
specfile |
is a character string containing the name of the specification file |
void qpb_import( struct qpb_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 qpb_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 qpb_reset_control( struct qpb_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 qpb_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 qpb_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 qpb_information(void **data, struct qpb_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see qpb_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void qpb_terminate( void **data, struct qpb_control_type* control, struct qpb_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see qpb_control_type) |
inform |
is a struct containing output information (see qpb_inform_type) |
available structures#
qpb_control_type structure#
#include <galahad_qpb.h> struct qpb_control_type { // components bool f_indexing; ipc_ error; ipc_ out; ipc_ print_level; ipc_ start_print; ipc_ stop_print; ipc_ maxit; ipc_ itref_max; ipc_ cg_maxit; ipc_ indicator_type; ipc_ restore_problem; ipc_ extrapolate; ipc_ path_history; ipc_ factor; ipc_ max_col; ipc_ indmin; ipc_ valmin; ipc_ infeas_max; ipc_ precon; ipc_ nsemib; ipc_ path_derivatives; ipc_ fit_order; ipc_ sif_file_device; rpc_ infinity; rpc_ stop_p; rpc_ stop_d; rpc_ stop_c; rpc_ theta_d; rpc_ theta_c; rpc_ beta; rpc_ prfeas; rpc_ dufeas; rpc_ muzero; rpc_ reduce_infeas; rpc_ obj_unbounded; rpc_ pivot_tol; rpc_ pivot_tol_for_dependencies; rpc_ zero_pivot; rpc_ identical_bounds_tol; rpc_ inner_stop_relative; rpc_ inner_stop_absolute; rpc_ initial_radius; rpc_ mu_min; rpc_ inner_fraction_opt; 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 center; bool primal; bool puiseux; bool feasol; bool array_syntax_worse_than_do_loop; bool space_critical; bool deallocate_error_fatal; bool generate_sif_file; char sif_file_name[31]; char prefix[31]; struct lsqp_control_type lsqp_control; struct fdc_control_type fdc_control; struct sbls_control_type sbls_control; struct gltr_control_type gltr_control; struct fit_control_type fit_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_ maxit
at most maxit inner iterations are allowed
ipc_ itref_max
the maximum number of iterative refinements allowed
ipc_ cg_maxit
the maximum number of CG iterations allowed. If cg_maxit < 0, this number will be reset to the dimension of the system + 1
ipc_ indicator_type
specifies the type of indicator function used. Pssible values are
1 primal indicator: constraint active <=> distance to nearest bound <= .indicator_p_tol
2 primal-dual indicator: constraint active <=> distance to nearest bound <= .indicator_tol_pd * size of corresponding multiplier
3 primal-dual indicator: constraint active <=> distance to nearest bound <= .indicator_tol_tapia * distance to same bound at previous iteration
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_ extrapolate
should extrapolation be used to track the central path? Possible values
0 never
1 after the final major iteration
2 at each major iteration
ipc_ path_history
the maximum number of previous path points to use when fitting the data
ipc_ factor
the factorization to be used. Possible values are
0 automatic
1 Schur-complement factorization
2 augmented-system factorization
ipc_ max_col
the maximum number of nonzeros in a column of A which is permitted with the Schur-complement factorization
ipc_ indmin
an initial guess as to the integer workspace required by SBLS
ipc_ valmin
an initial guess as to the real workspace required by SBLS
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_ precon
the preconditioner to be used for the CG is defined by precon. Possible values are
0 automatic
1 no preconditioner, i.e, the identity within full factorization
2 full factorization
3 band within full factorization
4 diagonal using the barrier terms within full factorization
ipc_ nsemib
the semi-bandwidth of a band preconditioner, if appropriate
ipc_ path_derivatives
the maximum order of path derivative to use
ipc_ fit_order
the order of (Puiseux) series to fit to the path data: <=0 to fit all data
ipc_ sif_file_device
specifies the unit number to write generated SIF file describing the current problem
rpc_ infinity
any bound larger than infinity in modulus will be regarded as infinite
rpc_ stop_p
the required accuracy for the primal infeasibility
rpc_ stop_d
the required accuracy for the dual infeasibility
rpc_ stop_c
the required accuracy for the complementarity
rpc_ theta_d
tolerances used to terminate the inner iteration (for given mu): dual feasibility <= MAX( theta_d * mu ** beta, 0.99 * stop_d ) complementarity <= MAX( theta_c * mu ** beta, 0.99 * stop_d )
rpc_ theta_c
see theta_d
rpc_ beta
see theta_d
rpc_ prfeas
initial primal variables will not be closer than prfeas from their bound
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_ 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_ pivot_tol
the threshold pivot used by the matrix factorization. See the documentation for SBLS for details
rpc_ pivot_tol_for_dependencies
the threshold pivot used by the matrix factorization when attempting to detect linearly dependent constraints. See the documentation for FDC for details
rpc_ zero_pivot
any pivots smaller than zero_pivot in absolute value will be regarded to zero when attempting to detect linearly dependent constraints
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_ inner_stop_relative
the search direction is considered as an acceptable approximation to the minimizer of the model if the gradient of the model in the preconditioning(inverse) norm is less than max( inner_stop_relative * initial preconditioning(inverse) gradient norm, inner_stop_absolute )
rpc_ inner_stop_absolute
see inner_stop_relative
rpc_ initial_radius
the initial trust-region radius
rpc_ mu_min
start terminal extrapolation when mu reaches mu_min
rpc_ inner_fraction_opt
a search direction which gives at least inner_fraction_opt times the optimal model decrease will be found
rpc_ indicator_tol_p
if .indicator_type = 1, a constraint/bound will be deemed to be active <=> distance to nearest bound <= .indicator_p_tol
rpc_ indicator_tol_pd
if .indicator_type = 2, a constraint/bound will be deemed to be active <=> distance to nearest bound <= .indicator_tol_pd * size of corresponding multiplier
rpc_ indicator_tol_tapia
if .indicator_type = 3, a constraint/bound will be deemed to be active <=> distance to nearest bound <= .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 center
if .center is true, the algorithm will use the analytic center of the feasible set as its initial feasible point. Otherwise, a feasible point as close as possible to the initial point will be used. We recommend using the analytic center
bool primal
if .primal, is true, a primal barrier method will be used in place of t primal-dual method
bool puiseux
If extrapolation is to be used, decide between Puiseux and Taylor series.
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 array_syntax_worse_than_do_loop
if .array_syntax_worse_than_do_loop is true, f77-style do loops will be used rather than f90-style array syntax for vector operations
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
char sif_file_name[31]
name of generated SIF 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 lsqp_control_type lsqp_control
control parameters for LSQP
struct fdc_control_type fdc_control
control parameters for FDC
struct sbls_control_type sbls_control
control parameters for SBLS
struct gltr_control_type gltr_control
control parameters for GLTR
struct fit_control_type fit_control
control parameters for FIT
qpb_time_type structure#
#include <galahad_qpb.h> struct qpb_time_type { // components rpc_ total; rpc_ preprocess; rpc_ find_dependent; rpc_ analyse; rpc_ factorize; rpc_ solve; rpc_ phase1_total; rpc_ phase1_analyse; rpc_ phase1_factorize; rpc_ phase1_solve; rpc_ clock_total; rpc_ clock_preprocess; rpc_ clock_find_dependent; rpc_ clock_analyse; rpc_ clock_factorize; rpc_ clock_solve; rpc_ clock_phase1_total; rpc_ clock_phase1_analyse; rpc_ clock_phase1_factorize; rpc_ clock_phase1_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 factorizatio
rpc_ factorize
the CPU time spent factorizing the required matrices
rpc_ solve
the CPU time spent computing the search direction
rpc_ phase1_total
the total CPU time spent in the initial-point phase of the package
rpc_ phase1_analyse
the CPU time spent analysing the required matrices prior to factorizatio in the inital-point phase
rpc_ phase1_factorize
the CPU time spent factorizing the required matrices in the inital-point phase
rpc_ phase1_solve
the CPU time spent computing the search direction in the inital-point ph
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 factorizat
rpc_ clock_factorize
the clock time spent factorizing the required matrices
rpc_ clock_solve
the clock time spent computing the search direction
rpc_ clock_phase1_total
the total clock time spent in the initial-point phase of the package
rpc_ clock_phase1_analyse
the clock time spent analysing the required matrices prior to factorizat in the inital-point phase
rpc_ clock_phase1_factorize
the clock time spent factorizing the required matrices in the inital-poi phase
rpc_ clock_phase1_solve
the clock time spent computing the search direction in the inital-point
qpb_inform_type structure#
#include <galahad_qpb.h> struct qpb_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_ nbacts; ipc_ nmods; rpc_ obj; rpc_ non_negligible_pivot; bool feasible; struct qpb_time_type time; struct lsqp_inform_type lsqp_inform; struct fdc_inform_type fdc_inform; struct sbls_inform_type sbls_inform; struct gltr_inform_type gltr_inform; struct fit_inform_type fit_inform; };
detailed documentation#
inform derived type as a C struct
components#
ipc_ status
return status. See QPB_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 conjugate gradient 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_ nmods
the total number of factorizations which were modified to ensure that th matrix was an appropriate preconditioner
rpc_ obj
the value of the objective function at the best estimate of the solution determined by QPB_solve
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?
struct qpb_time_type time
timings (see above)
struct lsqp_inform_type lsqp_inform
inform parameters for LSQP
struct fdc_inform_type fdc_inform
inform parameters for FDC
struct sbls_inform_type sbls_inform
inform parameters for SBLS
struct gltr_inform_type gltr_inform
return information from GLTR
struct fit_inform_type fit_inform
return information from FIT
example calls#
This is an example of how to use the package to solve a quadratic program; the code is available in $GALAHAD/src/qpb/C/qpbt.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.
/* qpbt.c */
/* Full test for the QPB C interface using C sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_qpb.h"
int main(void) {
// Derived types
void *data;
struct qpb_control_type control;
struct qpb_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 QPB
qpb_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';
qpb_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
qpb_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';
qpb_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 );
qpb_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};
qpb_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"zero", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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;
}
qpb_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: QPB_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
qpb_terminate( &data, &control, &inform );
}
}
This is the same example, but now fortran-style indexing is used; the code is available in $GALAHAD/src/qpb/C/qpbtf.c .
/* qpbtf.c */
/* Full test for the QPB C interface using Fortran sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_qpb.h"
int main(void) {
// Derived types
void *data;
struct qpb_control_type control;
struct qpb_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 QPB
qpb_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';
qpb_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
qpb_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';
qpb_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 );
qpb_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};
qpb_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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';
qpb_import( &control, &data, &status, n, m,
"zero", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
qpb_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;
}
qpb_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: QPB_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
qpb_terminate( &data, &control, &inform );
}
}