SBLS#

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#

The unsymmetric \(m\) by \(n\) matrix \(A\) may be presented and stored in a variety of convenient input formats.

Dense storage format: The matrix \(A\) is stored as a compact dense matrix by rows, that is, the values of the entries of each row in turn are stored in order within an appropriate real one-dimensional array. In this case, component \(n \ast i + j\) of the storage array A_val will hold the value \(A_{ij}\) for \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\). The string A_type = ‘dense’ should be specified.

Dense by columns storage format: The matrix \(A\) is stored as a compact dense matrix by columns, that is, the values of the entries of each column in turn are stored in order within an appropriate real one-dimensional array. In this case, component \(m \ast j + i\) of the storage array A_val will hold the value \(A_{ij}\) for \(0 \leq i \leq m-1\), $0 \leq j \leq n-1 The string A_type = ‘dense_by_columns’ should be specified.

Sparse co-ordinate storage format: Only the nonzero entries of the matrices are stored. For the \(l\)-th entry, \(0 \leq l \leq ne-1\), of \(A\), its row index i, column index j and value \(A_{ij}\), \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\), are stored as the \(l\)-th components of the integer arrays A_row and A_col and real array A_val, respectively, while the number of nonzeros is recorded as A_ne = \(ne\). The string A_type = ‘coordinate’should be specified.

Sparse row-wise storage format: Again only the nonzero entries are stored, but this time they are ordered so that those in row i appear directly before those in row i+1. For the i-th row of \(A\) the i-th component of the integer array A_ptr holds the position of the first entry in this row, while A_ptr(m) holds the total number of entries. The column indices j, \(0 \leq j \leq n-1\), and values \(A_{ij}\) of the nonzero entries in the i-th row are stored in components l = A_ptr(i), \(\ldots\), A_ptr(i+1)-1, \(0 \leq i \leq m-1\), of the integer array A_col, and real array A_val, respectively. For sparse matrices, this scheme almost always requires less storage than its predecessor. The string A_type = ‘sparse_by_rows’ should be specified.

Sparse column-wise storage format: Once again only the nonzero entries are stored, but this time they are ordered so that those in column j appear directly before those in column j+1. For the j-th column of \(A\) the j-th component of the integer array A_ptr holds the position of the first entry in this column, while A_ptr(n) holds the total number of entries. The row indices i, \(0 \leq i \leq m-1\), and values \(A_{ij}\) of the nonzero entries in the j-th columnsare stored in components l = A_ptr(j), \(\ldots\), A_ptr(j+1)-1, \(0 \leq j \leq n-1\), of the integer array A_row, and real array A_val, respectively. As before, for sparse matrices, this scheme almost always requires less storage than the co-ordinate format. The string A_type = ‘sparse_by_columns’ should be specified.

The symmetric \(n\) by \(n\) matrix \(H\), 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 \(0 \leq j \leq i \leq n-1\)) need be held. In this case the lower triangle should be stored by rows, that is component \(i * i / 2 + j\) of the storage array H_val will hold the value \(H_{ij}\) (and, by symmetry, \(H_{ji}\)) for \(0 \leq j \leq i \leq n-1\). The string H_type = ‘dense’ should be specified.

Sparse co-ordinate storage format: Only the nonzero entries of the matrices are stored. For the \(l\)-th entry, \(0 \leq l \leq ne-1\), of \(H\), its row index i, column index j and value \(H_{ij}\), \(0 \leq j \leq i \leq n-1\), are stored as the \(l\)-th components of the integer arrays H_row and H_col and real array H_val, respectively, while the number of nonzeros is recorded as H_ne = \(ne\). Note that only the entries in the lower triangle should be stored. The string H_type = ‘coordinate’ should be specified.

Sparse row-wise storage format: Again only the nonzero entries are stored, but this time they are ordered so that those in row i appear directly before those in row i+1. For the i-th row of \(H\) the i-th component of the integer array H_ptr holds the position of the first entry in this row, while H_ptr(n) holds the total number of entries. The column indices j, \(0 \leq j \leq i\), and values \(H_{ij}\) of the entries in the i-th row are stored in components l = H_ptr(i), …, H_ptr(i+1)-1 of the integer array H_col, and real array H_val, respectively. Note that as before only the entries in the lower triangle should be stored. For sparse matrices, this scheme almost always requires less storage than its predecessor. The string H_type = ‘sparse_by_rows’ should be specified.

Diagonal storage format: If \(H\) is diagonal (i.e., \(H_{ij} = 0\) for all \(0 \leq i \neq j \leq n-1\)) only the diagonals entries \(H_{ii}\), \(0 \leq i \leq n-1\) need be stored, and the first n components of the array H_val may be used for the purpose. The string H_type = ‘diagonal’ should be specified.

