GALAHAD CRO package#
purpose#
The cro
package provides a crossover from a
primal-dual interior-point
solution to given convex quadratic program to a basic one in which
the matrix of defining active constraints/variables is of full rank.
This applies to the problem of minimizing the quadratic objective function
See Section 4 of $GALAHAD/doc/cro.pdf for additional details.
terminology#
Any required solution \(x\) necessarily satisfies the primal optimality conditions
method#
Denote the active constraints by \(A_A x = c_A\) and the active bounds by \(I_A x = x_A\). Then any optimal solution satisfies the linear system
ULS
. If \(K\) is
non singular, the solution is unique and the solution input by the user
provides a linearly independent active set. Otherwise \(K\) is singular,
and partitions \(A_A^T = ( A_{AB}^T \;\; A_{AN}^T)\) and
\(I_A^T = ( I_{AB}^T \;\; I_{AN}^T)\) are found so that
SLS
,
the alternative solution \((x + \alpha x, y + \alpha y, z + \alpha z)\),
featuring \((\Delta x, \Delta y_{AB}, \Delta z_{AB})\) from (3)
in which successively one of the components of \(\Delta y_{AN}\)
and \(\Delta z_{AN}\) in turn is non zero, is taken.
The scalar \(\alpha\) at each stage
is chosen to be the largest possible that guarantees (1);
this may happen when a non-basic multiplier/dual variable reaches zero,
in which case the corresponding constraint is disregarded, or when this
happens for a basic multiplier/dual variable, in which case this constraint is
exchanged with the non-basic one under consideration and disregarded.
The latter corresponds to changing the basic-non-basic partition
in (3), and subsequent solutions may be found by updating
the factorization of the coefficient matrix in (3)
following the basic-non-basic swap using the package SCU
.
matrix storage#
The unsymmetric \(m\) by \(n\) matrix \(A\) must be presented and stored in sparse row-wise storage format. For this, only the nonzero entries are stored, and 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+1) holds the total number of entries plus one. The column indices j, \(1 \leq j \leq n\), 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, \(1 \leq i \leq m\), of the integer array A_col, and real array A_val, respectively.
The symmetric \(n\) by \(n\) matrix \(H\) must also be presented and stored in sparse row-wise storage format. But, crucially, now symmetry is exploited by only storing values from the lower triangular part (i.e, those entries that lie on or below the leading diagonal). As before, only the nonzero entries of the matrices are stored. Only the nonzero entries from the lower triangle are stored, and these 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+1) holds the total number of entries plus one. The column indices j, \(1 \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.
introduction to function calls#
To solve a given problem, functions from the cro package must be called in the following order:
To solve a given problem, functions from the cro package must be called in the following order:
cro_initialize - provide default control parameters and set up initial data structures
cro_read_specfile (optional) - override control values by reading replacement values from a file
cro_crossover_solution - move from a primal-dual soution to a full rank one
cro_terminate - deallocate data structures
See the examples section for illustrations of use.
parametric real type T#
Below, the symbol T refers to a parametric real type that may be Float32 (single precision) or Float64 (double precision).
callable functions#
function cro_initialize(T, data, control, status)
Set default control values and initialize private data
Parameters:
data |
holds private internal data |
control |
is a structure containing control information (see cro_control_type) |
status |
is a scalar variable of type Int32 that gives the exit status from the package. Possible values are (currently):
|
function cro_read_specfile(T, control, 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/cro/CRO.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/cro.pdf for a list of how these keywords relate to the components of the control structure.
Parameters:
control |
is a structure containing control information (see cro_control_type) |
specfile |
is a one-dimensional array of type Vararg{Cchar} that must give the name of the specification file |
function cro_crossover_solution(T, data, control, inform, n, m, m_equal, h_ne, H_val, H_col, H_ptr, a_ne, A_val, A_col, A_ptr, g, c_l, c_u, x_l, x_u, x, c, y, z, x_stat, c_stat)
Crosover the solution from a primal-dual to a basic one.
Parameters:
control |
is a structure whose members provide control parameters for the remaining procedures (see cro_control_type). The parameter .status is as follows: |
data |
holds private internal data. |
inform |
is a structure containing output information (see cro_inform_type). The component .status gives the exit status from the package. Possible values are:
|
n |
is a scalar variable of type Int32 that holds the number of variables. |
m |
is a scalar variable of type Int32 that holds the number of general linear constraints. |
m_equal |
is a scalar variable of type Int32 that holds the number of general linear equality constraints. Such constraints must occur first in \(A\). |
h_ne |
is a scalar variable of type Int32 that holds the number of entries in the lower triangular part of the Hessian matrix \(H\). |
H_val |
is a one-dimensional array of type T that holds the values of the entries of the lower triangular part of the Hessian matrix \(H\). The entries are stored by consecutive rows, the order within each row is unimportant. |
H_col |
is a one-dimensional array of type Int32 that holds the column indices of the lower triangular part of \(H\), in the same order as those in H_val. |
H_ptr |
is a one-dimensional array of size n+1 and type Int32 that holds the starting position of each row of the lower triangular part of \(H\). The n+1-st component holds the total number of entries (plus one if fortran indexing is used). |
a_ne |
is a scalar variable of type Int32 that holds the number of entries in the constraint Jacobian matrix \(A\). |
A_val |
is a one-dimensional array of type T that holds the values of the entries of the constraint Jacobian matrix \(A\). The entries are stored by consecutive rows, the order within each row is unimportant. Equality constraints must be ordered first. |
A_col |
is a one-dimensional array of size A_ne and type Int32 that holds the column indices of \(A\) in the same order as those in A_val. |
A_ptr |
is a one-dimensional array of size m+1 and type Int32 that holds the starting position of each row of \(A\). The m+1-st component holds the total number of entries (plus one if fortran indexing is used). |
g |
is a one-dimensional array of size n and type T that holds the linear term \(g\) of the objective function. The j-th component of |
c_l |
is a one-dimensional array of size m and type T that holds the lower bounds \(c^l\) on the constraints \(A x\). The i-th component of |
c_u |
is a one-dimensional array of size m and type T that holds the upper bounds \(c^l\) on the constraints \(A x\). The i-th component of |
x_l |
is a one-dimensional array of size n and type T that holds the lower bounds \(x^l\) on the variables \(x\). The j-th component of |
x_u |
is a one-dimensional array of size n and type T that holds the upper bounds \(x^l\) on the variables \(x\). The j-th component of |
x |
is a one-dimensional array of size n and type T that holds the values \(x\) of the optimization variables. The j-th component of |
c |
is a one-dimensional array of size m and type T that holds the residual \(c(x) = A x\). The i-th component of |
y |
is a one-dimensional array of size n and type T that holds the values \(y\) of the Lagrange multipliers for the general linear constraints. The j-th component of |
z |
is a one-dimensional array of size n and type T that holds the values \(z\) of the dual variables. The j-th component of |
x_stat |
is a one-dimensional array of size n and type Int32 that must be set on entry to give the status of the problem variables. If x_stat(j) is negative, the variable \(x_j\) is active on its lower bound, if it is positive, it is active and lies on its upper bound, and if it is zero, it is inactiive and lies between its bounds. On exit, the \(j\) -th component of x_stat is -1 if the variable is basic and active on its lower bound, -2 it is non-basic but active on its lower bound, 1 if it is basic and active on its upper bound, 2 it is non-basic but active on its upper bound, and 0 if it is inactive. |
c_stat |
is a one-dimensional array of size m and type Int32 that must be set on entry to give the status of the general linear constraints. If c_stat(i) is negative, the constraint value \(a_i^Tx\) is active on its lower bound, if it is positive, it is active and lies on its upper bound, and if it is zero, it is inactiive and lies between its bounds. On exit, the \(i\) -th component of x_stat is -1 if the constraint is basic and active on its lower bound, -2 it is non-basic but active on its lower bound, 1 if it is basic and active on its upper bound, 2 it is non-basic but active on its upper bound, and 0 if it is inactive. |
function cro_terminate(T, data, control, inform)
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a structure containing control information (see cro_control_type) |
inform |
is a structure containing output information (see cro_inform_type) |
available structures#
cro_control_type structure#
struct cro_control_type{T} f_indexing::Bool error::Int32 out::Int32 print_level::Int32 max_schur_complement::Int32 infinity::T feasibility_tolerance::T check_io::Bool refine_solution::Bool space_critical::Bool deallocate_error_fatal::Bool symmetric_linear_solver::NTuple{31,Cchar} unsymmetric_linear_solver::NTuple{31,Cchar} prefix::NTuple{31,Cchar} sls_control::sls_control_type{T} sbls_control::sbls_control_type{T} uls_control::uls_control_type{T} ir_control::ir_control_type{T}
detailed documentation#
control derived type as a Julia structure
components#
Bool f_indexing
use C or Fortran sparse matrix indexing
Int32 error
error and warning diagnostics occur on stream error
Int32 out
general output occurs on stream out
Int32 print_level
the level of output required is specified by print_level
Int32 max_schur_complement
the maximum permitted size of the Schur complement before a refactorization is performed
T infinity
any bound larger than infinity in modulus will be regarded as infinite
T feasibility_tolerance
feasibility tolerance for KKT violation
Bool check_io
if .check_io is true, the input (x,y,z) will be fully tested for consistency
Bool refine_solution
if .refine solution is true, attempt to satisfy the KKT conditions as accurately as possible
Bool space_critical
if .space_critical is 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
char symmetric_linear_solver[31]
indefinite linear equation solver
char unsymmetric_linear_solver[31]
unsymmetric linear equation solver
NTuple{31,Cchar} prefix
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 sls_control_type sls_control
control parameters for SLS
struct sbls_control_type sbls_control
control parameters for SBLS
struct uls_control_type uls_control
control parameters for ULS
struct ir_control_type ir_control
control parameters for iterative refinement
cro_time_type structure#
struct cro_time_type{T} total::Float32 analyse::Float32 factorize::Float32 solve::Float32 clock_total::T clock_analyse::T clock_factorize::T clock_solve::T
detailed documentation#
time derived type as a Julia structure
components#
Float32 total
the total CPU time spent in the package
Float32 analyse
the CPU time spent reordering the matrix prior to factorization
Float32 factorize
the CPU time spent factorizing the required matrices
Float32 solve
the CPU time spent computing corrections
T clock_total
the total clock time spent in the package
T clock_analyse
the clock time spent analysing the required matrices prior to factorizat
T clock_factorize
the clock time spent factorizing the required matrices
T clock_solve
the clock time spent computing corrections
cro_inform_type structure#
struct cro_inform_type{T} status::Int32 alloc_status::Int32 bad_alloc::NTuple{81,Cchar} dependent::Int32 time::cro_time_type{T} sls_inform::sls_inform_type{T} sbls_inform::sbls_inform_type{T} uls_inform::uls_inform_type{T} scu_status::Int32 scu_inform::scu_inform_type ir_inform::ir_inform_type{T}
detailed documentation#
inform derived type as a Julia structure
components#
Int32 status
return status. See CRO_solve for details
Int32 alloc_status
the status of the last attempted allocation/deallocation
NTuple{81,Cchar} bad_alloc
the name of the array for which an allocation/deallocation error occurred
Int32 dependent
the number of dependent active constraints
struct cro_time_type time
timings (see above)
struct sls_inform_type sls_inform
information from SLS
struct sbls_inform_type sbls_inform
information from SBLS
struct uls_inform_type uls_inform
information from ULS
Int32 scu_status
information from SCU
struct scu_inform_type scu_inform
see scu_status
struct ir_inform_type ir_inform
information from IR
example calls#
This is an example of how to use the package to crossover from a primal-dual QP solution to a basic one; the code is available in $GALAHAD/src/cro/Julia/test_cro.jl . A variety of supported Hessian and constraint matrix storage formats are shown.
# test_cro.jl
# Simple code to test the Julia interface to CRO
using GALAHAD
using Test
using Printf
using Accessors
function test_cro(::Type{T}) where T
# Derived types
data = Ref{Ptr{Cvoid}}()
control = Ref{cro_control_type{T}}()
inform = Ref{cro_inform_type{T}}()
# Set problem dimensions
n = 11 # dimension
m = 3 # number of general constraints
m_equal = 1 # number of equality constraints
# describe the objective function
H_ne = 21
H_val = T[1.0, 0.5, 1.0, 0.5, 1.0, 0.5, 1.0, 0.5, 1.0, 0.5, 1.0, 0.5, 1.0, 0.5, 1.0,
0.5, 1.0, 0.5, 1.0, 0.5, 1.0]
H_col = Cint[1, 1, 2, 2, 3, 3, 4, 4, 5, 5, 6, 6, 7, 7, 8, 8, 9, 9, 10, 10, 11]
H_ptr = Cint[1, 2, 4, 6, 8, 10, 12, 14, 16, 18, 20, 22]
g = T[0.5, -0.5, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -0.5]
# describe constraints
A_ne = 30
A_val = T[1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0]
A_col = Cint[1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 3, 4, 5, 6, 7, 8, 9, 10, 11, 2, 3, 4, 5,
6, 7, 8, 9, 10, 11]
A_ptr = Cint[1, 12, 21, 31]
c_l = T[10.0, 9.0, -Inf]
c_u = T[10.0, Inf, 10.0]
x_l = T[0.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0]
x_u = T[Inf, Inf, Inf, Inf, Inf, Inf, Inf, Inf, Inf, Inf, Inf]
# provide optimal variables, Lagrange multipliers and dual variables
x = T[0.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0]
c = T[10.0, 9.0, 10.0]
y = T[-1.0, 1.5, -2.0]
z = T[2.0, 4.0, 2.5, 2.5, 2.5, 2.5, 2.5, 2.5, 2.5, 2.5, 2.5]
# provide interior-point constraint and variable status
c_stat = Cint[-1, -1, 1]
x_stat = Cint[-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1]
# Set output storage
status = Ref{Cint}()
@printf(" Fortran sparse matrix indexing\n\n")
# Initialize CRO
cro_initialize(T, data, control, status)
# Set user-defined control options
@reset control[].f_indexing = true # Fortran sparse matrix indexing
# crossover the solution
cro_crossover_solution(T, data, control, inform, n, m, m_equal, H_ne, H_val, H_col, H_ptr,
A_ne, A_val, A_col, A_ptr, g, c_l, c_u, x_l, x_u, x, c, y, z,
x_stat, c_stat)
@printf(" CRO_crossover exit status = %1i\n", inform[].status)
# Delete internal workspace
cro_terminate(T, data, control, inform)
return 0
end
@testset "CRO" begin
@test test_cro(Float32) == 0
@test test_cro(Float64) == 0
end