RQS#
purpose#
The rqs
package uses matrix factorization to find the
global minimizer of a regularized quadratic objective function.
The aim is to minimize the regularized quadratic objective function
The matrix
Factorization of matrices of the form
glrt
may be preferred.
See Section 4 of $GALAHAD/doc/rqs.pdf for a brief description of the method employed and other details.
See Section 4 of $GALAHAD/doc/trs.pdf for additional details.
method#
The required solution
The method is iterative, and proceeds in two phases.
Firstly, lower and upper bounds,
The dominant cost is the requirement that we solve a sequence of
linear systems (2). In the absence of linear constraints, an
efficient sparse Cholesky factorization with precautions to
detect indefinite
reference#
The method is described in detail in
H. S. Dollar, N. I. M. Gould and D. P. Robinson. ``On solving trust-region and other regularised subproblems in optimization’’. Mathematical Programming Computation 2(1) (2010) 21–57.
matrix storage#
The unsymmetric
Dense storage format:
The matrix
Dense by columns storage format:
The matrix
Sparse co-ordinate storage format:
Only the nonzero entries of the matrices are stored.
For the
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
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
The symmetric
Dense storage format:
The matrix
Sparse co-ordinate storage format:
Only the nonzero entries of the matrices are stored.
For the
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
Diagonal storage format:
If
Multiples of the identity storage format:
If
The identity matrix format:
If
The zero matrix format:
The same is true if
functions#
- rqs.initialize()#
Set default option values and initialize private data
Returns:
- optionsdict
- dictionary containing default control options:
- errorint
error and warning diagnostics occur on stream error.
- outint
general output occurs on stream out.
- problemint
unit to write problem data into file problem_file.
- print_levelint
the level of output required is specified by print_level. Possible values are
<=0
gives no output,
1
gives a one-line summary for every iteration.
2
gives a summary of the inner iteration for each iteration.
>=3
gives increasingly verbose (debugging) output.
- dense_factorizationint
should the problem be solved by dense factorization? Possible values are
0
sparse factorization will be used
1
dense factorization will be used
other
the choice is made automatically depending on the
dimension and sparsity.
- new_hint
how much of
has changed since the previous call. Possible values are
0
unchanged
1
values but not indices have changed
2
values and indices have changed.
- new_mint
how much of
has changed since the previous call. Possible values are
0
unchanged
1
values but not indices have changed
2
values and indices have changed.
- new_aint
how much of
has changed since the previous call. Possible values are 0 unchanged 1 values but not indices have changed 2 values and indices have changed. - max_factorizationsint
the maximum number of factorizations (=iterations) allowed. -ve implies no limit.
- inverse_itmaxint
the number of inverse iterations performed in the “maybe hard” case.
- taylor_max_degreeint
maximum degree of Taylor approximant allowed.
- initial_multiplierfloat
initial estimate of the Lagrange multipler.
- lowerfloat
lower and upper bounds on the multiplier, if known.
- upperfloat
see lower.
- stop_normalfloat
stop when
stop_normal * max - stop_hardfloat
stop when bracket on optimal multiplier <= stop_hard * max( bracket ends ).
- start_invit_tolfloat
start inverse iteration when bracket on optimal multiplier <= stop_start_invit_tol * max( bracket ends ).
- start_invitmax_tolfloat
start full inverse iteration when bracket on multiplier <= stop_start_invitmax_tol * max( bracket ends).
- use_initial_multiplierbool
ignore initial_multiplier?.
- initialize_approx_eigenvectorbool
should a suitable initial eigenvector should be chosen or should a previous eigenvector may be used?.
- space_criticalbool
if
space_critical
is True, every effort will be made to use as little space as possible. This may result in longer computation time.- deallocate_error_fatalbool
if
deallocate_error_fatal
is True, any array/pointer deallocation error will terminate execution. Otherwise, computation will continue.- problem_filestr
name of file into which to write problem data.
- symmetric_linear_solverstr
symmetric (indefinite) linear equation solver.
- definite_linear_solverstr
definite linear equation solver.
- prefixstr
all output lines will be prefixed by the string contained in quotes within
prefix
, e.g. ‘word’ (note the qutoes) will result in the prefix word.- sls_optionsdict
default control options for SLS (see
sls.initialize
).- ir_optionsdict
default control options for IR (see
ir.initialize
).
- rqs.load(n, H_type, H_ne, H_row, H_col, H_ptr, options=None)#
Import problem data into internal storage prior to solution.
Parameters:
- nint
holds the number of variables.
- H_typestring
specifies the symmetric storage scheme used for the Hessian
. It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, ‘diagonal’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’; lower or upper case variants are allowed. - H_neint
holds the number of entries in the lower triangular part of
in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. - H_rowndarray(H_ne)
holds the row indices of the lower triangular part of
in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes, and in this case can be None. - H_colndarray(H_ne)
holds the column indices of the lower triangular part of
in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the other storage schemes are used, and in this case can be None. - H_ptrndarray(n+1)
holds the starting position of each row of the lower triangular part of
, 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 None. - optionsdict, optional
dictionary of control options (see
rqs.initialize
).
- rqs.load_m(n, M_type, M_ne, M_row, M_col, M_ptr, options=None)#
Import problem data for the scaling matrix
, if needed, into internal storage prior to solution. Parameters:
- nint
holds the number of variables.
- M_typestring
specifies the symmetric storage scheme used for the Hessian
. It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, ‘diagonal’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’; lower or upper case variants are allowed. - M_neint
holds the number of entries in the lower triangular part of
in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. - M_rowndarray(M_ne)
holds the row indices of the lower triangular part of
in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes, and in this case can be None. - M_colndarray(M_ne)
holds the column indices of the lower triangular part of
in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the other storage schemes are used, and in this case can be None. - M_ptrndarray(n+1)
holds the starting position of each row of the lower triangular part of
, 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 None. - optionsdict, optional
dictionary of control options (see
rqs.initialize
).
- rqs.load_a(m, A_type, A_ne, A_row, A_col, A_ptr, options=None)#
Import problem data for the constraint matrix
, if needed, into internal storage prior to solution. Parameters:
- mint
holds the number of constraints.
- A_typestring
specifies the unsymmetric storage scheme used for the Hessian
. It should be one of ‘coordinate’, ‘sparse_by_rows’ or ‘dense’; lower or upper case variants are allowed. - A_neint
holds the number of entries in
in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. - A_rowndarray(A_ne)
holds the row indices of
in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes, and in this case can be None. - A_colndarray(A_ne)
holds the column indices of
in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the other storage schemes are used, and in this case can be None. - A_ptrndarray(m+1)
holds the starting position of
, 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 None. - optionsdict, optional
dictionary of control options (see
rqs.initialize
).
- rqs.solve_problem(n, power, weight, f, g, H_ne, H_val, M_ne, M_val, m, A_ne, A_val)#
Find the global minimizer of the regularized quadratic objective function
subject to the affine constraints. Parameters:
- nint
holds the number of variables.
- powerfloat
holds the regularization power,
. - weightfloat
holds the strictly positive regularization weight,
. - ffloat
holds the constant term
in the objective function. - gndarray(n)
holds the values of the linear term
in the objective function. - H_neint
holds the number of entries in the lower triangular part of the Hessian
. - H_valndarray(H_ne)
holds the values of the nonzeros in the lower triangle of the Hessian
in the same order as specified in the sparsity pattern in rqs.load
.- M_neint
holds the number of entries in the lower triangular part of the scaling matrix
if it is not the identity matrix. Otherwise it should be None. - M_valndarray(M_ne)
holds the values of the nonzeros in the lower triangle of the scaling matrix
in the same order as specified in the sparsity pattern in rqs.load_m
if needed. Otherwise it should be None.- mint
holds the number of constraints.
- A_neint
holds the number of entries in the lower triangular part of the constraint matrix
if . Otherwise it should be None. - A_valndarray(A_ne)
holds the values of the nonzeros in the lower triangle of the constraint matrix
in the same order as specified in the sparsity pattern in rqs.load_a
if needed. Otherwise it should be None.Returns:
- xndarray(n)
holds the values of the approximate minimizer
after a successful call. - yndarray(m)
holds the values of the Lagrange multipliers associated with the affine constraints, if any. Absent if
trs.load_a
has not been called.
- [optional] rqs.information()
Provide optional output information
Returns:
- informdict
- dictionary containing output information:
- statusint
return status. Possible values are:
0
The run was successful.
-1
An allocation error occurred. A message indicating the offending array is written on unit options[‘error’], and the returned allocation status and a string containing the name of the offending array are held in inform[‘alloc_status’] and inform[‘bad_alloc’] respectively.
-2
A deallocation error occurred. A message indicating the offending array is written on unit options[‘error’] and the returned allocation status and a string containing the name of the offending array are held in inform[‘alloc_status’] and inform[‘bad_alloc’] respectively.
-3
The restriction n > 0, m > 0, weight > 0 or power
2, or requirement that type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, ‘diagonal’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’ has been violated. -9
The analysis phase of the factorization failed; the return status from the factorization package is given by inform[‘factor_status’].
-10
The factorization failed; the return status from the factorization package is given by inform[‘factor_status’].
-11
The solution of a set of linear equations using factors from the factorization package failed; the return status from the factorization package is given by inform[‘factor_status’].
-15
does not appear to be strictly diagonally dominant -16
The problem is so ill-conditioned that further progress is impossible.
-19
The CPU time limit has been reached. This may happen if options[‘cpu_time_limit’] is too small, but may also be symptomatic of a badly scaled problem.
-23
An entry from the strict upper triangle of
has been specified. - alloc_statusint
the status of the last attempted allocation/deallocation.
- bad_allocstr
the name of the array for which an allocation/deallocation error occurred.
- factorizationsint
the number of factorizations performed.
- max_entries_factorslong
the maximum number of entries in the factors.
- len_historyint
the number of
pairs in the history. - objfloat
the value of the quadratic function.
- obj_regularizedfloat
the value of the regularized quadratic function.
- x_normfloat
the
-norm of , . - multiplierfloat
the Lagrange multiplier corresponding to the regularization.
- polefloat
a lower bound
, where is the left-most eigenvalue of . - dense_factorizationbool
was a dense factorization used?.
- hard_casebool
has the hard case occurred?.
- timedict
- dictionary containing timing information:
- totalfloat
total CPU time spent in the package.
- assemblefloat
CPU time spent building
. - analysefloat
CPU time spent reordering
prior to factorization. - factorizefloat
CPU time spent factorizing
. - solvefloat
CPU time spent solving linear systems inolving
. - clock_totalfloat
total clock time spent in the package.
- clock_assemblefloat
clock time spent building
. - clock_analysefloat
clock time spent reordering
prior to factorization. - clock_factorizefloat
clock time spent factorizing
. - clock_solvefloat
clock time spent solving linear systems inolving
. - historydict
- dictionary recording the history of the iterates:
- lambdandarray(100)
the values of
for the first min(100, len_history
) iterations.- x_normndarray(100)
the corresponding values of
. - sls_informdict
inform parameters for SLS (see
sls.information
).- ir_informdict
inform parameters for IR (see
ir.information
).
- rqs.terminate()#
Deallocate all internal private storage.
example code#
from galahad import rqs
import numpy as np
np.set_printoptions(precision=4,suppress=True,floatmode='fixed')
print("\n** python test: rqs")
# set parameters
p = 1.0
n = 3
m = 1
infinity = float("inf")
# describe objective function
f = 0.96
g = np.array([0.0,2.0,0.0])
H_type = 'coordinate'
H_ne = 4
H_row = np.array([0,1,2,2])
H_col = np.array([0,1,2,0])
H_ptr = None
H_val = np.array([1.0,2.0,3.0,4.0])
# describe norm
M_type = 'coordinate'
M_ne = 3
M_row = np.array([0,1,2])
M_col = np.array([0,1,2])
M_ptr = None
M_val = np.array([1.0,2.0,1.0])
# describe constraint
A_type = 'coordinate'
A_ne = 3
A_row = np.array([0,0,0])
A_col = np.array([0,1,2])
A_ptr = None
A_val = np.array([1.0,1.0,1.0])
# set regularization parameters
power = 3.0
weight = 1.0
# allocate internal data and set default options
options = rqs.initialize()
# set some non-default options
options['print_level'] = 0
#print("options:", options)
# load data (and optionally non-default options)
rqs.load(n, H_type, H_ne, H_row, H_col, H_ptr, options)
# find minimum of quadratic within the trust region
print("\n solve problem 1")
x = rqs.solve_problem(n, power, weight, f, g, H_ne, H_val)
print(" x:",x)
# get information
inform = rqs.information()
print(" f: %.4f" % inform['obj'])
# load data (and optionally non-default options)
rqs.load_m(n, M_type, M_ne, M_row, M_col, M_ptr)
# find minimum of quadratic within the trust region
print("\n solve problem 2 with additional non-unit norm")
x = rqs.solve_problem(n, power, weight, f, g, H_ne, H_val, M_ne, M_val)
print(" x:",x)
# get information
inform = rqs.information()
print(" f: %.4f" % inform['obj'])
# load data (and optionally non-default options)
rqs.load_a(m, A_type, A_ne, A_row, A_col, A_ptr)
# find minimum of quadratic within the trust region
print("\n solve problem 3 with additional linear constraint")
x, y = rqs.solve_problem(n, power, weight, f, g, H_ne, H_val,
M_ne, M_val, m, A_ne, A_val)
print(" x:",x)
print(" y:",y)
# get information
inform = rqs.information()
print(" f: %.4f" % inform['obj'])
print('** rqs exit status:', inform['status'])
# deallocate internal data
rqs.terminate()
This example code is available in $GALAHAD/src/rqs/Python/test_rqs.py .