Multiples of the identity storage format: If \(H\) is a multiple of the identity matrix, (i.e., \(H = \alpha I\) where \(I\) is the n by n identity matrix and \(\alpha\) is a scalar), it suffices to store \(\alpha\) as the first component of H_val. The string H_type = ‘scaled_identity’ should be specified.

The identity matrix format: If \(H\) is the identity matrix, no values need be stored. The string H_type = ‘identity’ should be specified.

The zero matrix format: The same is true if \(H\) is the zero matrix, but now the string H_type = ‘zero’ or ‘none’ should be specified.

functions#

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

print_levelint

the level of output required is specified by print_level. Possible values are

  • <=0

    gives no output,

  • 1

    gives a summary of the progress of the method.

  • >=2

    gives increasingly verbose (debugging) output.

indminint

initial estimate of integer workspace for SLS (obsolete).

valminint

initial estimate of real workspace for SLS (obsolete).

len_ulsminint

initial estimate of workspace for ULS (obsolete).

itref_maxint

maximum number of iterative refinements with preconditioner allowed.

maxit_pcgint

maximum number of projected CG iterations allowed.

new_aint

how much has \(A\) changed since last factorization. Possible values are

  • 0

    unchanged.

  • 1

    values but not indices have changed.

  • 2

    values and indices have changed.

new_hint

how much has \(H\) changed since last factorization. Possible values are

  • 0

    unchanged.

  • 1

    values but not indices have changed.

  • 2

    values and indices have changed.

new_cint

how much has \(C\) changed since last factorization. Possible values are

  • 0

    unchanged.

  • 1

    values but not indices have changed.

  • 2

    values and indices have changed.

preconditionerint

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}\).

semi_bandwidthint

the semi-bandwidth for band(H).

factorizationint

the explicit factorization used:

  • 0

    selected automatically

  • 1

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

  • 2

    augmented system

  • 3

    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.

max_colint

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

scalingint

not used at present.

orderingint

see scaling.

pivot_tolfloat

the relative pivot tolerance used by ULS (obsolete).

pivot_tol_for_basisfloat

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

zero_pivotfloat

the absolute pivot tolerance used by ULS (obsolete).

static_tolerancefloat

not used at present.

static_levelfloat

see static_tolerance.

min_diagonalfloat

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

stop_absolutefloat

the required absolute and relative accuracies.

stop_relativefloat

see stop_absolute.

remove_dependenciesbool

preprocess equality constraints to remove linear dependencies.

find_basis_by_transposebool

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

affinebool

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

allow_singularbool

do we tolerate “singular” preconditioners?.

perturb_to_make_definitebool

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

get_norm_residualbool

compute the residual when applying the preconditioner?.

check_basisbool

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

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.

symmetric_linear_solverstr

indefinite linear equation solver used.

definite_linear_solverstr

definite linear equation solver used.

unsymmetric_linear_solverstr

unsymmetric linear equation solver used.

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

uls_optionsdict

default control options for ULS (see uls.initialize).

sbls.load(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, options=None)#

Import problem data into internal storage prior to solution.

Parameters:

nint

holds the dimension of \(H\) (equivalently, the number of columns of \(A\)).

mint

holds the dimension of \(C\) (equivalently, the number of rows of \(A\)).

H_typestring

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’; lower or upper case variants are allowed.

H_neint

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_rowndarray(H_ne)

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 schemes, and in this case can be None.

