GALAHAD TREK package#

purpose#

The trek package uses an extended-Krylov-subspace iteration 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 + c^T x + \frac{1}{2} x^T H x,\]
where the vector \(x\) is required to satisfy the ellipsoidal trust-region constraint \(\|x\|_S \leq \Delta\), where the \(S\)-norm of \(x\) is defined to be \(\|x\|_S = \sqrt{x^T S x}\), and where the radius \(\Delta > 0\). The matrix \(S\) need not be provided in the commonly-occurring \(\ell_2\)-trust-region case for which \(S = I\), the \(n\) by \(n\) identity matrix.

Factorization of the matrices \(H\) and, if present, \(S\) 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/trek.pdf for additional details.

method#

The required solution \(x_*\) necessarily satisfies the optimality condition \(H x_* + \lambda_* S x_* + c = 0\), where \(\lambda_* \geq 0\) is a Lagrange multiplier corresponding to the constraint \(\|x\|_M \leq \Delta\). In addition in all cases, the matrix \(H + \lambda_* S\) will be positive semi-definite; in most instances it will actually be positive definite, but in special “hard” cases singularity is a possibility.

The method is iterative, and is based upon building a solution approximation from an orthogonal basis of the evolving extended Krylov subspaces \({\cal K}_{2m+1}(H,c) = \mbox{span}\{c,H^{-1}c,H c,H^{-2}c,H^2c,\ldots,\) \(H^{-m}c,H^{m}c\}\) as \(m\) increases. The key observations are (i) the manifold of solutions to the optimality system \[ ( H + \lambda I ) x(\lambda) = - c\] as a function of \(\sigma\) is of approximately very low rank, (ii) the subspace \({\cal K}_{2m+1}(H,c)\) rapidly gives a very good approximation to this manifold, (iii) it is straightforward to build an orthogonal basis of \({\cal K}_{2m+1}(H,c)\) using short-term recurrences and a single factorization of \(H\), and (iv) solutions to the trust-region subproblem restricted to elements of the orthogonal subspace may be found very efficiently using effective high-order root-finding methods. The fact that the second element in the subspace is \(H^{-1} c\) means that it is easy to check for the interior-solution possibility \(x = - H^{-1} c\) that occurs when such a \(x\) satisfies \(\|x\| \leq \Delta\). Coping with general scalings \(S\) is a straightforward extension so long as factorization of \(S\) is also possible.

reference#

The method is described in detail in

H. Al Daas and N. I. M. Gould. Extended-Krylov-subspace methods for trust-region and norm-regularization subproblems. Preprint STFC-P-2025-002, Rutherford Appleton Laboratory, Oxfordshire, England.

matrix storage#

symmetric storage#

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

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

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

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

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

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

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

introduction to function calls#

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

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

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

  • trek_import - set up problem data structures and fixed values

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

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

  • trek_solve_problem - solve the trust-region problem

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

  • trek_terminate - deallocate data structures

See the examples section for illustrations of use.

parametric real type T and integer type INT#

Below, the symbol T refers to a parametric real type that may be Float32 (single precision), Float64 (double precision) or, if supported, Float128 (quadruple precision). The symbol INT refers to a parametric integer type that may be Int32 (32-bit integer) or Int64 (64-bit integer).

callable functions#

    function trek_initialize(T, INT, data, control, status)

Set default control values and initialize private data

Parameters:

data

holds private internal data

control

is a structure containing control information (see trek_control_type)

status

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

  • 0

    The initialization was successful.

    function trek_read_specfile(T, INT, control, specfile)

Read the content of a specification file, and assign values associated with given keywords to the corresponding control parameters. An in-depth discussion of specification files is available, and a detailed list of keywords with associated default values is provided in $GALAHAD/src/trek/TREK.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/trek.pdf for a list of how these keywords relate to the components of the control structure.

Parameters:

control

is a structure containing control information (see trek_control_type)

specfile

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

    function trek_import(T, INT, control, data, status, n,
                         H_type, H_ne, H_row, H_col, H_ptr)

Import problem data into internal storage prior to solution.

Parameters:

control

is a structure whose members provide control parameters for the remaining procedures (see trek_control_type)

data

holds private internal data

status

is a scalar variable of type INT 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 restriction n > 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 INT that holds the number of rows (and columns) of H.

