GALAHAD RQS package#

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

\[r(x) = f + g^T x + \frac{1}{2} x^T H x + \frac{\sigma}{p} \|x\|_{M}^p,\]
where the weight \(\sigma > 0\), the power \(p \geq 2\), the vector \(x\) may optionally be required to satisfy affine constraints \(A x = 0,\) and where the \(M\)-norm of \(x\) is defined to be \(\|x\|_{M} = \sqrt{x^T M x}\).

The matrix \(M\) need not be provided in the commonly-occurring \(\ell_2\)-regularization case for which \(M = I\), the \(n\) by \(n\) identity matrix.

Factorization of matrices of the form \(H + \lambda M\), or

\[\begin{split}\left(\begin{array}{cc} H + \lambda M & A^T \\ A & 0 \end{array}\right)\end{split}\]
in cases where \(A x = 0\) is imposed, for a succession of scalars \(\lambda\) will be required, so this package is most suited for the case where such a factorization may be found efficiently. If this is not the case, the package 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 \(x_*\) necessarily satisfies the optimality condition \(H x_* + \lambda_* M x_* + A^T y_* + g = 0\) and \(A x_* = 0\), where \(\lambda_* = \sigma \|x_*\|^{p-2}\) is a Lagrange multiplier corresponding to the regularization, and \(y_*\) are Lagrange multipliers for the linear constraints \(A x = 0\), if any. In addition in all cases, the matrix \(H + \lambda_* M\) will be positive semi-definite on the null-space of \(A\); in most instances it will actually be positive definite, but in special “hard” cases singularity is a possibility.

The method is iterative, and proceeds in two phases. Firstly, lower and upper bounds, \(\lambda_L\) and \(\lambda_U\), on \(\lambda_*\) are computed using Gershgorin’s theorems and other eigenvalue bounds. The first phase of the computation proceeds by progressively shrinking the bound interval \([\lambda_L,\lambda_U]\) until a value \(\lambda\) for which \(\|x(\lambda)\|_{M} \geq \sigma \|x(\lambda)\|_M^{p-2}\) is found. Here \(x(\lambda)\) and its companion \(y(\lambda)\) are defined to be a solution of

\[(H + \lambda M)x(\lambda) + A^T y(\lambda) = - g \;\;\mbox{and}\;\; A x(\lambda) = 0.\;\;\mbox{(2)}\]
Once the terminating \(\lambda\) from the first phase has been discovered, the second phase consists of applying Newton or higher-order iterations to the nonlinear “secular” equation \(\lambda = \sigma \|x(\lambda)\|_M^{p-2}\) with the knowledge that such iterations are both globally and ultimately rapidly convergent. It is possible in the “hard” case that the interval in the first-phase will shrink to the single point \(\lambda_*\), and precautions are taken, using inverse iteration with Rayleigh-quotient acceleration to ensure that this too happens rapidly.

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 \(H + \lambda M\) is used. If \(A x = 0\) is required, a sparse symmetric, indefinite factorization of (1) is used rather than a Cholesky factorization.

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#

unsymmetric storage#

The unsymmetric \(m\) by \(n\) matrix \(A\), if it is needed, 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 \(1 \leq i \leq m\), \(1 \leq j \leq n\). 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 \(1 \leq i \leq m\), \(1 \leq j \leq n\). 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, \(1 \leq l \leq ne\), of \(A\), its row index i, column index j and value \(A_{ij}\), \(1 \leq i \leq m\), \(1 \leq j \leq n\), 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+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. 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+1) holds the total number of entries plus one. The row indices i, \(1 \leq i \leq m\), 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, \(1 \leq j \leq n\), 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.

symmetric storage#

The symmetric \(n\) by \(n\) matrices \(H\) and, optionally. \(M\) 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 \(1 \leq j \leq i \leq n\)) need be held. In this case the lower triangle should be stored by rows, that is component \((i-1) * i / 2 + j\) of the storage array H_val will hold the value \(H_{ij}\) (and, by symmetry, \(H_{ji}\)) for \(1 \leq j \leq i \leq n\). 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, \(1 \leq l \leq ne\), of \(H\), its row index i, column index j and value \(H_{ij}\), \(1 \leq j \leq i \leq n\), 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+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. 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 \(1 \leq i \neq j \leq n\)) only the diagonals entries \(H_{ii}\), \(1 \leq i \leq n\) 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 rqs package must be called in the following order:

  • rqs_initialize - provide default control parameters and set up initial data structures

  • rqs_read_specfile (optional) - override control values by reading replacement values from a file

  • rqs_import - set up problem data structures and fixed values

  • rqs_import_m - (optional) set up problem data structures and fixed values for the scaling matrix \(M\), if any

  • rqs_import_a - (optional) set up problem data structures and fixed values for the constraint matrix \(A\), if any

  • rqs_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved

  • rqs_solve_problem - solve the regularised quadratic problem

  • rqs_information (optional) - recover information about the solution and solution process

  • rqs_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 rqs_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 rqs_control_type)

