GALAHAD SBLS package#

purpose#

Given a block, symmetric matrix

\[\begin{split}K_{H} = \begin{pmatrix}H & A^T \\ A & - C\end{pmatrix},\end{split}\]
the sbls package constructs a variety of preconditioners of the form
\[\begin{split}K_{G} = \begin{pmatrix}G & A^T \\ A & - C\end{pmatrix}.\end{split}\]
Here, the leading-block matrix \(G\) is a suitably-chosen approximation to \(H\); it may either be prescribed explicitly, in which case a symmetric indefinite factorization of \(K_{G}\) will be formed using the GALAHAD package SLS, or implicitly, by requiring certain sub-blocks of \(G\) be zero. In the latter case, a factorization of \(K_{G}\) will be obtained implicitly (and more efficiently) without recourse to sls. In particular, for implicit preconditioners, a reordering
\[\begin{split}K_{G} = P \begin{pmatrix} G_{11}^{} & G_{21}^T & A_1^T \\ G_{21}^{} & G_{22}^{} & A_2^T \\ A_{1}^{} & A_{2}^{} & - C \end{pmatrix} P^T\end{split}\]
involving a suitable permutation \(P\) will be found, for some invertible sub-block (“basis”) \(A_1\) of the columns of \(A\); the selection and factorization of \(A_1\) uses the package uls. Once the preconditioner has been constructed, solutions to the preconditioning system
\[\begin{split}\begin{pmatrix}G & A^T \\ A & - C\end{pmatrix} \begin{pmatrix}x \\ y\end{pmatrix} = \begin{pmatrix}a \\ b\end{pmatrix}\end{split}\]
may be obtained by the package. Full advantage is taken of any zero coefficients in the matrices \(H\), \(A\) and \(C\).

See Section 4 of $GALAHAD/doc/sbls.pdf for additional details.

method#

The method used depends on whether an explicit or implicit factorization is required. In the explicit case, the package is really little more than a wrapper for the symmetric, indefinite linear solver SLS in which the system matrix \(K_G\) is assembled from its constituents \(A\), \(C\) and whichever \(G\) is requested by the user. Implicit-factorization preconditioners are more involved, and there is a large variety of different possibilities. The essential ideas 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.

The range-space factorization is based upon the decomposition

\[\begin{split}K_G = \begin{pmatrix}G & 0 \\ A & I\end{pmatrix} \begin{pmatrix}G^{-1} & 0 \\ 0 & - S\end{pmatrix} \begin{pmatrix}G & A^T \\ 0 & I\end{pmatrix},\end{split}\]
where the ``Schur complement’’ \(S = C + A G^{-1} A^T\). Such a method requires that \(S\) is easily invertible. This is often the case when \(G\) is a diagonal matrix, in which case \(S\) is frequently sparse, or when \(m \ll n\) in which case \(S\) is small and a dense Cholesky factorization may be used.

When \(C = 0\), the null-space factorization is based upon the decomposition

\[\begin{split}K_{G} = P \begin{pmatrix} G_{11}^{} & 0 & I \\ G_{21}^{} & I & A_{2}^{T} A_{1}^{-T} \\ A_{1}^{} & 0 & 0 \end{pmatrix} \begin{pmatrix} 0 & 0 & I \\ \;\;\; 0 \;\; & \;\; R \;\; & 0 \\ I & 0 & - G_{11}^{} \end{pmatrix} \begin{pmatrix} G_{11}^{} & G_{21}^T & A_{1}^T \\ 0 & I & 0 \\ I & A_{1}^{-1} A_{2}^{} & 0 \end{pmatrix} P^T,\end{split}\]
where the ``reduced Hessian’’
\[\begin{split}R = ( - A_{2}^{T} A_1^{-T} \;\; I ) \begin{pmatrix}G_{11}^{} & G_{21}^{T} \\ G_{21}^{} & G_{22}^{}\end{pmatrix} \begin{pmatrix}- A_1^{-1} A_2^{} \\ I\end{pmatrix}\end{split}\]
and \(P\) is a suitably-chosen permutation for which \(A_1\) is invertible. The method is most useful when \(m \approx n\) as then the dimension of \(R\) is small and a dense Cholesky factorization may be used.

