GALAHAD EQP package#
purpose#
The eqp
package uses an iterative method to solve a
given equality-constrained quadratic program.
The aim is to minimize the quadratic objective function
See Section 4 of $GALAHAD/doc/eqp.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#
A solution to the problem is found in two phases. In the first, a point \(x_F\) satisfying (1) is found. In the second, the required solution \(x = x_F + s\) is determined by finding \(s\) to minimize \(q(s) = \frac{1}{2} s^T H s + g_F^T s + f_F^{}\) subject to the homogeneous constraints \(A s = 0\), where \(g_F^{} = H x_F^{} + g\) and \(f_F^{} = \frac{1}{2} x_F^T H x_F^{} + g^T x_F^{} + f\). The required constrained minimizer of \(q(s)\) is obtained by implictly applying the preconditioned conjugate-gradient method in the null space of \(A\). Any preconditioner of the form
SBLS
provides a number of possibilities. In order to ensure that the
minimizer obtained is finite, an additional, precautionary trust-region
constraint \(\|s\| \leq \Delta\) for some suitable positive radius
\(\Delta\) is imposed, and the package GLTR
is used to solve
this additionally-constrained problem.
references#
The preconditioning aspcets are described in detail in
H. S. Dollar, N. I. M. Gould and A. J. Wathen. ``On implicit-factorization constraint preconditioners’’. In Large Scale Nonlinear Optimization (G. Di Pillo and M. Roma, eds.) Springer Series on Nonconvex Optimization and Its Applications, Vol. 83, Springer Verlag (2006) 61–82
and
H. S. Dollar, N. I. M. Gould, W. H. A. Schilders and A. J. Wathen ``On iterative methods and implicit-factorization preconditioners for regularized saddle-point systems’’. SIAM Journal on Matrix Analysis and Applications 28(1) (2006) 170–189,
while the constrained conjugate-gradient method is discussed in
N. I. M. Gould, S. Lucidi, M. Roma and Ph. L. Toint, ``Solving the trust-region subproblem using the Lanczos method’’. SIAM Journal on Optimization 9(2) (1999), 504–525.
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 eqp package must be called in the following order:
eqp_initialize - provide default control parameters and set up initial data structures
eqp_read_specfile (optional) - override control values by reading replacement values from a file
eqp_import - set up problem data structures and fixed values
eqp_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved
solve the problem by calling one of
eqp_solve_qp - solve the quadratic program
eqp_solve_sldqp - solve the shifted least-distance problem
eqp_resolve_qp (optional) - resolve the problem with the same Hessian and Jacobian, but different \(g\), \(f\) and/or \(c\)
eqp_information (optional) - recover information about the solution and solution process
eqp_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 eqp_control_type; struct eqp_inform_type; struct eqp_time_type; // function calls void eqp_initialize(void **data, struct eqp_control_type* control, ipc_ *status); void eqp_read_specfile(struct eqp_control_type* control, const char specfile[]); void eqp_import( struct eqp_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 eqp_reset_control( struct eqp_control_type* control, void **data, ipc_ *status ); void eqp_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[], rpc_ c[], rpc_ x[], rpc_ y[] ); void eqp_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[], rpc_ c[], rpc_ x[], rpc_ y[] ); void eqp_resolve_qp( void **data, ipc_ *status, ipc_ n, ipc_ m, const rpc_ g[], const rpc_ f, rpc_ c[], rpc_ x[], rpc_ y[] ); void eqp_information(void **data, struct eqp_inform_type* inform, ipc_ *status); void eqp_terminate( void **data, struct eqp_control_type* control, struct eqp_inform_type* inform );
typedefs#
typedef float spc_
spc_
is real single precision
typedef double rpc_
rpc_
is the real working precision used, but may be changed to float
by
defining the preprocessor variable REAL_32
or (if supported) to
__real128
using the variable REAL_128
.
typedef int ipc_
ipc_
is the default integer word length used, but may be changed to
int64_t
by defining the preprocessor variable INTEGER_64
.
function calls#
void eqp_initialize(void **data, struct eqp_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 eqp_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void eqp_read_specfile(struct eqp_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/eqp/EQP.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/eqp.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 eqp_control_type) |
specfile |
is a character string containing the name of the specification file |
void eqp_import( struct eqp_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 eqp_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 eqp_reset_control( struct eqp_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 eqp_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 eqp_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[], rpc_ c[], rpc_ x[], rpc_ y[] )
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 |
is a one-dimensional array of size m and type rpc_, that holds the linear term \(c\) in the constraints. The i-th component of c, i = 0, … , m-1, contains \(c_i\). |
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\). |
y |
is a one-dimensional array of size n and type rpc_, that holds the values \(y\) of the Lagrange multipliers for the linear constraints. The j-th component of y, i = 0, … , m-1, contains \(y_i\). |
void eqp_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[], rpc_ c[], rpc_ x[], rpc_ y[] )
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 |
is a one-dimensional array of size m and type rpc_, that holds the linear term \(c\) in the constraints. The i-th component of c, i = 0, … , m-1, contains \(c_i\). |
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\). |
y |
is a one-dimensional array of size n and type rpc_, that holds the values \(y\) of the Lagrange multipliers for the linear constraints. The j-th component of y, i = 0, … , m-1, contains \(y_i\). |
void eqp_resolve_qp( void **data, ipc_ *status, ipc_ n, ipc_ m, const rpc_ g[], const rpc_ f, rpc_ c[], rpc_ x[], rpc_ y[] )
Resolve the quadratic program or shifted least-distance quadratic program when some or all of the data \(g\), \(f\) and \(c\) has changed
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. |
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. |
c |
is a one-dimensional array of size m and type rpc_, that holds the linear term \(c\) in the constraints. The i-th component of c, i = 0, … , m-1, contains \(c_i\). |
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\). |
y |
is a one-dimensional array of size n and type rpc_, that holds the values \(y\) of the Lagrange multipliers for the linear constraints. The j-th component of y, i = 0, … , m-1, contains \(y_i\). |
void eqp_information(void **data, struct eqp_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see eqp_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void eqp_terminate( void **data, struct eqp_control_type* control, struct eqp_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see eqp_control_type) |
inform |
is a struct containing output information (see eqp_inform_type) |
available structures#
eqp_control_type structure#
#include <galahad_eqp.h> struct eqp_control_type { // components bool f_indexing; ipc_ error; ipc_ out; ipc_ print_level; ipc_ factorization; ipc_ max_col; ipc_ indmin; ipc_ valmin; ipc_ len_ulsmin; ipc_ itref_max; ipc_ cg_maxit; ipc_ preconditioner; ipc_ semi_bandwidth; ipc_ new_a; ipc_ new_h; ipc_ sif_file_device; rpc_ pivot_tol; rpc_ pivot_tol_for_basis; rpc_ zero_pivot; rpc_ inner_fraction_opt; rpc_ radius; rpc_ min_diagonal; rpc_ max_infeasibility_relative; rpc_ max_infeasibility_absolute; rpc_ inner_stop_relative; rpc_ inner_stop_absolute; rpc_ inner_stop_inter; bool find_basis_by_transpose; bool remove_dependencies; bool space_critical; bool deallocate_error_fatal; bool generate_sif_file; char sif_file_name[31]; char prefix[31]; struct fdc_control_type fdc_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_ factorization
the factorization to be used. Possible values are /li 0 automatic /li 1 Schur-complement factorization /li 2 augmented-system factorization (OBSOLETE)
ipc_ max_col
the maximum number of nonzeros in a column of A which is permitted with the Schur-complement factorization (OBSOLETE)
ipc_ indmin
an initial guess as to the integer workspace required by SBLS (OBSOLETE)
ipc_ valmin
an initial guess as to the real workspace required by SBLS (OBSOLETE)
ipc_ len_ulsmin
an initial guess as to the workspace required by ULS (OBSOLETE)
ipc_ itref_max
the maximum number of iterative refinements allowed (OBSOLETE)
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_ preconditioner
the preconditioner to be used for the CG. 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 (OBSOLETE)
5 optionally supplied diagonal, G = D
ipc_ semi_bandwidth
the semi-bandwidth of a band preconditioner, if appropriate (OBSOLETE)
ipc_ new_a
how much has A changed since last problem solved: 0 = not changed, 1 = values changed, 2 = structure changed
ipc_ new_h
how much has H changed since last problem solved: 0 = not changed, 1 = values changed, 2 = structure changed
ipc_ sif_file_device
specifies the unit number to write generated SIF file describing the current problem
rpc_ pivot_tol
the threshold pivot used by the matrix factorization. See the documentation for SBLS for details (OBSOLETE)
rpc_ pivot_tol_for_basis
the threshold pivot used by the matrix factorization when finding the ba See the documentation for ULS for details (OBSOLETE)
rpc_ zero_pivot
any pivots smaller than zero_pivot in absolute value will be regarded to zero when attempting to detect linearly dependent constraints (OBSOLETE)
rpc_ inner_fraction_opt
the computed solution which gives at least inner_fraction_opt times the optimal value will be found (OBSOLETE)
rpc_ radius
an upper bound on the permitted step (-ve will be reset to an appropriat large value by eqp_solve)
rpc_ min_diagonal
diagonal preconditioners will have diagonals no smaller than min_diagonal (OBSOLETE)
rpc_ max_infeasibility_relative
if the constraints are believed to be rank defficient and the residual at a “typical” feasible point is larger than max( max_infeasibility_relative * norm A, max_infeasibility_absolute ) the problem will be marked as infeasible
rpc_ max_infeasibility_absolute
see max_infeasibility_relative
rpc_ inner_stop_relative
the computed solution is considered as an acceptable approximation to th minimizer of the problem if the gradient of the objective 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_ inner_stop_inter
see inner_stop_relative
bool find_basis_by_transpose
if .find_basis_by_transpose is true, implicit factorization precondition will be based on a basis of A found by examining A’s transpose (OBSOLETE)
bool remove_dependencies
if .remove_dependencies is true, the equality constraints will be preprocessed to remove any linear dependencies
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 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
eqp_time_type structure#
#include <galahad_eqp.h> struct eqp_time_type { // components rpc_ total; rpc_ find_dependent; rpc_ factorize; rpc_ solve; rpc_ solve_inter; rpc_ clock_total; rpc_ clock_find_dependent; rpc_ clock_factorize; rpc_ clock_solve; };
detailed documentation#
time derived type as a C struct
components#
rpc_ total
the total CPU time spent in the package
rpc_ find_dependent
the CPU time spent detecting linear dependencies
rpc_ factorize
the CPU time spent factorizing the required matrices
rpc_ solve
the CPU time spent computing the search direction
rpc_ solve_inter
see solve
rpc_ clock_total
the total clock time spent in the package
rpc_ clock_find_dependent
the clock time spent detecting linear dependencies
rpc_ clock_factorize
the clock time spent factorizing the required matrices
rpc_ clock_solve
the clock time spent computing the search direction
eqp_inform_type structure#
#include <galahad_eqp.h> struct eqp_inform_type { // components ipc_ status; ipc_ alloc_status; char bad_alloc[81]; ipc_ cg_iter; ipc_ cg_iter_inter; int64_t factorization_integer; int64_t factorization_real; rpc_ obj; struct eqp_time_type time; struct fdc_inform_type fdc_inform; struct sbls_inform_type sbls_inform; struct gltr_inform_type gltr_inform; };
detailed documentation#
inform derived type as a C struct
components#
ipc_ status
return status. See EQP_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_ cg_iter
the total number of conjugate gradient iterations required
ipc_ cg_iter_inter
see cg_iter
int64_t factorization_integer
the total integer workspace required for the factorization
int64_t factorization_real
the total real workspace required for the factorization
rpc_ obj
the value of the objective function at the best estimate of the solution determined by QPB_solve
struct eqp_time_type time
timings (see above)
struct fdc_inform_type fdc_inform
inform parameters for FDC
struct sbls_inform_type sbls_inform
inform parameters for SBLS
struct gltr_inform_type gltr_inform
return information from GLTR
example calls#
This is an example of how to use the package to solve an equality-constrained quadratic program; the code is available in $GALAHAD/src/eqp/C/eqpt.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.
/* eqpt.c */
/* Full test for the EQP 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_eqp.h"
#ifdef REAL_128
#include <quadmath.h>
#endif
int main(void) {
// Derived types
void *data;
struct eqp_control_type control;
struct eqp_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[] = {3.0, 0.0}; // rhs of the constraints
// Set output storage
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 EQP
eqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
control.fdc_control.use_sls = true ;
strcpy(control.fdc_control.symmetric_linear_solver, "sytr ") ;
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};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
eqp_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 2: // sparse by rows
st = 'R';
eqp_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 );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
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};
eqp_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
eqp_solve_qp( &data, &status, n, m, H_dense_ne, H_dense, g, f,
A_dense_ne, A_dense, c, x, y );
break;
case 4: // diagonal
st = 'L';
eqp_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 5: // scaled identity
st = 'S';
eqp_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 6: // identity
st = 'I';
eqp_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 7: // zero
st = 'Z';
eqp_import( &control, &data, &status, n, m,
"zero", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
}
eqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_if.h
#include "galahad_pquad_if.h"
#else
printf("%c:%6" i_ipc_ " cg iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.cg_iter, inform.obj, inform.status);
#endif
}else{
printf("%c: EQP_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
eqp_terminate( &data, &control, &inform );
}
// test shifted least-distance interface
for( ipc_ d=1; d <= 1; d++){
// Initialize EQP
eqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
control.fdc_control.use_sls = true ;
strcpy(control.fdc_control.symmetric_linear_solver, "sytr ") ;
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};
// 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';
eqp_import( &control, &data, &status, n, m,
"shifted_least_distance", H_ne, NULL, NULL, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
eqp_solve_sldqp( &data, &status, n, m, w, x_0, g, f,
A_ne, A_val, c, x, y );
break;
}
eqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_if.h
#include "galahad_pquad_if.h"
#else
printf("%c:%6" i_ipc_ " cg iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.cg_iter, inform.obj, inform.status);
#endif
}else{
printf("%c: EQP_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
eqp_terminate( &data, &control, &inform );
}
}
This is the same example, but now fortran-style indexing is used; the code is available in $GALAHAD/src/eqp/C/eqptf.c .
/* eqptf.c */
/* Full test for the EQP C interface using Fortran sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include <string.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_eqp.h"
#ifdef REAL_128
#include <quadmath.h>
#endif
int main(void) {
// Derived types
void *data;
struct eqp_control_type control;
struct eqp_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[] = {3.0, 0.0}; // rhs of the constraints
// Set output storage
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 EQP
eqp_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = true; // Fortran sparse matrix indexing
control.fdc_control.use_sls = true ;
strcpy(control.fdc_control.symmetric_linear_solver, "sytr ") ;
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};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
eqp_import( &control, &data, &status, n, m,
"coordinate", H_ne, H_row, H_col, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 2: // sparse by rows
st = 'R';
eqp_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 );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
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};
eqp_import( &control, &data, &status, n, m,
"dense", H_ne, NULL, NULL, NULL,
"dense", A_ne, NULL, NULL, NULL );
eqp_solve_qp( &data, &status, n, m, H_dense_ne, H_dense, g, f,
A_dense_ne, A_dense, c, x, y );
break;
case 4: // diagonal
st = 'L';
eqp_import( &control, &data, &status, n, m,
"diagonal", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 5: // scaled identity
st = 'S';
eqp_import( &control, &data, &status, n, m,
"scaled_identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 6: // identity
st = 'I';
eqp_import( &control, &data, &status, n, m,
"identity", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
case 7: // zero
st = 'Z';
eqp_import( &control, &data, &status, n, m,
"zero", H_ne, NULL, NULL, NULL,
"sparse_by_rows", A_ne, NULL, A_col, A_ptr );
eqp_solve_qp( &data, &status, n, m, H_ne, H_val, g, f,
A_ne, A_val, c, x, y );
break;
}
eqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
#include "galahad_pquad_if.h"
#else
printf("%c:%6" i_ipc_ " cg iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.cg_iter, inform.obj, inform.status);
#endif
}else{
printf("%c: EQP_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
eqp_terminate( &data, &control, &inform );
}
// test shifted least-distance interface
for( ipc_ d=1; d <= 1; d++){
// Initialize EQP
eqp_initialize( &data, &control, &status );
control.fdc_control.use_sls = true ;
strcpy(control.fdc_control.symmetric_linear_solver, "sytr ") ;
strcpy(control.sbls_control.symmetric_linear_solver, "sytr ") ;
strcpy(control.sbls_control.definite_linear_solver, "sytr ") ;
// 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};
// 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';
eqp_import( &control, &data, &status, n, m,
"shifted_least_distance", H_ne, NULL, NULL, NULL,
"coordinate", A_ne, A_row, A_col, NULL );
eqp_solve_sldqp( &data, &status, n, m, w, x_0, g, f,
A_ne, A_val, c, x, y );
break;
}
eqp_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
#include "galahad_pquad_if.h"
#else
printf("%c:%6" i_ipc_ " cg iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.cg_iter, inform.obj, inform.status);
#endif
}else{
printf("%c: EQP_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
eqp_terminate( &data, &control, &inform );
}
}