status

is a scalar variable of type Int32 that gives the exit status from the package. Possible values are (currently):

  • 0

    The initialization was successful.

    function rqs_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/rqs/RQS.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/rqs.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 rqs_control_type)

specfile

is a one-dimensional array of type Vararg{Cchar} that must give the name of the specification file

    function rqs_import(T, control, data, status, n,
                        H_type, H_ne, H_row, H_col, H_ptr)

Import problem data into internal storage prior to solution.

Parameters:

control

is a structure whose members provide control parameters for the remaining procedures (see rqs_control_type)

data

holds private internal data

status

is a scalar variable of type Int32 that gives the exit status from the package. Possible values are:

  • 0

    The import was successful

  • -1

    An allocation error occurred. A message indicating the offending array is written on unit control.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 control.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 restrictions n > 0 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, diagonal’ or ‘identity’ has been violated.

n

is a scalar variable of type Int32 that holds the number of rows (and columns) of H.

H_type

is a one-dimensional array of type Vararg{Cchar} that specifies the symmetric storage scheme used for the Hessian, \(H\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, or ‘diagonal’; lower or upper case variants are allowed.

H_ne

is a scalar variable of type Int32 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 Int32 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 C_NULL.

H_col

is a one-dimensional array of size H_ne and type Int32 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 or diagonal storage schemes are used, and in this case can be C_NULL.

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\), 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 C_NULL.

    function rqs_import_m(T, data, status, n, M_type, M_ne, M_row, M_col, M_ptr)

Import data for the scaling matrix M into internal storage prior to solution.

Parameters:

data

holds private internal data

status

is a scalar variable of type Int32 that gives the exit status from the package. Possible values are:

  • 0

    The import was successful

  • -1

    An allocation error occurred. A message indicating the offending array is written on unit control.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 control.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 restrictions n > 0 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, diagonal’ or ‘identity’ has been violated.

n

is a scalar variable of type Int32 that holds the number of rows (and columns) of M.

M_type