matrix storage#

unsymmetric 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 \(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\) matrix \(H\), as well as the \(m\) by \(m\) matrix \(C\), 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). We focus on \(H\), but everything we say applied equally to \(C\).

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 sbls package must be called in the following order:

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), Float64 (double precision) or, if supported, Float128 (quadruple precision).

callable functions#

    function sbls_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 sbls_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 sbls_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/sbls/SBLS.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/sbls.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 sbls_control_type)

specfile

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

    function sbls_import(T, control, data, status, n, m,
                         H_type, H_ne, H_row, H_col, H_ptr,
                         A_type, A_ne, A_row, A_col, A_ptr,
                         C_type, C_ne, C_row, C_col, C_ptr)

Import structural matrix data into internal storage prior to solution.

Parameters:

control

is a structure whose members provide control parameters for the remaining procedures (see sbls_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 or m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, ‘diagonal’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’ has been violated.

n

is a scalar variable of type Int32 that holds the number of rows in the symmetric matrix \(H\).

m

is a scalar variable of type Int32 that holds the number of rows in the symmetric matrix \(C\).

H_type

is a one-dimensional array of type Vararg{Cchar} that specifies the symmetric storage scheme used for the matrix \(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 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, diagonal or (scaled) identity 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.

A_type

is a one-dimensional array of type Vararg{Cchar} that specifies the unsymmetric storage scheme used for the matrix \(A\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’ or ‘absent’, the latter if access to the Jacobian is via matrix-vector products; lower or upper case variants are allowed.

A_ne

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

C_type

is a one-dimensional array of type Vararg{Cchar} that specifies the symmetric storage scheme used for the matrix \(C\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, ‘diagonal’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’, the latter pair if \(C=0\); lower or upper case variants are allowed.

C_ne

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

C_row

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

C_col

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

C_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 \(C\), 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 sbls_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 sbls_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 sbls_factorize_matrix(T, data, status, n, h_ne, H_val,
                               a_ne, A_val, c_ne, C_val, D)

Form and factorize the block matrix

\[\begin{split}K_{G} = \begin{pmatrix}G & A^T \\ A & - C\end{pmatrix}\end{split}\]
for some appropriate matrix \(G\).

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 factors were generated successfully.

  • -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’, ‘scaled_identity’, ‘identity’, ‘zero’ or ‘none’ has been violated.

  • -9

    An error was reported by SLS analyse. The return status from SLS analyse is given in inform.sls_inform.status. See the documentation for the GALAHAD package SLS for further details.

  • -10

    An error was reported by SLS_factorize. The return status from SLS factorize is given in inform.sls_inform.status. See the documentation for the GALAHAD package SLS for further details.

  • -13

    An error was reported by ULS_factorize. The return status from ULS_factorize is given in inform.uls_factorize_status. See the documentation for the GALAHAD package ULS for further details.

  • -15

    The computed preconditioner \(K_G\) is singular and is thus unsuitable.

  • -20

    The computed preconditioner \(K_G\) has the wrong inertia and is thus unsuitable.

  • -24

    An error was reported by the GALAHAD package SORT_reorder_by_rows. The return status from SORT_reorder_by_rows is given in inform.sort_status. See the documentation for the GALAHAD package SORT for further details.

n

is a scalar variable of type Int32 that holds the number of rows in the symmetric matrix \(H\).

h_ne

is a scalar variable of type Int32 that holds the number of entries in the lower triangular part of the symmetric 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 symmetric matrix \(H\) in any of the available storage schemes

a_ne

is a scalar variable of type Int32 that holds the number of entries in the unsymmetric matrix \(A\).

A_val

is a one-dimensional array of size a_ne and type T that holds the values of the entries of the unsymmetric matrix \(A\) in any of the available storage schemes.

c_ne

is a scalar variable of type Int32 that holds the number of entries in the lower triangular part of the symmetric matrix \(C\).

C_val

is a one-dimensional array of size c_ne and type T that holds the values of the entries of the lower triangular part of the symmetric matrix \(C\) in any of the available storage schemes

D

is a one-dimensional array of size n and type T that holds the values of the entries of the diagonal matrix \(D\) that is required if the user has specified control.preconditioner = 5. It need not be set otherwise.

    function sbls_solve_system(T, data, status, n, m, sol)

Solve the block linear system

\[\begin{split}\begin{pmatrix}G & A^T \\ A & - C\end{pmatrix} \begin{pmatrix}x \\ y\end{pmatrix} = \begin{pmatrix}a \\ b\end{pmatrix}.\end{split}\]

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 required solution was obtained.

  • -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.

  • -11

    An error was reported by SLS_solve. The return status from SLS solve is given in inform.sls_inform.status. See the documentation for the GALAHAD package SLS for further details.

  • -14

    An error was reported by ULS_solve. The return status from ULS_solve is given in inform.uls_solve_status. See the documentation for the GALAHAD package ULS for further details.

n

is a scalar variable of type Int32 that holds the number of entries in the vector \(a\).

m

is a scalar variable of type Int32 that holds the number of entries in the vector \(b\).

sol

is a one-dimensional array of size n + m and type T. on entry, its first n entries must hold the vector \(a\), and the following entries must hold the vector \(b\). On a successful exit, its first n entries contain the solution components \(x\), and the following entries contain the components \(y\).

    function sbls_information(T, data, inform, status)

Provides output information

Parameters:

data

holds private internal data

inform

is a structure containing output information (see sbls_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 sbls_terminate(T, data, control, inform)

Deallocate all internal private storage

Parameters:

data

holds private internal data

control

is a structure containing control information (see sbls_control_type)

inform

is a structure containing output information (see sbls_inform_type)

available structures#

sbls_control_type structure#

    struct sbls_control_type{T}
      f_indexing::Bool
      error::Int32
      out::Int32
      print_level::Int32
      indmin::Int32
      valmin::Int32
      len_ulsmin::Int32
      itref_max::Int32
      maxit_pcg::Int32
      new_a::Int32
      new_h::Int32
      new_c::Int32
      preconditioner::Int32
      semi_bandwidth::Int32
      factorization::Int32
      max_col::Int32
      scaling::Int32
      ordering::Int32
      pivot_tol::T
      pivot_tol_for_basis::T
      zero_pivot::T
      static_tolerance::T
      static_level::T
      min_diagonal::T
      stop_absolute::T
      stop_relative::T
      remove_dependencies::Bool
      find_basis_by_transpose::Bool
      affine::Bool
      allow_singular::Bool
      perturb_to_make_definite::Bool
      get_norm_residual::Bool
      check_basis::Bool
      space_critical::Bool
      deallocate_error_fatal::Bool
      symmetric_linear_solver::NTuple{31,Cchar}
      definite_linear_solver::NTuple{31,Cchar}
      unsymmetric_linear_solver::NTuple{31,Cchar}
      prefix::NTuple{31,Cchar}
      sls_control::sls_control_type{T}
      uls_control::uls_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 print_level

controls level of diagnostic output

Int32 indmin

initial estimate of integer workspace for SLS (obsolete)

Int32 valmin

initial estimate of real workspace for SLS (obsolete)

Int32 len_ulsmin

initial estimate of workspace for ULS (obsolete)

Int32 itref_max

maximum number of iterative refinements with preconditioner allowed

Int32 maxit_pcg

maximum number of projected CG iterations allowed

Int32 new_a

how much has \(A\) changed since last factorization: 0 = not changed, 1 = values changed, 2 = structure changed

Int32 new_h

how much has \(H\) changed since last factorization: 0 = not changed, 1 = values changed, 2 = structure changed

Int32 new_c

how much has \(C\) changed since last factorization: 0 = not changed, 1 = values changed, 2 = structure changed

Int32 preconditioner

which preconditioner to use:

  • 0 selected automatically

  • 1 explicit with \(G = I\)

  • 2 explicit with \(G = H\)

  • 3 explicit with \(G =\) diag(max(\(H\),min_diag))

  • 4 explicit with \(G =\) band \((H)\)

  • 5 explicit with \(G =\) (optional, diagonal) \(D\)

  • 11 explicit with \(G_{11} = 0\), \(G_{21} = 0\), \(G_{22} = H_{22}\)

  • 12 explicit with \(G_{11} = 0\), \(G_{21} = H_{21}\), \(G_{22} = H_{22}\)

  • -1 implicit with \(G_{11} = 0\), \(G_{21} = 0\), \(G_{22} = I\)

  • -2 implicit with \(G_{11} = 0\), \(G_{21} = 0\), \(G_{22} = H_{22}\)

Int32 semi_bandwidth

the semi-bandwidth for band(H)

Int32 factorization

the explicit factorization used:

  • 0 selected automatically

  • 1 Schur-complement if \(G\) is diagonal and successful otherwise augmented system

  • 2 augmented system

  • 3 C_NULL-space

  • 4 Schur-complement if \(G\) is diagonal and successful otherwise failure

  • 5 Schur-complement with pivoting if \(G\) is diagonal and successful otherwise failure

Int32 max_col

maximum number of nonzeros in a column of \(A\) for Schur-complement factorization

Int32 scaling

not used at present

Int32 ordering

see scaling

T pivot_tol

the relative pivot tolerance used by ULS (obsolete)

T pivot_tol_for_basis

the relative pivot tolerance used by ULS when determining the basis matrix

T zero_pivot

the absolute pivot tolerance used by ULS (obsolete)

T static_tolerance

not used at present

T static_level

see static_tolerance

T min_diagonal

the minimum permitted diagonal in diag(max(\(H\),min_diag))

T stop_absolute

the required absolute and relative accuracies

T stop_relative

see stop_absolute

Bool remove_dependencies

preprocess equality constraints to remove linear dependencies

Bool find_basis_by_transpose

determine implicit factorization preconditioners using a basis of A found by examining A’s transpose

Bool affine

can the right-hand side \(c\) be assumed to be zero?

Bool allow_singular

do we tolerate “singular” preconditioners?

Bool perturb_to_make_definite

if the initial attempt at finding a preconditioner is unsuccessful, should the diagonal be perturbed so that a second attempt succeeds?

Bool get_norm_residual

compute the residual when applying the preconditioner?

Bool check_basis

if an implicit or C_NULL-space preconditioner is used, assess and correct for ill conditioned basis matrices

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 symmetric_linear_solver[31]

indefinite linear equation solver

char definite_linear_solver[31]

definite 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 uls_control_type uls_control

control parameters for ULS

sbls_time_type structure#

    struct sbls_time_type{T}
      total::T
      form::T
      factorize::T
      apply::T
      clock_total::T
      clock_form::T
      clock_factorize::T
      clock_apply::T

detailed documentation#

time derived type as a Julia structure

components#

T total

total cpu time spent in the package

T form

cpu time spent forming the preconditioner \(K_G\)

T factorize

cpu time spent factorizing \(K_G\)

T apply

cpu time spent solving linear systems inolving \(K_G\)

T clock_total

total clock time spent in the package

T clock_form

clock time spent forming the preconditioner \(K_G\)

T clock_factorize

clock time spent factorizing \(K_G\)

T clock_apply

clock time spent solving linear systems inolving \(K_G\)

sbls_inform_type structure#

    struct sbls_inform_type{T}
      status::Int32
      alloc_status::Int32
      bad_alloc::NTuple{81,Cchar}
      sort_status::Int32
      factorization_integer::Int64
      factorization_real::Int64
      preconditioner::Int32
      factorization::Int32
      d_plus::Int32
      rank::Int32
      rank_def::Bool
      perturbed::Bool
      iter_pcg::Int32
      norm_residual::T
      alternative::Bool
      time::sbls_time_type{T}
      sls_inform::sls_inform_type{T}
      uls_inform::uls_inform_type{T}

detailed documentation#

inform derived type as a Julia structure

components#

Int32 status

return status. See SBLS_form_and_factorize 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 sort_status

the return status from the sorting routines

Int64 factorization_integer

the total integer workspace required for the factorization

Int64 factorization_real

the total real workspace required for the factorization

Int32 preconditioner

the preconditioner used

Int32 factorization

the factorization used

Int32 d_plus

how many of the diagonals in the factorization are positive

Int32 rank

the computed rank of \(A\)

Bool rank_def

is the matrix A rank defficient?

Bool perturbed

has the used preconditioner been perturbed to guarantee correct inertia?

Int32 iter_pcg

the total number of projected CG iterations required

T norm_residual

the norm of the residual

Bool alternative

has an “alternative” \(y\) : \(K y = 0\) and \(y^T c > 0\) been found when trying to solve \(K y = c\) for generic \(K\)?

struct sbls_time_type time

timings (see above)

struct sls_inform_type sls_inform

inform parameters from the GALAHAD package SLS used

struct uls_inform_type uls_inform

inform parameters from the GALAHAD package ULS used

example calls#

This is an example of how to use the package to solve a block system of linear equations; the code is available in $GALAHAD/src/sbls/Julia/test_sbls.jl . A variety of supported matrix storage formats are shown.

# test_sbls.jl
# Simple code to test the Julia interface to SBLS

using GALAHAD
using Test
using Printf
using Accessors
using Quadmath

function test_sbls(::Type{T}) where T
  # Derived types
  data = Ref{Ptr{Cvoid}}()
  control = Ref{sbls_control_type{T}}()
  inform = Ref{sbls_inform_type{T}}()

  # Set problem data
  n = 3 # dimension of H
  m = 2 # dimension of C
  H_ne = 4 # number of elements of H
  A_ne = 3 # number of elements of A
  C_ne = 3 # number of elements of C
  H_dense_ne = 6 # number of elements of H
  A_dense_ne = 6 # number of elements of A
  C_dense_ne = 3 # number of elements of C
  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]
  A_row = Cint[1, 1, 2]
  A_col = Cint[1, 2, 3]
  A_ptr = Cint[1, 3, 4]
  C_row = Cint[1, 2, 2]  # row indices, NB lower triangle
  C_col = Cint[1, 1, 2]
  C_ptr = Cint[1, 2, 4]
  H_val = T[1.0, 2.0, 3.0, 1.0]
  A_val = T[2.0, 1.0, 1.0]
  C_val = T[4.0, 1.0, 2.0]
  H_dense = T[1.0, 0.0, 2.0, 1.0, 0.0, 3.0]
  A_dense = T[2.0, 1.0, 0.0, 0.0, 0.0, 1.0]
  C_dense = T[4.0, 1.0, 2.0]
  H_diag = T[1.0, 1.0, 2.0]
  C_diag = T[4.0, 2.0]
  H_scid = T[2.0]
  C_scid = T[2.0]

  st = ' '
  status = Ref{Cint}()

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

  for d in 1:7

    # Initialize SBLS
    sbls_initialize(T, data, control, status)
    @reset control[].preconditioner = Cint(2)
    @reset control[].factorization = Cint(2)
    @reset control[].get_norm_residual = true

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

    # sparse co-ordinate storage
    if d == 1
      st = 'C'
      sbls_import(T, control, data, status, n, m,
                  "coordinate", H_ne, H_row, H_col, C_NULL,
                  "coordinate", A_ne, A_row, A_col, C_NULL,
                  "coordinate", C_ne, C_row, C_col, C_NULL)

      sbls_factorize_matrix(T, data, status, n,
                            H_ne, H_val,
                            A_ne, A_val,
                            C_ne, C_val, C_NULL)
    end

    # sparse by rows
    if d == 2
      st = 'R'
      sbls_import(T, control, data, status, n, m,
                  "sparse_by_rows", H_ne, C_NULL, H_col, H_ptr,
                  "sparse_by_rows", A_ne, C_NULL, A_col, A_ptr,
                  "sparse_by_rows", C_ne, C_NULL, C_col, C_ptr)

      sbls_factorize_matrix(T, data, status, n,
                            H_ne, H_val,
                            A_ne, A_val,
                            C_ne, C_val, C_NULL)
    end

    # dense
    if d == 3
      st = 'D'
      sbls_import(T, control, data, status, n, m,
                  "dense", H_ne, C_NULL, C_NULL, C_NULL,
                  "dense", A_ne, C_NULL, C_NULL, C_NULL,
                  "dense", C_ne, C_NULL, C_NULL, C_NULL)

      sbls_factorize_matrix(T, data, status, n,
                            H_dense_ne, H_dense,
                            A_dense_ne, A_dense,
                            C_dense_ne, C_dense,
                            C_NULL)
    end

    # diagonal
    if d == 4
      st = 'L'
      sbls_import(T, control, data, status, n, m,
                  "diagonal", H_ne, C_NULL, C_NULL, C_NULL,
                  "dense", A_ne, C_NULL, C_NULL, C_NULL,
                  "diagonal", C_ne, C_NULL, C_NULL, C_NULL)

      sbls_factorize_matrix(T, data, status, n,
                            n, H_diag,
                            A_dense_ne, A_dense,
                            m, C_diag,
                            C_NULL)
    end

    # scaled identity
    if d == 5
      st = 'S'
      sbls_import(T, control, data, status, n, m,
                  "scaled_identity", H_ne, C_NULL, C_NULL, C_NULL,
                  "dense", A_ne, C_NULL, C_NULL, C_NULL,
                  "scaled_identity", C_ne, C_NULL, C_NULL, C_NULL)

      sbls_factorize_matrix(T, data, status, n,
                            1, H_scid,
                            A_dense_ne, A_dense,
                            1, C_scid,
                            C_NULL)
    end

    # identity
    if d == 6
      st = 'I'
      sbls_import(T, control, data, status, n, m,
                  "identity", H_ne, C_NULL, C_NULL, C_NULL,
                  "dense", A_ne, C_NULL, C_NULL, C_NULL,
                  "identity", C_ne, C_NULL, C_NULL, C_NULL)

      sbls_factorize_matrix(T, data, status, n,
                            0, H_val,
                            A_dense_ne, A_dense,
                            0, C_val, C_NULL)
    end

    # zero
    if d == 7
      st = 'Z'
      sbls_import(T, control, data, status, n, m,
                  "identity", H_ne, C_NULL, C_NULL, C_NULL,
                  "dense", A_ne, C_NULL, C_NULL, C_NULL,
                  "zero", C_ne, C_NULL, C_NULL, C_NULL)

      sbls_factorize_matrix(T, data, status, n,
                            0, H_val,
                            A_dense_ne, A_dense,
                            0, C_NULL, C_NULL)
    end

    # Set right-hand side (a, b)
    sol = T[3.0, 2.0, 4.0, 2.0, 0.0]  # values

    sbls_solve_system(T, data, status, n, m, sol)

    sbls_information(T, data, inform, status)

    if inform[].status == 0
      @printf("%c: residual = %9.1e status = %1i\n",
              st, inform[].norm_residual, inform[].status)
    else
      @printf("%c: SBLS_solve exit status = %1i\n", st, inform[].status)
    end

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

    # Delete internal workspace
    sbls_terminate(T, data, control, inform)
  end

  return 0
end

@testset "SBLS" begin
  @test test_sbls(Float32) == 0
  @test test_sbls(Float64) == 0
  @test test_sbls(Float128) == 0
end