GALAHAD TRS package#

purpose#

The trs package uses matrix factorization to find the global minimizer of a quadratic objective function within an ellipsoidal region; this is commonly known as the trust-region subproblem. The aim is to minimize the quadratic objective function

\[q(x) = f + g^T x + \frac{1}{2} x^T H x,\]
where the vector \(x\) is required to satisfy the ellipsoidal trust-region constraint \(\|x\|_M \leq \Delta\), and optionally affine constraints \(A x = 0\), where the \(M\)-norm of \(x\) is defined to be \(\|x\|_M = \sqrt{x^T M x}\), and where the radius \(\Delta > 0\).

The package may also be used to solve the related problem in which \(x\) is instead required to satisfy the equality constraint \(\|x\|_M = \Delta\). The matrix \(M\) need not be provided in the commonly-occurring \(\ell_2\)-trust-region case for which \(M = I\), the \(n\) by \(n\) identity matrix.

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

\[\begin{split}\left(\begin{array}{cc} H + \lambda M & A^T \\ A & 0 \end{array}\right)\;\;\mbox{(1)}\end{split}\]
in cases where \(A x = 0\) is imposed, for a succession of scalars \(\lambda\) will be required, so this package is most suited for the case where such a factorization may be found efficiently. If this is not the case, the package gltr may be preferred.

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

method#

The required solution \(x_*\) necessarily satisfies the optimality condition \(H x_* + \lambda_* M x_* + A^T y_* + g = 0\) and \(A x_* = 0\), where \(\lambda_* \geq 0\) is a Lagrange multiplier corresponding to the constraint \(\|x\|_M \leq \Delta\) and \(y_*\) are Lagrange multipliers for the linear constraints \(A x = 0\), if any; for the equality-constrained problem \(\|x\|_M = \Delta\), the multiplier is unconstrained. In addition in all cases, the matrix \(H + \lambda_* M\) will be positive semi-definite on the null-space of \(A\); in most instances it will actually be positive definite, but in special ``hard’’ cases singularity is a possibility.

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

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

The dominant cost is the requirement that we solve a sequence of linear systems (2). In the absence of linear constraints, an efficient sparse Cholesky factorization with precautions to detect indefinite \(H + \lambda M\) is used. If \(A x = 0\) is required, a sparse symmetric, indefinite factorization of (1) is used rather than a Cholesky factorization.

reference#

The method is described in detail in

H. S. Dollar, N. I. M. Gould and D. P. Robinson. ``On solving trust-region and other regularised subproblems in optimization’’. Mathematical Programming Computation 2(1) (2010) 21–57.

matrix storage#

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

Dense storage format: The matrix \(A\) is stored as a compact dense matrix by rows, that is, the values of the entries of each row in turn are stored in order within an appropriate real one-dimensional array. In this case, component \(n \ast i + j\) of the storage array A_val will hold the value \(A_{ij}\) for \(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\) matrices \(H\) and, optionally. \(M\) may also be presented and stored in a variety of formats. But crucially symmetry is exploited by only storing values from the lower triangular part (i.e, those entries that lie on or below the leading diagonal).

