SSLS#

purpose#

Given a (possibly rectangular) matrix \(A\) and symmetric matrices \(H\) and \(C\), form and factorize the block, symmetric matrix

\[\begin{split}K = \begin{pmatrix}H & A^T \\ A & - C\end{pmatrix},\end{split}\]
and susequently solve systems
\[\begin{split}\begin{pmatrix}H & A^T \\ A & - C\end{pmatrix} \begin{pmatrix}x \\ y\end{pmatrix} = \begin{pmatrix}a \\ b\end{pmatrix},\end{split}\]
using the GALAHAD symmetric-indefinite factorization package SLS Full advantage is taken of any zero coefficients in the matrices \(H\), \(A\) and \(C\).

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

method#

The method simply assembles \(K\) from its components, and then relies on SLS for analysis, factorization and solves.

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#

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

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.

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

ssls.analyse(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)#

Assmeble the structure of the block matrix

\[\begin{split}K = \begin{pmatrix}H & A^T \\ A & - C\end{pmatrix}.\end{split}\]
and analyse the resulting matrix prior to factorization.

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

ssls.factorize_matrix(n, m, H_ne, H_val, A_ne, A_val, C_ne, C_val)#

Factorize the block matrix

\[\begin{split}K = \begin{pmatrix}H & A^T \\ A & - C\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\)).

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 ssls.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 ssls.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 ssls.load.

ssls.solve_system(n, m, rhs)#

Solve the block linear system

\[\begin{split}\begin{pmatrix}H & 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] ssls.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’].

alloc_statusint

the status of the last attempted allocation/deallocation.

bad_allocstr

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

factorization_integerlong

the total integer workspace required for the factorization.

factorization_reallong

the total real workspace required for the factorization.

rankint

the computed rank of \(A\).

rank_defbool

is the matrix A rank defficient?.

timedict
dictionary containing timing information:
totalfloat

total cpu time spent in the package.

analysefloat

cpu time spent forming and analysing \(K\).

factorizefloat

cpu time spent factorizing \(K\).

solvefloat

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

clock_totalfloat

total clock time spent in the package.

clock_analysefloat

clock time spent forming and analysing \(K\).

clock_factorizefloat

clock time spent factorizing \(K\).

clock_solvefloat

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

sls_informdict

inform parameters for SLS (see sls.information).

ssls.terminate()#

Deallocate all internal private storage.

example code#

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

#  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 = ssls.initialize()

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

# load data (and optionally non-default options), and analyse matrix structure
ssls.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
ssls.factorize_matrix(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 = ssls.solve_system(n, m, rhs)
print(" x:",x)

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

# deallocate internal data

ssls.terminate()

This example code is available in $GALAHAD/src/ssls/Python/test_ssls.py .