H_type

is a one-dimensional array of type Vararg{Cchar} that specifies the symmetric storage scheme used for the Hessian, \(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_ne

is a scalar variable of type INT 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 INT that holds the row indices of the lower triangular part of \(H\) in the sparse co-ordinate storage scheme. It need not be set for any of the other three schemes, and in this case can be C_NULL.

H_col

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

H_ptr

is a one-dimensional array of size n+1 and type INT that holds the starting position of each row of the lower triangular part of \(H\), as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the other schemes are used, and in this case can be C_NULL.

    function trek_import_s(T, INT, data, status, n, S_type, S_ne, S_row, S_col, S_ptr)

Import data for the scaling matrix \(S\) into internal storage prior to solution.

Parameters:

data

holds private internal data

status

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

n

is a scalar variable of type INT that holds the number of rows (and columns) of \(S\).

S_type

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

S_ne

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

S_row

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

S_col

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

S_ptr

is a one-dimensional array of size n+1 and type INT that holds the starting position of each row of the lower triangular part of \(S\), as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the other schemes are used, and in this case can be C_NULL.

    function trek_reset_control(T, INT, control, data, status)

Reset control parameters after import if required.

Parameters:

control

is a structure whose members provide control parameters for the remaining procedures (see trek_control_type)

data

holds private internal data

status

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

  • 0

    The import was successful.

    function trek_solve_problem(T, INT, data, status, n, H_ne, H_val, c, radius,
                                x, S_ne, S_val)

Solve the trust-region subproblem.

Parameters:

data

holds private internal data

status

is a scalar variable of type INT 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, or radius > 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

    The analysis phase of the factorization of the matrix \(H\) or \(S\) failed.

  • -10

    The factorization of the matrix \(H\) or \(S\) failed.

  • -11

    A solve involving the matrix \(H\) or \(S\) failed.

  • -15

    The matrix \(S\) appears not to be diagonally dominant.

  • -16

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

  • -18

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

    • -31 A resolve call has been made before an initial call (see control.new_radius and control.new_values).

    • -38 An error occurred in a call to an LAPACK subroutine.

n

is a scalar variable of type INT that holds the number of variables

H_ne

is a scalar variable of type INT 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 T that holds the values of the entries of the lower triangular part of the Hessian matrix \(H\) in any of the available storage schemes.

c

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

radius

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

x

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

S_ne

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

S_val

is a one-dimensional array of size S_ne and type T that holds the values of the entries of the scaling matrix \(S\), if it is not the identity matrix, in any of the available storage schemes. If S_val is C_NULL, \(S\) will be taken to be the identity matrix.

    function trek_information(T, INT, data, inform, status)

Provides output information

Parameters:

data

holds private internal data

inform

is a structure containing output information (see trek_inform_type)

status

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

  • 0

    The values were recorded successfully

    function trek_terminate(T, INT, data, control, inform)

Deallocate all internal private storage

Parameters:

data

holds private internal data

control

is a structure containing control information (see trek_control_type)

inform

is a structure containing output information (see trek_inform_type)

available structures#

trek_control_type structure#

    struct trek_control_type{T,INT}
      f_indexing::Bool
      error::INT
      out::INT
      print_level::INT
      eks_max::INT
      it_max::INT
      f::T
      reduction::T
      stop_residual::T
      reorthogonalize::Bool
      s_version_52::Bool
      perturb_c::Bool
      stop_check_all_orders::Bool
      new_radius::Bool
      new_values::Bool
      space_critical::Bool
      deallocate_error_fatal::Bool
      linear_solver::NTuple{31,Cchar}
      linear_solver_for_s::NTuple{31,Cchar}
      prefix::NTuple{31,Cchar}
      sls_control::sls_control_type{T,INT}
      sls_s_control::sls_control_type{T,INT}
      trs_control::ir_control_type{T,INT}

detailed documentation#

control derived type as a Julia structure

components#

Bool f_indexing

use C or Fortran sparse matrix indexing

INT error

unit for error messages

INT out

unit for monitor output

INT print_level

controls level of diagnostic output

INT eks_max

maximum dimension of the extended Krylov space employed. If a negative value is given, the value 100 will be used instead

INT it_max

the maximum number of iterations allowed. If a negative value is given, the value 100 will be used instead

T f

the value of \(f\) in the objective function. This value has no effect on the computed \(x\), and takes the value 0.0 by default

T reduction

the value of the reduction factor for a suggested subsequent trust-region radius, see control[‘next_radius’]. The suggested radius will be reduction times the smaller of the current radius and \(\|x\|_S\) at the output \(x\).

T stop_residual

the value of the stopping tolerance used by the algorithm. The iteration stops as soon as \(x\) and \(\lambda\) are found to satisfy \(\| ( H + \lambda S ) x + c \| <\) stop_residual \(\times \max( 1, \|c\| )\).

Bool reorthogonalize

should be set to true if the generated basis of the extended-Krylov subspace is to be reorthogonalized at every iteration. This can be very expensive, and is generally not warranted

Bool s_version_52

should be set to true if Algorithm 5.2 in the paper is used to generate the extended Krylov space recurrences when a non-unit \(S\) is given, and false if those from Algorithm B.3 ares used instead. In practice, there is very little difference in performance and accuracy

Bool perturb_c

should be set to true if the user wishes to make a tiny pseudo-random perturbations to the components of the term \(c\) to try to protect from the so-called (probability zero) “hard” case. Perturbations are generally not needed, and should only be used in very exceptional cases

Bool stop_check_all_orders

should be set to true if the algorithm checks for termination for each new member of the extended Krylov space. Such checks incur some extra cost, and experience shows that testing every second member is sufficient

Bool new_radius

should be set to true if the call retains the previous \(H\), \(S\) and \(c\), but with a new, smaller radius

Bool new_values

should be set to true if the any of the values of \(H\), \(S\) and \(c\) has changed since a previous call

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

the name of the linear equation solver used to solve any symmetric positive-definite linear system involving \(H\) that might arise. Possible choices are currently: ‘sils’, ‘ma27’, ‘ma57’, ‘ma77’, ‘ma86’, ‘ma87’, ‘ma97’, ssids, ‘pardiso’, ‘wsmp’, ‘sytr’, ‘potr’ and ‘pbtr’ although only ‘sytr’, ‘potr’, ‘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 linear_solver_for_s[31]

the name of the linear equation solver used to solve any symmetric positive-definite linear system involving the optional \(S\) that might arise. Possible choices are currently: ‘sils’, ‘ma27’, ‘ma57’, ‘ma77’, ‘ma86’, ‘ma87’, ‘ma97’, ssids, ‘pardiso’, ‘wsmp’, ‘sytr’, ‘potr’ and ‘pbtr’ although only ‘sytr’, ‘potr’, ‘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.

NTuple{31,Cchar} prefix

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

struct sls_control_type sls_control

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

struct sls_control_type sls_s_control

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

struct trs_control_type trs_control

control parameters for the solution of diagonal trust-region subproblems (see trs_c documentation)

trek_time_type structure#

    struct trek_time_type{T}
      total::T
      assemble::T
      analyse::T
      factorize::T
      solve::T
      clock_total::T
      clock_assemble::T
      clock_analyse::T
      clock_factorize::T
      clock_solve::T

detailed documentation#

time derived type as a Julia structure

components#

T total

total CPU time spent in the package

T assemble

CPU time spent building \(H\) and \(S\).

T analyse

CPU time spent reordering \(H\) and \(S\) prior to factorization.

T factorize

CPU time spent factorizing \(H\) and \(S\).

T solve

CPU time spent solving linear systems inolving \(H\) and \(S\).

T clock_total

total clock time spent in the package

T clock_assemble

clock time spent building \(H\) and \(S\)

T clock_analyse

clock time spent reordering \(H\) and \(S\) prior to factorization

T clock_factorize

clock time spent factorizing \(H\) and \(S\)

T clock_solve

clock time spent solving linear systems inolving \(H\) and \(S\)

trek_inform_type structure#

    struct trek_inform_type{T,INT}
      status::INT
      alloc_status::INT
      iter::INT
      n_vec::INT
      obj::T
      x_norm::T
      multiplier::T
      radius::T
      next_radius::T
      error::T
      bad_alloc::NTuple{81,Cchar}
      time::trek_time_type{T}
      sls_inform::sls_inform_type{T,INT}
      sls_s_inform::sls_inform_type{T,INT}
      trs_inform::trs_inform_type{T,INT}

detailed documentation#

inform derived type as a Julia structure

components#

INT 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, or the 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.sls_inform.status or inform.sls_s_inform.status as appropriate

  • -10

    the factorization failed; the return status from the factorization package is given by inform.sls_inform.status or inform.sls_s_inform.status as appropriate

  • -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.sls_inform.status or inform.sls_s_inform.status as appropriate

  • -15

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

  • -16

    ill-conditioning has prevented further progress

  • -18

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

  • -31

    a resolve call has been made before an initial call (see control.new_radius and control.new_values)

  • -38

    an error occurred in a call to an LAPACK subroutine

INT alloc_status

STAT value after allocate failure.

INT iter

the total number of iterations required

INT n_vec

the number of orthogonal vectors required (the dimension of the extended-Krylov subspace)

T obj

the value of the quadratic function

T x_norm

the \(S\) -norm of \(x\), \(||x||_S\)

T multiplier

the Lagrange multiplier corresponding to the trust-region constraint

T radius

the value of the current radius

T next_radius

the value of the the proposed next radius to be used if the current radius proves to be too large (see inform.reduction).

T error

the value of the norm of the maximum relative residual error, \(\|(H+\lambda S) x + c\|/\max(1,\|c\|)\)

NTuple{81,Cchar} bad_alloc

name of array that provoked an allocate failure

struct trek_time_type time

time information

struct sls_inform_type sls_inform

Cholesky information for factorization and solves with \(H\) (see sls_c documentation)

struct sls_inform_type sls_s_inform

Cholesky information for factorization and solves with \(S\) (see sls_c documentation)

struct trs_inform_type trs_inform

diagonal trust-region solve information (see trs_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/trek/Julia/test_trek.jl . A variety of supported Hessian and scaling matrix storage formats are shown.

# test_trek.jl
# Simple code to test the Julia interface to TREK

using GALAHAD
using Test
using Printf
using Accessors
using Quadmath

function test_trek(::Type{T}, ::Type{INT}; dls::String="pbtr") where {T,INT}
  # Derived types
  data = Ref{Ptr{Cvoid}}()
  control = Ref{trek_control_type{T,INT}}()
  inform = Ref{trek_inform_type{T,INT}}()

  # Set problem data
  n = INT(3)  # dimension of H
  m = INT(1)  # dimension of A
  H_ne = INT(4)  # number of elements of H
  S_ne = INT(3)  # number of elements of M
  H_dense_ne = INT(6)  # number of elements of H
  S_dense_ne = INT(6)  # number of elements of M
  H_row = INT[1, 2, 3, 3]  # row indices, NB lower triangle
  H_col = INT[1, 2, 3, 1]
  H_ptr = INT[1, 2, 3, 5]
  S_row = INT[1, 2, 3]  # row indices, NB lower triangle
  S_col = INT[1, 2, 3]
  S_ptr = INT[1, 2, 3, 4]
  H_val = T[1.0, 2.0, 3.0, 4.0]
  S_val = T[1.0, 2.0, 1.0]
  H_dense = T[1.0, 0.0, 2.0, 4.0, 0.0, 3.0]
  S_dense = T[1.0, 0.0, 2.0, 0.0, 0.0, 1.0]
  H_diag = T[1.0, 0.0, 2.0]
  S_diag = T[1.0, 2.0, 1.0]
  c = T[0.0, 2.0, 0.0]

  st = ' '
  status = Ref{INT}()
  x = zeros(T, n)
  sr = ""

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

  for s_is in 0:1 # include a scaling matrix?

    for storage_type in 1:4
      # Initialize TREK
      trek_initialize(T, INT, data, control, status)

      # Linear solvers
      @reset control[].linear_solver = galahad_linear_solver(dls)
      @reset control[].linear_solver_for_s = galahad_linear_solver(dls)

      # sparse co-ordinate storage
      if storage_type == 1
        st = 'C'
        # import the control parameters and structural data
        trek_import(T, INT, control, data, status, n,
                   "coordinate", H_ne, H_row, H_col, C_NULL)
        if s_is == 1
          trek_s_import(T, INT, data, status, n,
                        "coordinate", S_ne, S_row, S_col, C_NULL)
        end
      end

      # sparse by rows
      if storage_type == 2
        st = 'R'
        # import the control parameters and structural data
        trek_import(T, INT, control, data, status, n,
                   "sparse_by_rows", H_ne, C_NULL, H_col, H_ptr)
        if s_is == 1
          trek_s_import(T, INT, data, status, n,
                        "sparse_by_rows", S_ne, C_NULL, S_col, S_ptr)
        end
      end

      # dense
      if storage_type == 3
        st = 'D'
        # import the control parameters and structural data
        trek_import(T, INT, control, data, status, n,
                   "dense", H_ne, C_NULL, C_NULL, C_NULL)
        if s_is == 1
          trek_s_import(T, INT, data, status, n,
                        "dense", S_ne, C_NULL, C_NULL, C_NULL)
        end
      end

      # diagonal
      if storage_type == 4
        st = 'L'
        # import the control parameters and structural data
        trek_import(T, INT, control, data, status, n,
                   "diagonal", H_ne, C_NULL, C_NULL, C_NULL)
        if s_is == 1
          trek_s_import(T, INT, data, status, n,
                        "diagonal", S_ne, C_NULL, C_NULL, C_NULL)
        end
      end

      for r_is in 1:2 # original or smaller radius

        if (r_is == 1)
          radius = one(T)
        else
          radius = inform[].next_radius
          @reset control[].new_radius = true
          trek_reset_control(T, INT, control, data, status)
        end

        if (r_is == 2) && (s_is == 1)
          sr = "S-"
        elseif r_is == 2
          sr = "- "
        elseif s_is == 1
          sr = "S "
        else
          sr = "  "
        end

        # solve the problem

        # sparse co-ordinate storage
        if storage_type == 1
          if s_is == 1
            trek_solve_problem(T, INT, data, status, n,
                               H_ne, H_val, c, radius, x,
                               S_ne, S_val, 0, 0, C_NULL, C_NULL)
          else
            trek_solve_problem(T, INT, data, status, n,
                               H_ne, H_val, c, radius, x,
                               0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        # sparse by rows
        if storage_type == 2
          if s_is == 1
            trek_solve_problem(T, INT, data, status, n,
                               H_ne, H_val, c, radius, x,
                               S_ne, S_val, 0, 0, C_NULL, C_NULL)
          else
            trek_solve_problem(T, INT, data, status, n,
                               H_ne, H_val, c, radius, x,
                               0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        # dense
        if storage_type == 3
          if s_is == 1
            trek_solve_problem(T, INT, data, status, n,
                               H_dense_ne, H_dense_val, c, radius, x,
                               S_dense_ne, S_dense, 0, 0, C_NULL, C_NULL)
          else
            trek_solve_problem(T, INT, data, status, n,
                               H_dense_ne, H_dense_val, c, radius, x,
                               0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        # diagonal
        if storage_type == 4
          if s_is == 1
            trek_solve_problem(T, INT, data, status, n,
                               n, H_diah_val, c, radius, x,
                               n, S_diag, 0, 0, C_NULL, C_NULL)
          else
            trek_solve_problem(T, INT, data, status, n,
                               n, H_diah_val, c, radius, x,
                               0, C_NULL, 0, 0, C_NULL, C_NULL)
          end
        end

        trek_information(T, INT, data, inform, status)

        @printf("format %c%s: TREK_solve_problem exit status = %1i,
                f = %.2f\n", st, sr, inform[].status, inform[].obj)

      end

      # Delete internal workspace
      trek_terminate(T, INT, data, control, inform)

    end
  end

  return 0
end

for (T, INT, libgalahad) in ((Float32 , Int32, GALAHAD.libgalahad_single      ),
                             (Float32 , Int64, GALAHAD.libgalahad_single_64   ),
                             (Float64 , Int32, GALAHAD.libgalahad_double      ),
                             (Float64 , Int64, GALAHAD.libgalahad_double_64   ),
                             (Float128, Int32, GALAHAD.libgalahad_quadruple   ),
                             (Float128, Int64, GALAHAD.libgalahad_quadruple_64))
  if isfile(libgalahad)
    @testset "TREK -- $T -- $INT" begin
      @test test_trek(T, INT) == 0
    end
end