Dense storage format: The matrix \(H\) is stored as a compact dense matrix by rows, that is, the values of the entries of each row in turn are stored in order within an appropriate real one-dimensional array. Since \(H\) is symmetric, only the lower triangular part (that is the part \(H_{ij}\) for \(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.

introduction to function calls#

To solve a given problem, functions from the trs package must be called in the following order:

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

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

  • trs_import - set up problem data structures and fixed values

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

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

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

  • trs_solve_problem - solve the trust-region problem

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

  • trs_terminate - deallocate data structures

See the examples section for illustrations of use.

callable functions#

overview of functions provided#

// typedefs

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

// structs

struct trs_control_type;
struct trs_time_type;
struct trs_history_type;
struct trs_inform_type;

// global functions

void trs_initialize(void **data, struct trs_control_type* control, ipc_ *status);
void trs_read_specfile(struct trs_control_type* control, const char specfile[]);

void trs_import(
    struct trs_control_type* control,
    void **data,
    ipc_ *status,
    ipc_ n,
    const char H_type[],
    ipc_ H_ne,
    const ipc_ H_row[],
    const ipc_ H_col[],
    const ipc_ H_ptr[]
);

void trs_import_m(
    void **data,
    ipc_ *status,
    ipc_ n,
    const char M_type[],
    ipc_ M_ne,
    const ipc_ M_row[],
    const ipc_ M_col[],
    const ipc_ M_ptr[]
);

void trs_import_a(
    void **data,
    ipc_ *status,
    ipc_ m,
    const char A_type[],
    ipc_ A_ne,
    const ipc_ A_row[],
    const ipc_ A_col[],
    const ipc_ A_ptr[]
);

void trs_reset_control(
    struct trs_control_type* control,
    void **data,
    ipc_ *status
);

void trs_solve_problem(
    void **data,
    ipc_ *status,
    ipc_ n,
    const rpc_ radius,
    const rpc_ f,
    const rpc_ c[],
    ipc_ H_ne,
    const rpc_ H_val[],
    rpc_ x[],
    ipc_ M_ne,
    const rpc_ M_val[],
    ipc_ m,
    ipc_ A_ne,
    const rpc_ A_val[],
    rpc_ y[]
);

void trs_information(void **data, struct trs_inform_type* inform, ipc_ *status);

void trs_terminate(
    void **data,
    struct trs_control_type* control,
    struct trs_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 SINGLE.

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 trs_initialize(void **data, struct trs_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 trs_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 trs_read_specfile(struct trs_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/trs/TRS.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/trs.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 trs_control_type)

specfile

is a character string containing the name of the specification file

void trs_import(
    struct trs_control_type* control,
    void **data,
    ipc_ *status,
    ipc_ n,
    const char H_type[],
    ipc_ H_ne,
    const ipc_ H_row[],
    const ipc_ H_col[],
    const ipc_ H_ptr[]
)

Import problem data into internal storage prior to solution.

Parameters:

control

is a struct whose members provide control paramters for the remaining prcedures (see trs_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 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, diagonal’ or ‘identity’ has been violated.

n

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

H_type

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

H_ne

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

void trs_import_m(
    void **data,
    ipc_ *status,
    ipc_ n,
    const char M_type[],
    ipc_ M_ne,
    const ipc_ M_row[],
    const ipc_ M_col[],
    const ipc_ M_ptr[]
)

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

Parameters:

data

holds private internal data

status

is a scalar variable of type 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 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, diagonal’ or ‘identity’ has been violated.

n

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

M_type

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

M_ne

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

M_row

is a one-dimensional array of size M_ne and type ipc_, that holds the row indices of the lower triangular part of \(M\) in the sparse co-ordinate storage scheme. It need not be set for any of the other three schemes, and in this case can be NULL.

M_col

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

M_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 \(M\), as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the other schemes are used, and in this case can be NULL.

void trs_import_a(
    void **data,
    ipc_ *status,
    ipc_ m,
    const char A_type[],
    ipc_ A_ne,
    const ipc_ A_row[],
    const ipc_ A_col[],
    const ipc_ A_ptr[]
)

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

Parameters:

data

holds private internal data

status

is a scalar variable of type 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 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’ or ‘sparse_by_rows’ has been violated.

m

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

A_type

is a one-dimensional array of type char that specifies the unsymmetric storage scheme used for the constraint Jacobian, \(A\) if any. It should be one of ‘coordinate’, ‘sparse_by_rows’ or ‘dense’; lower or upper case variants are allowed.

A_ne

is a scalar variable of type ipc_, that holds the number of entries in \(A\), if used, in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes.

A_row

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

void trs_reset_control(
    struct trs_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 trs_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:

    1. The import was successful.

void trs_solve_problem(
    void **data,
    ipc_ *status,
    ipc_ n,
    const rpc_ radius,
    const rpc_ f,
    const rpc_ c[],
    ipc_ H_ne,
    const rpc_ H_val[],
    rpc_ x[],
    ipc_ M_ne,
    const rpc_ M_val[],
    ipc_ m,
    ipc_ A_ne,
    const rpc_ A_val[],
    rpc_ y[]
)

Solve the trust-region problem.

Parameters:

data

holds private internal data

status

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

On initial entry, status must be set to 1.

Possible exit values are:

  • 0

    The run was successful.

  • -1

    An allocation error occurred. A message indicating the offending array is written on unit control.error, and the returned allocation status and a string containing the name of the offending array are held in inform.alloc_status and inform.bad_alloc respectively.

  • -2

    A deallocation error occurred. A message indicating the offending array is written on unit control.error and the returned allocation status and a string containing the name of the offending array are held in inform.alloc_status and inform.bad_alloc respectively.

  • -3

    The restrictions n > 0, radius > 0 and m > 0 or requirement that a type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, ‘diagonal’ or ‘identity’ has been violated.

  • -9

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

  • -10

    The factorization of the matrix (1) failed.

  • -15

    The matrix M appears not to be diagonally dominant.

  • -16

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

  • -18

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

n

is a scalar variable of type ipc_, that holds the number of variables

radius

is a scalar of type rpc_, that holds the trust-region radius, \(\Delta\), used. radius must be strictly positive

f

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

c

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

H_ne

is a scalar variable of type ipc_, that holds the number of entries in the lower triangular part of the Hessian matrix \(H\).

H_val

is a one-dimensional array of size h_ne and type rpc_, that holds the values of the entries of the lower triangular part of the Hessian matrix \(H\) in any of the available storage schemes.

x

is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\).

M_ne

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

M_val

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

m

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

A_ne

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

A_val

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

y

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

void trs_information(void **data, struct trs_inform_type* inform, ipc_ *status)

Provides output information

Parameters:

data

holds private internal data

inform

is a struct containing output information (see trs_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 trs_terminate(
    void **data,
    struct trs_control_type* control,
    struct trs_inform_type* inform
)

Deallocate all internal private storage

Parameters:

data

holds private internal data

control

is a struct containing control information (see trs_control_type)

inform

is a struct containing output information (see trs_inform_type)

available structures#

trs_control_type structure#

#include <galahad_trs.h>

struct trs_control_type {
    // fields

    bool f_indexing;
    ipc_ error;
    ipc_ out;
    ipc_ problem;
    ipc_ print_level;
    ipc_ dense_factorization;
    ipc_ new_h;
    ipc_ new_m;
    ipc_ new_a;
    ipc_ max_factorizations;
    ipc_ inverse_itmax;
    ipc_ taylor_max_degree;
    rpc_ initial_multiplier;
    rpc_ lower;
    rpc_ upper;
    rpc_ stop_normal;
    rpc_ stop_absolute_normal;
    rpc_ stop_hard;
    rpc_ start_invit_tol;
    rpc_ start_invitmax_tol;
    bool equality_problem;
    bool use_initial_multiplier;
    bool initialize_approx_eigenvector;
    bool force_Newton;
    bool space_critical;
    bool deallocate_error_fatal;
    char problem_file[31];
    char symmetric_linear_solver[31];
    char definite_linear_solver[31];
    char prefix[31];
    struct sls_control_type sls_control;
    struct ir_control_type ir_control;
};

detailed documentation#

control derived type as a C struct

components#

bool f_indexing

use C or Fortran sparse matrix indexing

ipc_ error

unit for error messages

ipc_ out

unit for monitor output

ipc_ problem

unit to write problem data into file problem_file

ipc_ print_level

controls level of diagnostic output

ipc_ dense_factorization

should the problem be solved by dense factorization? Possible values are

  • 0 sparse factorization will be used

  • 1 dense factorization will be used

  • other the choice is made automatically depending on the dimension & sparsity

ipc_ new_h

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

  • 0 unchanged

  • 1 values but not indices have changed

  • 2 values and indices have changed

ipc_ new_m

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

  • 0 unchanged

  • 1 values but not indices have changed

  • 2 values and indices have changed

ipc_ new_a

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

  • 0 unchanged

  • 1 values but not indices have changed

  • 2 values and indices have changed

ipc_ max_factorizations

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

ipc_ inverse_itmax

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

ipc_ taylor_max_degree

maximum degree of Taylor approximant allowed

rpc_ initial_multiplier

initial estimate of the Lagrange multipler

rpc_ lower

lower and upper bounds on the multiplier, if known

rpc_ upper

see lower

rpc_ stop_normal

stop when \(| ||x|| - radius | \leq\) max( stop_normal \* radius, stop_absolute_normal )

rpc_ stop_absolute_normal

see stop_normal

rpc_ stop_hard

stop when bracket on optimal multiplier \(\leq\) stop_hard \* max( bracket ends )

rpc_ start_invit_tol

start inverse iteration when bracket on optimal multiplier \(\leq\) stop_start_invit_tol \* max( bracket ends )

rpc_ start_invitmax_tol

start full inverse iteration when bracket on multiplier \(\leq\) stop_start_invitmax_tol \* max( bracket ends)

bool equality_problem

is the solution is <b<required to lie on the boundary (i.e., is the constraint an equality)?

bool use_initial_multiplier

ignore initial_multiplier?

bool initialize_approx_eigenvector

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

bool force_Newton

ignore the trust-region if \(H\) is positive definite

bool space_critical

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

bool deallocate_error_fatal

exit if any deallocation fails

char problem_file[31]

name of file into which to write problem data

char symmetric_linear_solver[31]

the name of the symmetric-indefinite linear equation solver used. Possible choices are currently: ‘sils’, ‘ma27’, ‘ma57’, ‘ma77’, ‘ma86’, ‘ma97’, ‘ssids’, ‘mumps’, ‘pardiso’, ‘mkl_pardiso’, ‘pastix’, ‘wsmp’, and ‘sytr’, although only ‘sytr’ and, for OMP 4.0-compliant compilers, ‘ssids’ are installed by default; others are easily installed (see README.external). More details of the capabilities of each solver are provided in the documentation for galahad_sls.

char definite_linear_solver[31]

the name of the definite linear equation solver used. Possible choices are currently: ‘sils’, ‘ma27’, ‘ma57’, ‘ma77’, ‘ma86’, ‘ma87’, ‘ma97’, ‘ssids’, ‘mumps’, ‘pardiso’, ‘mkl_pardiso’, ‘pastix’, ‘wsmp’, ‘potr’, ‘sytr’ and ‘pbtr’, although only ‘potr’, ‘sytr’, ‘pbtr’ and, for OMP 4.0-compliant compilers, ‘ssids’ are installed by default; others are easily installed (see README.external). More details of the capabilities of each solver are provided in the documentation for galahad_sls.

char prefix[31]

all output lines will be prefixed by prefix(2:LEN(TRIM(.prefix))-1) where prefix contains the required string enclosed in quotes, e.g. “string” or ‘string’

struct sls_control_type sls_control

control parameters for the Cholesky factorization and solution (see sls_c documentation)

struct ir_control_type ir_control

control parameters for iterative refinement (see ir_c documentation)

trs_time_type structure#

#include <galahad_trs.h>

struct trs_time_type {
    // fields

    rpc_ total;
    rpc_ assemble;
    rpc_ analyse;
    rpc_ factorize;
    rpc_ solve;
    rpc_ clock_total;
    rpc_ clock_assemble;
    rpc_ clock_analyse;
    rpc_ clock_factorize;
    rpc_ clock_solve;
};

detailed documentation#

time derived type as a C struct

components#

rpc_ total

total CPU time spent in the package

rpc_ assemble

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

rpc_ analyse

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

rpc_ factorize

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

rpc_ solve

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

rpc_ clock_total

total clock time spent in the package

rpc_ clock_assemble

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

rpc_ clock_analyse

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

rpc_ clock_factorize

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

rpc_ clock_solve

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

trs_history_type structure#

#include <galahad_trs.h>

struct trs_history_type {
    // fields

    rpc_ lambda;
    rpc_ x_norm;
};

detailed documentation#

history derived type as a C struct

components#

rpc_ lambda

the value of \(\lambda\)

rpc_ x_norm

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

trs_inform_type structure#

#include <galahad_trs.h>

struct trs_inform_type {
    // fields

    ipc_ status;
    ipc_ alloc_status;
    ipc_ factorizations;
    int64_t max_entries_factors;
    ipc_ len_history;
    rpc_ obj;
    rpc_ x_norm;
    rpc_ multiplier;
    rpc_ pole;
    bool dense_factorization;
    bool hard_case;
    char bad_alloc[81];
    struct trs_time_type time;
    struct trs_history_type history[100];
    struct sls_inform_type sls_inform;
    struct ir_inform_type ir_inform;
};

detailed documentation#

inform derived type as a C struct

components#

ipc_ status

reported return status:

  • 0

    the solution has been found

  • -1

    an array allocation has failed

  • -2

    an array deallocation has failed

  • -3

    n and/or Delta is not positive

  • -9

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

  • -10

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

  • -15

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

  • -16

    ill-conditioning has prevented further progress

ipc_ alloc_status

STAT value after allocate failure.

ipc_ factorizations

the number of factorizations performed

int64_t max_entries_factors

the maximum number of entries in the factors

ipc_ len_history

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

rpc_ obj

the value of the quadratic function

rpc_ x_norm

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

rpc_ multiplier

the Lagrange multiplier corresponding to the trust-region constraint

rpc_ pole

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

bool dense_factorization

was a dense factorization used?

bool hard_case

has the hard case occurred?

char bad_alloc[81]

name of array that provoked an allocate failure

struct trs_time_type time

time information

struct trs_history_type history[100]

history information

struct sls_inform_type sls_inform

cholesky information (see sls_c documentation)

struct ir_inform_type ir_inform

iterative_refinement information (see ir_c documentation)

example calls#

This is an example of how to use the package to solve a trust-region subproblem; the code is available in $GALAHAD/src/trs/C/trst.c . A variety of supported Hessian and constraint matrix storage formats are shown.

Notice that C-style indexing is used, and that this is flagged by setting control.f_indexing to false. The floating-point type rpc_ is set in galahad_precision.h to double by default, but to float if the preprocessor variable SINGLE is defined. Similarly, the integer type ipc_ from galahad_precision.h is set to int by default, but to int64_t if the preprocessor variable INTEGER_64 is defined.

/* trst.c */
/* Full test for the TRS C interface using C sparse matrix indexing */

#include <stdio.h>
#include <string.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_trs.h"

int main(void) {

    // Derived types
    void *data;
    struct trs_control_type control;
    struct trs_inform_type inform;

    // Set problem data
    ipc_ n = 3; // dimension of H
    ipc_ m = 1; // dimension of A
    ipc_ H_ne = 4; // number of elements of H
    ipc_ M_ne = 3; // number of elements of M
    ipc_ A_ne = 3; // number of elements of A
    ipc_ H_dense_ne = 6; // number of elements of H
    ipc_ M_dense_ne = 6; // number of elements of M
    ipc_ H_row[] = {0, 1, 2, 2}; // row indices, NB lower triangle
    ipc_ H_col[] = {0, 1, 2, 0};
    ipc_ H_ptr[] = {0, 1, 2, 4};
    ipc_ M_row[] = {0, 1, 2}; // row indices, NB lower triangle
    ipc_ M_col[] = {0, 1, 2};
    ipc_ M_ptr[] = {0, 1, 2, 3};
    ipc_ A_row[] = {0, 0, 0} ;
    ipc_ A_col[] = {0, 1, 2};
    ipc_ A_ptr[] = {0, 3};
    rpc_ H_val[] = {1.0, 2.0, 3.0, 4.0};
    rpc_ M_val[] = {1.0, 2.0, 1.0};
    rpc_ A_val[] = {1.0, 1.0, 1.0};
    rpc_ H_dense[] = {1.0, 0.0, 2.0, 4.0, 0.0, 3.0};
    rpc_ M_dense[] = {1.0, 0.0, 2.0, 0.0, 0.0, 1.0};
    rpc_ H_diag[] = {1.0, 0.0, 2.0};
    rpc_ M_diag[] = {1.0, 2.0, 1.0};
    rpc_ f = 0.96;
    rpc_ radius = 1.0;
    rpc_ c[] = {0.0, 2.0, 0.0};

    char st = ' ';
    ipc_ status;
    rpc_ x[n];
    char ma[3];

    printf(" C sparse matrix indexing\n\n");

    printf(" basic tests of storage formats\n\n");

    for( ipc_ a_is=0; a_is <= 1; a_is++){ // add a linear constraint?
      for( ipc_ m_is=0; m_is <= 1; m_is++){ // include a scaling matrix?

        if (a_is == 1 && m_is == 1 ) {
          strcpy(ma, "MA");
        }
        else if (a_is == 1) {
          strcpy(ma, "A ");
        }
        else if (m_is == 1) {
          strcpy(ma, "M ");
        }
        else {
          strcpy(ma, "  ");
        }

        for( ipc_ storage_type=1; storage_type <= 4; storage_type++){

          // Initialize TRS
          trs_initialize( &data, &control, &status );

          // Set user-defined control options
          control.f_indexing = false; // C sparse matrix indexing

          switch(storage_type){
              case 1: // sparse co-ordinate storage
                  st = 'C';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                             "coordinate", H_ne, H_row, H_col, NULL );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                  "coordinate", M_ne, M_row, M_col, NULL );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                  "coordinate", A_ne, A_row, A_col, NULL );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, m, A_ne, A_val, NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              case 2: // sparse by rows
                  st = 'R';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                              "sparse_by_rows", H_ne, NULL, H_col, H_ptr );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                  "sparse_by_rows", M_ne, NULL, M_col, M_ptr );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                 "sparse_by_rows", A_ne, NULL, A_col, A_ptr );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, m, A_ne, A_val, NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              case 3: // dense
                  st = 'D';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                              "dense", H_ne, NULL, NULL, NULL );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                 "dense", M_ne, NULL, NULL, NULL );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                 "dense", A_ne, NULL, NULL, NULL );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       M_dense_ne, M_dense, m, A_ne, A_val,
                                       NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       M_dense_ne, M_dense, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              case 4: // diagonal
                  st = 'L';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                              "diagonal", H_ne, NULL, NULL, NULL );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                 "diagonal", M_ne, NULL, NULL, NULL );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                 "dense", A_ne, NULL, NULL, NULL );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       n, M_diag, m, A_ne, A_val, NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       n, M_diag, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              }

          trs_information( &data, &inform, &status );

          printf("format %c%s: TRS_solve_problem exit status = %1" i_ipc_ ", f = %.2f\n",
                 st, ma, inform.status, inform.obj );
          //printf("x: ");
          //for( ipc_ i = 0; i < n+m; i++) printf("%f ", x[i]);

          // Delete internal workspace
          trs_terminate( &data, &control, &inform );
       }
     }
   }
}

