overview of functions provided#

// namespaces

namespace conf;

// typedefs

typedef float spc_;
typedef double rpc_;
typedef int ipc_;

// structs

struct sbls_control_type;
struct sbls_inform_type;
struct sbls_time_type;

// global functions

void sbls_initialize(
    void **data,
    struct sbls_control_type* control,
    ipc_ *status
);

void sbls_read_specfile(
    struct sbls_control_type* control,
    const char specfile[]
);

void sbls_import(
    struct sbls_control_type* control,
    void **data,
    ipc_ *status,
    ipc_ n,
    ipc_ m,
    const char H_type[],
    ipc_ H_ne,
    const ipc_ H_row[],
    const ipc_ H_col[],
    const ipc_ H_ptr[],
    const char A_type[],
    ipc_ A_ne,
    const ipc_ A_row[],
    const ipc_ A_col[],
    const ipc_ A_ptr[],
    const char C_type[],
    ipc_ C_ne,
    const ipc_ C_row[],
    const ipc_ C_col[],
    const ipc_ C_ptr[]
);

void sbls_reset_control(
    struct sbls_control_type* control,
    void **data,
    ipc_ *status
);

void sbls_factorize_matrix(
    void **data,
    ipc_ *status,
    ipc_ n,
    ipc_ h_ne,
    const rpc_ H_val[],
    ipc_ a_ne,
    const rpc_ A_val[],
    ipc_ c_ne,
    const rpc_ C_val[],
    const rpc_ D[]
);

void sbls_solve_system(void **data, ipc_ *status, ipc_ n, ipc_ m, rpc_ sol[]);
void sbls_information(void **data, struct sbls_inform_type* inform, ipc_ *status);

void sbls_terminate(
    void **data,
    struct sbls_control_type* control,
    struct sbls_inform_type* inform
);

typedefs#

typedef float spc_

spc_ is real single precision

typedef double rpc_

rpc_ is the real working precision used, but may be changed to float by defining the preprocessor variable REAL_32 or (if supported) to __real128 using the variable REAL_128.

typedef int ipc_

ipc_ is the default integer word length used, but may be changed to int64_t by defining the preprocessor variable INTEGER_64.

function calls#

void sbls_initialize(
    void **data,
    struct sbls_control_type* control,
    ipc_ *status
)

Set default control values and initialize private data

Parameters:

data

holds private internal data

control

is a struct containing control information (see sbls_control_type)

status

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

  • 0

    The initialization was successful.

void sbls_read_specfile(
    struct sbls_control_type* control,
    const char 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 struct containing control information (see sbls_control_type)

specfile

is a character string containing the name of the specification file

void sbls_import(
    struct sbls_control_type* control,
    void **data,
    ipc_ *status,
    ipc_ n,
    ipc_ m,
    const char H_type[],
    ipc_ H_ne,
    const ipc_ H_row[],
    const ipc_ H_col[],
    const ipc_ H_ptr[],
    const char A_type[],
    ipc_ A_ne,
    const ipc_ A_row[],
    const ipc_ A_col[],
    const ipc_ A_ptr[],
    const char C_type[],
    ipc_ C_ne,
    const ipc_ C_row[],
    const ipc_ C_col[],
    const ipc_ C_ptr[]
)

Import structural matrix data into internal storage prior to solution.

Parameters:

control

is a struct whose members provide control paramters for the remaining prcedures (see sbls_control_type)

data

holds private internal data

status

is a scalar variable of type ipc_, 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 ipc_, that holds the number of rows in the symmetric matrix \(H\).

m

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

H_type

is a one-dimensional array of type char 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 ipc_, 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 ipc_, 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 NULL.

H_col

is a one-dimensional array of size H_ne and type ipc_, 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 NULL.

H_ptr

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

A_type

is a one-dimensional array of type char that specifies the symmetric 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 ipc_, 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 ipc_, 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 NULL.

A_col

is a one-dimensional array of size A_ne and type ipc_, 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 NULL.

A_ptr

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

C_type

is a one-dimensional array of type char 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 ipc_, 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 ipc_, 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 NULL.

C_col

is a one-dimensional array of size C_ne and type ipc_, 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 NULL.

C_ptr

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

void sbls_reset_control(
    struct sbls_control_type* control,
    void **data,
    ipc_ *status
)

Reset control parameters after import if required.

Parameters:

control

is a struct whose members provide control paramters for the remaining prcedures (see sbls_control_type)

data

holds private internal data

status

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

  • 0

    The import was successful.

void sbls_factorize_matrix(
    void **data,
    ipc_ *status,
    ipc_ n,
    ipc_ h_ne,
    const rpc_ H_val[],
    ipc_ a_ne,
    const rpc_ A_val[],
    ipc_ c_ne,
    const rpc_ C_val[],
    const rpc_ 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 ipc_, 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 ipc_, that holds the number of rows in the symmetric matrix \(H\).

h_ne

is a scalar variable of type ipc_, 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 rpc_, 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 ipc_, that holds the number of entries in the unsymmetric matrix \(A\).

A_val

is a one-dimensional array of size a_ne and type rpc_, 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 ipc_, 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 rpc_, 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 rpc_, 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.

void sbls_solve_system(void **data, ipc_ *status, ipc_ n, ipc_ m, rpc_ 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 ipc_, 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 ipc_, that holds the number of entries in the vector \(a\).

m

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

sol

is a one-dimensional array of size n + m and type double. 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\).

void sbls_information(void **data, struct sbls_inform_type* inform, ipc_ *status)

Provides output information

Parameters:

data

holds private internal data

inform

is a struct containing output information (see sbls_inform_type)

status

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

  • 0

    The values were recorded successfully

void sbls_terminate(
    void **data,
    struct sbls_control_type* control,
    struct sbls_inform_type* inform
)

Deallocate all internal private storage

Parameters:

data

holds private internal data

control

is a struct containing control information (see sbls_control_type)

inform

is a struct containing output information (see sbls_inform_type)