H_colndarray(H_ne)

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 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 \(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 None.

A_typestring

specifies the unsymmetric storage scheme used for the matrix \(A\), 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 \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other two schemes.

A_rowndarray(A_ne)

holds the row indices of \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other two schemes, and in this case can be None.

A_colndarray(A_ne)

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 storage scheme is used, and in this case can be None.

A_ptrndarray(m+1)

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

C_typestring

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’; lower or upper case variants are allowed.

C_neint

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_rowndarray(C_ne)

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 schemes, and in this case can be None.

C_colndarray(C_ne)

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 other storage schemes are used, and in this case can be None.

C_ptrndarray(m+1)

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

optionsdict, optional

dictionary of control options (see sbls.initialize).

sbls.factorize_matrix(n, m, 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:

nint

holds the dimension of \(H\) (equivalently, the number of columns of \(A\)).

mint

holds the dimension of \(C\) (equivalently, the number of rows of \(A\)).

H_neint

holds the number of entries in the lower triangular part of the matrix \(H\).

H_valndarray(H_ne)

holds the values of the nonzeros in the lower triangle of the matrix \(H\) in the same order as specified in the sparsity pattern in sbls.load.

A_neint

holds the number of entries in the matrix \(A\).

A_valndarray(A_ne)

holds the values of the nonzeros in the matrix \(A\) in the same order as specified in the sparsity pattern in sbls.load.

C_neint

holds the number of entries in the lower triangular part of the matrix \(C\).

C_valndarray(C_ne)

holds the values of the nonzeros in the lower triangle of the matrix \(C\) in the same order as specified in the sparsity pattern in sbls.load.

Dndarray(n)

holds the values of the diagonals of the matrix \(D\) that is required if options[preconditioner]=5 has been specified. Otherwise it shuld be set to None.

sbls.solve_system(n, m, rhs)#

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:

nint

holds the dimension of \(H\) (equivalently, the number of columns of \(A\)).

mint

holds the dimension of \(C\) (equivalently, the number of rows of \(A\)).

solndarray(n+m)

holds the values of the right-hand side vector \((a,b)\).

Returns:

solndarray(n+m)

holds the values of the solution vector \((x,y)\) after a successful call.

[optional] sbls.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 or m > 0 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’].

  • -12

    The analysis phase of an unsymmetric factorization failed; the return status from the factorization package is given by inform[‘factor_status’].

  • -13

    An unsymmetric factorization failed; the return status from the factorization package is given by inform[‘factor_status’].

  • -15

    The computed preconditioner \(P_G\) is singular, and is thus unsuitable

  • -20

    The computed preconditioner \(P_G\) has the wrong inertia, and is thus unsuitable

  • -24

    An error was reported by the sort routine; the return status is returned in sort_status.

alloc_statusint

the status of the last attempted allocation/deallocation.

bad_allocstr

the name of the array for which an allocation/deallocation error occurred.

sort_statusint

the return status from the sorting routines.

factorization_integerlong

the total integer workspace required for the factorization.

factorization_reallong

the total real workspace required for the factorization.

preconditionerint

the preconditioner used.

factorizationint

the factorization used.

d_plusint

how many of the diagonals in the factorization are positive.

rankint

the computed rank of \(A\).

rank_defbool

is the matrix A rank defficient?.

perturbedbool

has the used preconditioner been perturbed to guarantee correct inertia?.

iter_pcgint

the total number of projected CG iterations required.

norm_residualfloat

the norm of the residual.

alternativebool

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

timedict
dictionary containing timing information:
totalfloat

total cpu time spent in the package.

formfloat

cpu time spent forming the preconditioner \(K_G\).

factorizefloat

cpu time spent factorizing \(K_G\).

applyfloat

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

clock_totalfloat

total clock time spent in the package.

clock_formfloat

clock time spent forming the preconditioner \(K_G\).

clock_factorizefloat

clock time spent factorizing \(K_G\).

clock_applyfloat

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

sls_informdict

inform parameters for SLS (see sls.information).

uls_informdict

inform parameters for ULS (see uls.information).

sbls.terminate()#

Deallocate all internal private storage.

example code#

from galahad import sbls
import numpy as np
np.set_printoptions(precision=4,suppress=True,floatmode='fixed')
print("\n** python test: sbls")

#  describe problem (only the lower triangle of matrix H is required)
#  ( 1     |  2   )     (   3   )
#  (   2 1 |  1 1 ) x = (   5   )
#  (   1 3 |    1 )     (   5   )
#  ( ------------ )     ( ----- )
#  ( 2 1   | -t   )     ( 3 - t )
#  (   1 1 |   -t )     ( 2 - t )

n = 3
m = 2

H_type = 'coordinate'
H_ne = 4
H_row = np.array([0,1,2,2])
H_col = np.array([0,1,1,2])
H_ptr = None
H_val = np.array([1.0,2.0,1.0,3.0])

A_type = 'coordinate'
A_ne = 4
A_row = np.array([0,0,1,1])
A_col = np.array([0,1,1,2])
A_ptr = None
A_val = np.array([2.0,1.0,1.0,1.0])

C_type = 'scaled_identity'
C_ne = 1
C_row = None
C_col = None
C_ptr = None
#t = 0.00000001
t = 0.0
C_val = np.array([t])

# allocate internal data and set default options
options = sbls.initialize()

# set some non-default options
#options['print_level'] = 1
options['symmetric_linear_solver'] = 'sytr '
options['definite_linear_solver'] = 'sytr '
options['preconditioner'] = 2
options['factorization'] = 2
#print(options['preconditioner'])
#print("options:", options)

# load data (and optionally non-default options), and analyse matrix structure
sbls.load(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, options)

# factorize matrix
sbls.factorize_matrix(n, m, H_ne, H_val, A_ne, A_val, C_ne, C_val)

# solve system
rhs = np.array([3.0,5.0,5.0,3.0-t,2.0-t])
x = sbls.solve_system(n, m, rhs)
print(" x:",x)

# get information
inform = sbls.information()
print(" rank:",inform['rank'])
print('** sbls exit status:', inform['status'])

# deallocate internal data

sbls.terminate()

This example code is available in $GALAHAD/src/sbls/Python/test_sbls.py .