This is the same example, but now fortran-style indexing is used; the code is available in $GALAHAD/src/trs/C/trstf.c .

/* trstf.c */
/* Full test for the TRS C interface using Fortran sparse matrix indexing */

#include <stdio.h>
#include <string.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_trs.h"

int main(void) {

    // Derived types
    void *data;
    struct trs_control_type control;
    struct trs_inform_type inform;

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

    char st = ' ';
    ipc_ status;
    rpc_ x[n];
    char ma[3];

    printf(" Fortran sparse matrix indexing\n\n");

    printf(" basic tests of storage formats\n\n");

    for( ipc_ a_is=0; a_is <= 1; a_is++){ // add a linear constraint?
      for( ipc_ m_is=0; m_is <= 1; m_is++){ // include a scaling matrix?

        if (a_is == 1 && m_is == 1 ) {
          strcpy(ma, "MA");
        }
        else if (a_is == 1) {
          strcpy(ma, "A ");
        }
        else if (m_is == 1) {
          strcpy(ma, "M ");
        }
        else {
          strcpy(ma, "  ");
        }

        for( ipc_ storage_type=1; storage_type <= 4; storage_type++){

          // Initialize TRS
          trs_initialize( &data, &control, &status );

          // Set user-defined control options
          control.f_indexing = true; // fortran sparse matrix indexing

          switch(storage_type){
              case 1: // sparse co-ordinate storage
                  st = 'C';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                             "coordinate", H_ne, H_row, H_col, NULL );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                  "coordinate", M_ne, M_row, M_col, NULL );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                  "coordinate", A_ne, A_row, A_col, NULL );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, m, A_ne, A_val, NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              case 2: // sparse by rows
                  st = 'R';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                              "sparse_by_rows", H_ne, NULL, H_col, H_ptr );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                  "sparse_by_rows", M_ne, NULL, M_col, M_ptr );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                 "sparse_by_rows", A_ne, NULL, A_col, A_ptr );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, m, A_ne, A_val, NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       M_ne, M_val, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_ne, H_val, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              case 3: // dense
                  st = 'D';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                              "dense", H_ne, NULL, NULL, NULL );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                 "dense", M_ne, NULL, NULL, NULL );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                 "dense", A_ne, NULL, NULL, NULL );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       M_dense_ne, M_dense, m, A_ne, A_val,
                                       NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       M_dense_ne, M_dense, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, H_dense_ne, H_dense, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              case 4: // diagonal
                  st = 'L';
                  // import the control parameters and structural data
                  trs_import( &control, &data, &status, n,
                              "diagonal", H_ne, NULL, NULL, NULL );
                  if (m_is == 1) {
                    trs_import_m( &data, &status, n,
                                 "diagonal", M_ne, NULL, NULL, NULL );
                  }
                  if (a_is == 1) {
                    trs_import_a( &data, &status, m,
                                 "dense", A_ne, NULL, NULL, NULL );
                  }
                  // solve the problem
                  if (a_is == 1 && m_is == 1 ) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       n, M_diag, m, A_ne, A_val, NULL );
                  }
                  else if (a_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       0, NULL, m, A_ne, A_val, NULL );
                  }
                  else if (m_is == 1) {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       n, M_diag, 0, 0, NULL, NULL );
                  }
                  else {
                    trs_solve_problem( &data, &status, n,
                                       radius, f, c, n, H_diag, x,
                                       0, NULL, 0, 0, NULL, NULL );
                  }
                  break;
              }

          trs_information( &data, &inform, &status );

          printf("format %c%s: TRS_solve_problem exit status = %1" i_ipc_ ", f = %.2f\n",
                 st, ma, inform.status, inform.obj );
          //printf("x: ");
          //for( ipc_ i = 0; i < n+m; i++) printf("%f ", x[i]);

          // Delete internal workspace
          trs_terminate( &data, &control, &inform );
       }
     }
   }
}