is a one-dimensional array of type Vararg{Cchar} that specifies the symmetric storage scheme used for the scaling matrix, \(M\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, or ‘diagonal’; lower or upper case variants are allowed.

M_ne

is a scalar variable of type Int32 that holds the number of entries in the lower triangular part of \(M\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes.

M_row

is a one-dimensional array of size M_ne and type Int32 that holds the row indices of the lower triangular part of \(M\) 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 C_NULL.

M_col

is a one-dimensional array of size M_ne and type Int32 that holds the column indices of the lower triangular part of \(M\) in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the dense, diagonal or identity storage schemes are used, and in this case can be C_NULL.

M_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 \(M\), 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 C_NULL.

    function rqs_import_a(T, data, status, m, A_type, A_ne, A_row, A_col, A_ptr)

Import data for the constraint matrix A into internal storage prior to solution.

Parameters:

data

holds private internal data

status

is a scalar variable of type Int32 that gives the exit status from the package. Possible values are:

  • 0

    The import was successful

  • -1

    An allocation error occurred. A message indicating the offending array is written on unit control.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 control.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 restrictions n > 0 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’ or ‘sparse_by_rows’ has been violated.

m

is a scalar variable of type Int32 that holds the number of general linear constraints, i.e., the number of rows of A, if any. m must be non-negative.

A_type

is a one-dimensional array of type Vararg{Cchar} that specifies the unsymmetric storage scheme used for the constraint Jacobian, \(A\) if any. 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 Int32 that holds the number of entries in \(A\), if used, 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 Int32 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 C_NULL.

A_col

is a one-dimensional array of size A_ne and type Int32 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 C_NULL.

A_ptr

is a one-dimensional array of size n+1 and type Int32 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 C_NULL.

    function rqs_reset_control(T, control, data, status)

Reset control parameters after import if required.

Parameters:

control

is a structure whose members provide control parameters for the remaining procedures (see rqs_control_type)

data

holds private internal data

status

is a scalar variable of type Int32 that gives the exit status from the package. Possible values are:

  • 0

    The import was successful.

    function rqs_solve_problem(T, data, status, n, power, weight, f, c,
                              H_ne, H_val, x, M_ne, M_val, m, A_ne, A_val, y)

Solve the regularised quadratic problem.

Parameters:

data

holds private internal data

status

is a scalar variable of type Int32 that gives the entry and exit status from the package.

On initial entry, status must be set to 1.

Possible exit values are:

  • 0

    The run was successful.

  • -1

    An allocation error occurred. A message indicating the offending array is written on unit control.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 control.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 restrictions n > 0, power > 2, weight > 0 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, ‘diagonal’ or ‘identity’ has been violated.

  • -9

    The analysis phase of the factorization of the matrix (1) failed.

  • -10

    The factorization of the matrix (1) failed.

  • -15

    The matrix M appears not to be diagonally dominant.

  • -16

    The problem is so ill-conditioned that further progress is impossible.

  • -18

    Too many factorizations have been required. This may happen if control.max factorizations is too small, but may also be symptomatic of a badly scaled problem.

n

is a scalar variable of type Int32 that holds the number of variables.

power

is a scalar of type T that holds the order of regularisation, \(p\), used. power must be no smaller than 2.

weight

is a scalar of type T that holds the regularisation weight, \(\sigma\), used. weight must be strictly positive.

c

is a one-dimensional array of size n and type T that holds the linear term \(c\) of the objective function. The j-th component of c, j = 1, … , m, contains \(c_j\).

f

is a scalar of type T that holds the constant term \(f\) of the objective function.

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 size h_ne and type T that holds the values of the entries of the lower triangular part of the Hessian matrix \(H\) in any of the available storage schemes.

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 x, j = 1, … , n, contains \(x_j\).

M_ne

is a scalar variable of type Int32 that holds the number of entries in the scaling matrix \(M\) if it not the identity matrix.

M_val

is a one-dimensional array of size M_ne and type T that holds the values of the entries of the scaling matrix \(M\), if it is not the identity matrix, in any of the available storage schemes. If M_val is C_NULL, M will be taken to be the identity matrix.

m

is a scalar variable of type Int32 that holds the number of general linear constraints, if any. m must be non-negative.

A_ne

is a scalar variable of type Int32 that holds the number of entries in the constraint Jacobian matrix \(A\) if used. A_ne must be non-negative.

A_val

is a one-dimensional array of size A_ne and type T that holds the values of the entries of the constraint Jacobian matrix \(A\), if used, in any of the available storage schemes. If A_val is C_NULL, no constraints will be enforced.

y

is a one-dimensional array of size m and type T that holds the values \(y\) of the Lagrange multipliers for the equality constraints \(A x = 0\) if used. The i-th component of y, i = 1, … , m, contains \(y_i\).

    function rqs_information(T, data, inform, status)

Provides output information

Parameters:

data

holds private internal data

inform

is a structure containing output information (see rqs_inform_type)

status

is a scalar variable of type Int32 that gives the exit status from the package. Possible values are (currently):

  • 0

    The values were recorded successfully

    function rqs_terminate(T, data, control, inform)

Deallocate all internal private storage

Parameters:

data

holds private internal data

control

is a structure containing control information (see rqs_control_type)

inform

is a structure containing output information (see rqs_inform_type)

available structures#

rqs_control_type structure#

    struct rqs_control_type{T}
      f_indexing::Bool
      error::Int32
      out::Int32
      problem::Int32
      print_level::Int32
      dense_factorization::Int32
      new_h::Int32
      new_m::Int32
      new_a::Int32
      max_factorizations::Int32
      inverse_itmax::Int32
      taylor_max_degree::Int32
      initial_multiplier::T
      lower::T
      upper::T
      stop_normal::T
      stop_hard::T
      start_invit_tol::T
      start_invitmax_tol::T
      use_initial_multiplier::Bool
      initialize_approx_eigenvector::Bool
      space_critical::Bool
      deallocate_error_fatal::Bool
      problem_file::NTuple{31,Cchar}
      symmetric_linear_solver::NTuple{31,Cchar}
      definite_linear_solver::NTuple{31,Cchar}
      prefix::NTuple{31,Cchar}
      sls_control::sls_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

unit for error messages

Int32 out

unit for monitor output

Int32 problem

unit to write problem data into file problem_file

Int32 print_level

controls level of diagnostic output

Int32 dense_factorization

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

Int32 new_h

how much of \(H\) has changed since the previous call. Possible values are

  • 0 unchanged

  • 1 values but not indices have changed

  • 2 values and indices have changed

Int32 new_m

how much of \(M\) has changed since the previous call. Possible values are

  • 0 unchanged

  • 1 values but not indices have changed

  • 2 values and indices have changed

Int32 new_a

how much of \(A\) has changed since the previous call. Possible values are 0 unchanged 1 values but not indices have changed 2 values and indices have changed

Int32 max_factorizations

the maximum number of factorizations (=iterations) allowed. -ve implies no limit

Int32 inverse_itmax

the number of inverse iterations performed in the “maybe hard” case

Int32 taylor_max_degree

maximum degree of Taylor approximant allowed

T initial_multiplier

initial estimate of the Lagrange multipler

T lower

lower and upper bounds on the multiplier, if known

T upper

see lower

T stop_normal

stop when \(| \|x\| - (multiplier/\sigma)^(1/(p-2)) | \leq\) stop_normal \* max \(( \|x\|, (multiplier/\sigma)^(1/(p-2)) )\) REAL ( KIND = wp ) :: stop_normal = epsmch \*\* 0.75

T stop_hard

stop when bracket on optimal multiplier <= stop_hard * max( bracket ends ) REAL ( KIND = wp ) :: stop_hard = epsmch ** 0.75

T start_invit_tol

start inverse iteration when bracket on optimal multiplier <= stop_start_invit_tol * max( bracket ends )

T start_invitmax_tol

start full inverse iteration when bracket on multiplier <= stop_start_invitmax_tol * max( bracket ends)

Bool use_initial_multiplier

ignore initial_multiplier?

Bool initialize_approx_eigenvector

should a suitable initial eigenvector should be chosen or should a previous eigenvector may be used?

Bool space_critical

if space is critical, ensure allocated arrays are no bigger than needed

Bool deallocate_error_fatal

exit if any deallocation fails

char problem_file[31]

name of file into which to write problem data

char symmetric_linear_solver[31]

symmetric (indefinite) linear equation solver

char definite_linear_solver[31]

definite 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 the Cholesky factorization and solution (see sls_c documentation)

struct ir_control_type ir_control

control parameters for iterative refinement (see ir_c documentation)

rqs_time_type structure#

    struct rqs_time_type{T}
      total::T
      assemble::T
      analyse::T
      factorize::T
      solve::T
      clock_total::T
      clock_assemble::T
      clock_analyse::T
      clock_factorize::T
      clock_solve::T

detailed documentation#

time derived type as a Julia structure

components#

T total

total CPU time spent in the package

T assemble

CPU time spent building \(H + \lambda M\).

T analyse

CPU time spent reordering \(H + \lambda M\) prior to factorization.

T factorize

CPU time spent factorizing \(H + \lambda M\).

T solve

CPU time spent solving linear systems inolving \(H + \lambda M\).

T clock_total

total clock time spent in the package

T clock_assemble

clock time spent building \(H + \lambda M\)

T clock_analyse

clock time spent reordering \(H + \lambda M\) prior to factorization

T clock_factorize

clock time spent factorizing \(H + \lambda M\)

T clock_solve

clock time spent solving linear systems inolving \(H + \lambda M\)

rqs_history_type structure#

    struct rqs_history_type{T}
      lambda::T
      x_norm::T

detailed documentation#

history derived type as a Julia structure

components#

T lambda

the value of \(\lambda\)

T x_norm

the corresponding value of \(\|x(\lambda)\|_M\)

rqs_inform_type structure#

    struct rqs_inform_type{T}
      status::Int32
      alloc_status::Int32
      factorizations::Int32
      max_entries_factors::Int64
      len_history::Int32
      obj::T
      obj_regularized::T
      x_norm::T
      multiplier::T
      pole::T
      dense_factorization::Bool
      hard_case::Bool
      bad_alloc::NTuple{81,Cchar}
      time::rqs_time_type{T}
      history::NTuple{100,rqs_history_type{T}}
      sls_inform::sls_inform_type{T}
      ir_inform::ir_inform_type{T}

detailed documentation#

inform derived type as a Julia structure

components#

Int32 status

reported return status:

  • 0

    the solution has been found

  • -1

    an array allocation has failed

  • -2

    an array deallocation has failed

  • -3

    n and/or sigma is not positive and/or p <= 2

  • -9

    the analysis phase of the factorization of \(H + \lambda M\) failed

  • -10

    the factorization of \(H + \lambda M\) failed

  • -15

    \(M\) does not appear to be strictly diagonally dominant

  • -16

    ill-conditioning has prevented furthr progress

Int32 alloc_status

STAT value after allocate failure.

Int32 factorizations

the number of factorizations performed

Int64 max_entries_factors

the maximum number of entries in the factors

Int32 len_history

the number of \((\|x\|_M,\lambda)\) pairs in the history

T obj

the value of the quadratic function

T obj_regularized

the value of the regularized quadratic function

T x_norm

the \(M\) -norm of \(x\), \(\|x\|_M\)

T multiplier

the Lagrange multiplier corresponding to the regularization

T pole

a lower bound max \((0,-\lambda_1)\), where \(\lambda_1\) is the left-most eigenvalue of \((H,M)\)

Bool dense_factorization

was a dense factorization used?

Bool hard_case

has the hard case occurred?

NTuple{81,Cchar} bad_alloc

name of array which provoked an allocate failure

struct rqs_time_type time

time information

struct rqs_history_type history[100]

history information

struct sls_inform_type sls_inform

cholesky information (see sls_c documentation)

struct ir_inform_type ir_inform

iterative_refinement information (see ir_c documentation)

example calls#

This is an example of how to use the package to solve a regularized quadratic subproblem; the code is available in $GALAHAD/src/rqs/Julia/test_rqs.jl . A variety of supported Hessian and constraint matrix storage formats are shown.

# test_rqs.jl
# Simple code to test the Julia interface to RQS

using GALAHAD
using Test
using Printf
using Accessors

function test_rqs(::Type{T}) where T
  # Derived types
  data = Ref{Ptr{Cvoid}}()
  control = Ref{rqs_control_type{T}}()
  inform = Ref{rqs_inform_type{T}}()

  # Set problem data
  n = 3 # dimension of H
  m = 1 # dimension of A
  H_ne = 4 # number of elements of H
  M_ne = 3 # number of elements of M
  A_ne = 3 # number of elements of A
  H_dense_ne = 6 # number of elements of H
  M_dense_ne = 6 # number of elements of M
  H_row = Cint[1, 2, 3, 3]  # row indices, NB lower triangle
  H_col = Cint[1, 2, 3, 1]
  H_ptr = Cint[1, 2, 3, 5]
  M_row = Cint[1, 2, 3]  # row indices, NB lower triangle
  M_col = Cint[1, 2, 3]
  M_ptr = Cint[1, 2, 3, 4]
  A_row = Cint[1, 1, 1]
  A_col = Cint[1, 2, 3]
  A_ptr = Cint[1, 4]
  H_val = T[1.0, 2.0, 3.0, 4.0]
  M_val = T[1.0, 2.0, 1.0]
  A_val = T[1.0, 1.0, 1.0]
  H_dense = T[1.0, 0.0, 2.0, 4.0, 0.0, 3.0]
  M_dense = T[1.0, 0.0, 2.0, 0.0, 0.0, 1.0]
  H_diag = T[1.0, 0.0, 2.0]
  M_diag = T[1.0, 2.0, 1.0]
  f = 0.96
  power = 3.0
  weight = 1.0
  c = T[0.0, 2.0, 0.0]

  st = ' '
  status = Ref{Cint}()
  x = zeros(T, n)

  @printf(" Fortran sparse matrix indexing\n\n")
  @printf(" basic tests of storage formats\n\n")

  for a_is in 0:1 # add a linear constraint?
    for m_is in 0:1 # include a scaling matrix?
      if (a_is == 1) && (m_is == 1)
        ma = "MA"
      elseif a_is == 1
        ma = "A"
      elseif m_is == 1
        ma = "M"
      else
        ma = ""
      end

      for storage_type in 1:4
        # Initialize RQS
        rqs_initialize(T, data, control, status)

        # Set user-defined control options
        @reset control[].f_indexing = true # fortran sparse matrix indexing

        # sparse co-ordinate storage
        if storage_type == 1
          st = 'C'
          # import the control parameters and structural data
          rqs_import(T, control, data, status, n,
                     "coordinate", H_ne, H_row, H_col, C_NULL)
          if m_is == 1
            rqs_import_m(T, data, status, n,
                         "coordinate", M_ne, M_row, M_col, C_NULL)
          end
          if a_is == 1
            rqs_import_a(T, data, status, m,
                         "coordinate", A_ne, A_row, A_col, C_NULL)
          end

          # solve the problem
          if (a_is == 1) && (m_is == 1)
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              M_ne, M_val, m, A_ne, A_val, C_NULL)
          elseif a_is == 1
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              0, C_NULL, m, A_ne, A_val, C_NULL)
          elseif m_is == 1
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              M_ne, M_val, 0, 0, C_NULL, C_NULL)
          else
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        # sparse by rows
        if storage_type == 2
          st = 'R'
          # import the control parameters and structural data
          rqs_import(T, control, data, status, n,
                     "sparse_by_rows", H_ne, C_NULL, H_col, H_ptr)
          if m_is == 1
            rqs_import_m(T, data, status, n,
                         "sparse_by_rows", M_ne, C_NULL, M_col, M_ptr)
          end
          if a_is == 1
            rqs_import_a(T, data, status, m,
                         "sparse_by_rows", A_ne, C_NULL, A_col, A_ptr)
          end

          # solve the problem
          if (a_is == 1) && (m_is == 1)
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              M_ne, M_val, m, A_ne, A_val, C_NULL)
          elseif a_is == 1
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              0, C_NULL, m, A_ne, A_val, C_NULL)
          elseif m_is == 1
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              M_ne, M_val, 0, 0, C_NULL, C_NULL)
          else
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, H_ne, H_val, x,
                              0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        # dense
        if storage_type == 3
          st = 'D'
          # import the control parameters and structural data
          rqs_import(T, control, data, status, n,
                     "dense", H_ne, C_NULL, C_NULL, C_NULL)
          if m_is == 1
            rqs_import_m(T, data, status, n,
                         "dense", M_ne, C_NULL, C_NULL, C_NULL)
          end
          if a_is == 1
            rqs_import_a(T, data, status, m,
                         "dense", A_ne, C_NULL, C_NULL, C_NULL)
          end

          # solve the problem
          if (a_is == 1) && (m_is == 1)
            rqs_solve_problem(T, data, status, n, power, weight,
                              f, c, H_dense_ne, H_dense, x,
                              M_dense_ne, M_dense, m, A_ne, A_val,
                              C_NULL)
          elseif a_is == 1
            rqs_solve_problem(T, data, status, n, power, weight,
                              f, c, H_dense_ne, H_dense, x,
                              0, C_NULL, m, A_ne, A_val, C_NULL)
          elseif m_is == 1
            rqs_solve_problem(T, data, status, n, power, weight,
                              f, c, H_dense_ne, H_dense, x,
                              M_dense_ne, M_dense, 0, 0, C_NULL, C_NULL)
          else
            rqs_solve_problem(T, data, status, n, power, weight,
                              f, c, H_dense_ne, H_dense, x,
                              0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        # diagonal
        if storage_type == 4
          st = 'L'
          # import the control parameters and structural data
          rqs_import(T, control, data, status, n,
                     "diagonal", H_ne, C_NULL, C_NULL, C_NULL)
          if m_is == 1
            rqs_import_m(T, data, status, n,
                         "diagonal", M_ne, C_NULL, C_NULL, C_NULL)
          end
          if a_is == 1
            rqs_import_a(T, data, status, m,
                         "dense", A_ne, C_NULL, C_NULL, C_NULL)
          end

          # solve the problem
          if (a_is == 1) && (m_is == 1)
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, n, H_diag, x,
                              n, M_diag, m, A_ne, A_val, C_NULL)
          elseif a_is == 1
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, n, H_diag, x,
                              0, C_NULL, m, A_ne, A_val, C_NULL)
          elseif m_is == 1
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, n, H_diag, x,
                              n, M_diag, 0, 0, C_NULL, C_NULL)
          else
            rqs_solve_problem(T, data, status, n,
                              power, weight, f, c, n, H_diag, x,
                              0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        rqs_information(T, data, inform, status)

        @printf("format %c%s: RQS_solve_problem exit status = %1i, f = %.2f\n", st, ma,
                inform[].status, inform[].obj_regularized)

        # @printf("x: ")
        # for i in 1:n+m
        #   @printf("%f ", x[i])
        # end

        # Delete internal workspace
        rqs_terminate(T, data, control, inform)
      end
    end
  end
  return 0
end

@testset "RQS" begin
  @test test_rqs(Float32) == 0
  @test test_rqs(Float64) == 0
end