GALAHAD EXPO package#

purpose#

The expo package uses an exponential-penalty function method to solve a given constrained optimization problem. The aim is to find a (local) minimizer of a differentiable objective function \(f(x)\) of \(n\) variables \(x\), subject to \(m\) general constraints \(c_l \leq c(x) \leq c_u\) and simple-bound constraints \(x_l \leq x \leq x_u\) on the variables. Here, any of the components of the vectors of bounds \(c_l\), \(c_u\), \(x_l\) and \(x_u\) may be infinite. The method offers the choice of direct and iterative solution of the key unconstrained-optimization subproblems, and is most suitable for large problems. First derivatives are required, and if second derivatives can be calculated, they will be exploited—if the product of second derivatives with a vector may be found but not the derivatives themselves, that may also be exploited.

N.B. This package is currently a beta release, and aspects may change before it is formally released

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

terminology#

The exponential penalty function is defined to be

\[\begin{split}\begin{array}{rl}\phi(x,w,\mu,v,\nu) \!\! & = f(x) + \sum_{i} \mu_{li} w_{li} \exp[(c_{li} - c_i(x))/\mu_{li}] \\ & \;\;\;\;\;\;\;\;\;\;\;\;\; + \sum_{i} \mu_{ui} w_{ui} \exp[(c_i(x) - c_{ui})/\mu_{ui}] \\ & \;\;\;\;\;\;\;\;\;\;\;\;\; + \sum_{j} \nu_{lj} v_{lj} \exp[(x_{lj} - x_j)/\nu_{lj}] \\ & \;\;\;\;\;\;\;\;\;\;\;\;\; + \sum_{j} \nu_{uj} v_{uj} \exp[(x_j - x_{uj})/\nu_{uj}], \end{array} \end{split}\]
where \(c_{li}\), \(c_{ui}\) and \(c_i(x)\) are the \(i\)-th components of \(c_l\), \(c_u\) and \(c(x)\), and \(c_{lj}\), \(c_{uj}\) and \(x_j\) are the \(j\)-th components of \(x_l\), \(x_u\) and \(x\), respectively. Here the components of \(\mu_l\), \(\mu_u\), \(\nu_l\) and \(\nu_u\) are separate penalty parameters for each lower and upper, general and simple-bound constraint, respectively, while those of \(w_l\), \(w_u\), \(v_l\), \(v_u\) are likewise separate weights for the same. The algorithm iterates by approximately minimizing \(\phi(x,w,\mu,v,\nu)\) for a fixed set of penalty parameters and weights, and then adjusting these parameters and weights. The adjustments are designed so the sequence of approximate minimizers of \(\phi\) converge to that of the specified constrained optimization problem.

Key constructs are the gradient of the objective function

\[g(x) := \nabla_x f(x),\]
the Jacobian of the vector of constraints,
\[J(x) := \nabla_x c(x),\]
and the gradient and Hessian of the Lagrangian function
\[g_L(x,y,z) := g(x) - J^T(x)y - z \;\;\mbox{and}\;\; H_L(x,y) := \nabla_{xx} \left[ f(x) - \sum_{i} y_i c_i(x)\right]\]
for given vectors \(y\) and \(z\).

Any required solution \(x\) necessarily satisfies the primal optimality conditions

\[c_l \leq c(x) \leq c_u, \;\; x_l \leq x \leq x_u,\;\;\mbox{(1)}\]
the dual optimality conditions
\[g(x) = J^{T}(x) y + z,\;\; y = y_l + y_u \;\;\mbox{and}\;\; z = z_l + z_u,\;\;\mbox{(2a)}\]
and
\[y_l \geq 0, \;\; y_u \leq 0, \;\; z_l \geq 0 \;\;\mbox{and}\;\; z_u \leq 0,\;\;\mbox{(2b)}\]
and the complementary slackness conditions
\[( c(x) - c_l )^{T} y_l = 0,\;\; ( c(x) - c_u )^{T} y_u = 0,\;\; (x -x_l )^{T} z_l = 0 \;\;\mbox{and}\;\;(x -x_u )^{T} z_u = 0,\;\;\mbox{(3)}\]
where the vectors \(y\) and \(z\) are known as the Lagrange multipliers for the general constraints, and the dual variables for the simple bounds, respectively, and where the vector inequalities hold component-wise.

method#

The method employed involves a sequential minimization of the exponential penalty function \(\phi(x,w,\mu,v,\nu)\) for a sequence of positive penalty parameters \((\mu_{lk}, \mu_{uk}, \nu_{lk}, \nu_{uk})\) and weights \((w_{lk}, w_{uk}, v_{lk}, v_{uk})\), for increasing \(k \geq 0\). Convergence is ensured if the penalty parameters are forced to zero, and may be accelerated by adjusting the weights. The minimization of \(\phi(x,w,\mu,v,\nu)\) is accomplished using the trust-region unconstrained solver TRU. Although critical points \(\{x_k\}\) of \(\phi(x,w_k,\mu_k,v_k,\nu_k)\) converge to a local solution \(x_*\) of the underlying problem, the reduction of the penalty parameters to zero often results in \(x_k\) being a poor starting point for the minimization of \(\phi(x,w_{k+1},\mu_{k+1},v_{k+1},\nu_{k+1})\). Consequently, a careful extrapolated starting point from \(x_k\) is used instead. Moreover, once the algorithm is confident that it is sufficiently close to \(x_*\), it switches to Newton’s method to accelerate the convergence. Both the extrapolation and the Newton iteration rely on the block-linear-system solver SSLS.

The iteration is terminated as soon as residuals to the optimality conditions (1)–(3) are sufficiently small. For infeasible problems, this will not be possible, and instead the residuals to (1) will be made as small as possible.

references#

The method is described in detail in

N.Gould, S.Leyffer, A.Montoison and C.Vanaret (2025) The exponential multiplier method in the 21st century. RAL Technical Report, in preparation.

matrix storage#

The unsymmetric \(m\) by \(n\) Jacobian matrix \(J = J(x)\) may be presentedand stored in a variety of convenient input formats.

Dense storage format: The matrix \(J\) 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 J_val will hold the value \(J_{ij}\) for \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\).

Dense by columns storage format: The matrix \(J\) 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 J_val will hold the value \(J_{ij}\) for \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\).

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 \(J\), its row index i, column index j and value \(J_{ij}\), \(0 \leq i \leq m-1\), \(0 \leq j \leq n-1\), are stored as the \(l\)-th components of the integer arrays J_row and J_col and real array J_val, respectively, while the number of nonzeros is recorded as J_ne = \(ne\).

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 \(J\) the i-th component of the integer array J_ptr holds the position of the first entry in this row, while J_ptr(m) holds the total number of entries. The column indices j, \(0 \leq j \leq n-1\), and values \(J_{ij}\) of the nonzero entries in the i-th row are stored in components l = J_ptr(i), \(\ldots\), J_ptr(i+1)-1, \(0 \leq i \leq m-1\), of the integer array J_col, and real array J_val, respectively. For sparse matrices, this scheme almost always requires less storage than its predecessor.

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 \(J\) the j-th component of the integer array J_ptr holds the position of the first entry in this column, while J_ptr(n) holds the total number of entries. The row indices i, \(0 \leq i \leq m-1\), and values \(J_{ij}\) of the nonzero entries in the j-th columnsare stored in components l = J_ptr(j), \(\ldots\), J_ptr(j+1)-1, \(0 \leq j \leq n-1\), of the integer array J_row, and real array J_val, respectively. As before, for sparse matrices, this scheme almost always requires less storage than the co-ordinate format.

The symmetric \(n\) by \(n\) matrix \(H = Hl(x,y)\) may 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\).

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.

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.

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.

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 identity matrix format: If \(H\) is the identity matrix, no values need be stored.

The zero matrix format: The same is true if \(H\) is the zero matrix.

introduction to function calls#

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

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

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

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

  • expo_import - set up problem data structures and fixed values

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

  • expo_solve_hessian_direct - solve the problem using function calls to evaluate function and derivative values

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

  • expo_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 expo_control_type;
struct expo_inform_type;
struct expo_time_type;

// function calls

void expo_initialize(
    void **data,
    struct expo_control_type* control,
    struct expo_inform_type* inform
);

void expo_read_specfile(struct expo_control_type* control, const char specfile[]);

void expo_import(
    struct expo_control_type* control,
    void **data,
    ipc_ *status,
    ipc_ n,
    ipc_ m,
    const char J_type[],
    ipc_ J_ne,
    const ipc_ J_row[],
    const ipc_ J_col[],
    const ipc_ J_ptr[],
    const char H_type[],
    ipc_ H_ne,
    const ipc_ H_row[],
    const ipc_ H_col[],
    const ipc_ H_ptr[]
);

void expo_reset_control(
    struct expo_control_type* control,
    void **data,
    ipc_ *status
);

void expo_solve_hessian_direct(
    void **data,
    void *userdata,
    ipc_ *status,
    ipc_ n,
    ipc_ m,
    ipc_ J_ne,
    ipc_ H_ne,
    const rpc_ c_l[],
    const rpc_ c_u[],
    const rpc_ x_l[],
    const rpc_ x_u[],
    rpc_ x[],
    rpc_ y[],
    rpc_ z[],
    rpc_ c[],
    rpc_ gl[],
    ipc_(*)(ipc_, ipc_, const rpc_[], rpc_, rpc_[],
             const void*) eval_fc,
    ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], rpc_[], , rpc_[],
             const void*) eval_gj,
    ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[],
             const void*) eval_hl
);

void expo_information(void **data, struct expo_inform_type* inform, ipc_ *status);

void expo_terminate(
    void **data,
    struct expo_control_type* control,
    struct expo_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 and structure names#

The function and structure names described below are appropriate for the default real working precision (double) and integer word length (int32_t). To use the functions and structures with different precisions and integer word lengths, an additional suffix must be added to their names (and the arguments set accordingly). The appropriate suffices are:

_s for single precision (float) reals and standard 32-bit (int32_t) integers;

_q for quadruple precision (__real128) reals (if supported) and standard 32-bit (int32_t) integers;

_64 for standard precision (double) reals and 64-bit (int64_t) integers;

_s_64 for single precision (float) reals and 64-bit (int64_t) integers; and

_q_64 for quadruple precision (__real128) reals (if supported) and 64-bit (int64_t) integers.

Thus a call to arc_initialize below will instead be

void arc_initialize_s_64(void **data, struct arc_control_type_s_64* control,
                         int64_t *status)

if single precision (float) reals and 64-bit (int64_t) integers are required. Thus it is possible to call functions for this package with more that one precision and/or integer word length at same time. An example illustrates this feature.

function calls#

void expo_initialize(
    void **data,
    struct expo_control_type* control,
    struct expo_inform_type* inform
)

Set default control values and initialize private data

Parameters:

data

holds private internal data

control

is a struct containing control information (see expo_control_type)

inform

is a struct containing output information (see expo_inform_type)

void expo_read_specfile(struct expo_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/expo/EXPO.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/expo.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 expo_control_type)

specfile

is a character string containing the name of the specification file

void expo_import(
    struct expo_control_type* control,
    void **data,
    ipc_ *status,
    ipc_ n,
    ipc_ m,
    const char J_type[],
    ipc_ J_ne,
    const ipc_ J_row[],
    const ipc_ J_col[],
    const ipc_ J_ptr[],
    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 expo_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, and the package is ready for the solve phase

  • -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, m \(\geq\) 0 or requirement that J/H_type contains its relevant string ‘dense’, ‘dense_by_columns’, ‘coordinate’, ‘sparse_by_rows’, ‘sparse_by_columns’, ‘diagonal’ or ‘absent’ has been violated.

n

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

m

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

J_type

is a one-dimensional array of type char that specifies the unsymmetric storage scheme used for the Jacobian, \(J\). 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.

J_ne

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

J_row

is a one-dimensional array of size J_ne and type ipc_, that holds the row indices of \(J\) 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.

J_col

is a one-dimensional array of size J_ne and type ipc_, that holds the column indices of \(J\) 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.

J_ptr

is a one-dimensional array of size m+1 and type ipc_, that holds the starting position of each row of \(J\), 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.

H_type

is a one-dimensional array of type char that specifies the symmetric storage scheme used for the Hessian, \(H_L\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, ‘diagonal’ or ‘absent’, the latter if access to \(H\) is via matrix-vector products; 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_L\) in the sparse co-ordinate storage scheme. It need not be set for any of the other three 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_L\) 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_L\), 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 expo_reset_control(
    struct expo_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 expo_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, and the package is ready for the solve phase

void expo_solve_hessian_direct(
    void **data,
    void *userdata,
    ipc_ *status,
    ipc_ n,
    ipc_ m,
    ipc_ J_ne,
    ipc_ H_ne,
    const rpc_ c_l[],
    const rpc_ c_u[],
    const rpc_ x_l[],
    const rpc_ x_u[],
    rpc_ x[],
    rpc_ y[],
    rpc_ z[],
    rpc_ c[],
    rpc_ gl[],
    ipc_(*)(ipc_, ipc_, const rpc_[], rpc_[],
            const void*) eval_c,
    ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], rpc_[],
            const void*) eval_j,
    ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[],
            const void*) eval_h,
)

Find a local minimizer of a given constrained optimization problem.

This call is for the case where \(H(x,y) = \nabla_{xx}f(x) - \sum_i y_i \nabla_{xx}c_i(x)\) is provided specifically, and all function/derivative information is available by function calls.

Parameters:

data

holds private internal data

userdata

is a structure that allows data to be passed into the function and derivative evaluation programs.

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 restriction n > 0, m \(\geq\) 0 or requirement that type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, or ‘diagonal’ has been violated.

  • -9

    The analysis phase of the factorization failed; the return status from the factorization package is given in the component inform.factor_status

  • -10

    The factorization failed; the return status from the factorization package is given in the component inform.factor_status.

  • -11

    The solution of a set of linear equations using factors from the factorization package failed; the return status from the factorization package is given in the component inform.factor_status.

  • -16

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

  • -17

    The step is too small to make further impact.

  • -18

    Too many iterations have been performed. This may happen if control.max_it or control.max_eval is too small, but may also be symptomatic of a badly scaled problem.

  • -19

    The CPU time limit has been reached. This may happen if control.cpu_time_limit is too small, but may also be symptomatic of a badly scaled problem.

  • -82

    The user has forced termination of solver by removing the file named control.alive_file from unit unit control.alive_unit.

n

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

m

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

J_ne

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

H_ne

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

c_l

is a one-dimensional array of size m and type rpc_, that holds the values \(c_l\) of the lower bounds on the constraint functions \(c(x)\). The i-th component of c_l, \(i = 0, \ldots, m-1\), contains \(c_{li}\).

c_u

is a one-dimensional array of size m and type rpc_, that holds the values \(c_u\) of the upper bounds on the constraint functions \(c(x)\). The i-th component of c_u, \(i = 0, \ldots, m-1\), contains \(c_{ui}\).

x_l

is a one-dimensional array of size n and type rpc_, that holds the values \(x_l\) of the lower bounds on the optimization variables \(x\). The j-th component of x_l, \(j = 0, \ldots, n-1\), contains \(x_{lj}\).

x_u

is a one-dimensional array of size n and type rpc_, that holds the values \(x_u\) of the upper bounds on the optimization variables \(x\). The j-th component of x_u, \(j = 0, \ldots, n-1\), contains \(x_{uj}\).

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\). This should be set on input to an estimate of the minimizer.

y

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

z

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

c

is a one-dimensional array of size m and type rpc_, that holds the constraints \(c(x)\). The i-th component of c, i = 0, … , n-1, contains \(c_i(x)\).

gl

is a one-dimensional array of size n and type rpc_, that holds the gradient \(g_L(x,y)\) of the Lagrangian function. The j-th component of gl, j = 0, … , n-1, contains \(g_{Lj}\).

eval_fc

is a user-supplied function that must have the following signature:

ipc_ eval_fc( ipc_ n, const rpc_ x[], rpc_ f, rpc_ c[],
              const void *userdata )

The value of the objective function \(f(x)\) and the components of the constraint function \(c(x)\) evaluated at x= \(x\) must be assigned to f and c, respectively, and the function return value set to 0. If the evaluation is impossible at x, return should be set to a nonzero value. Data may be passed into eval_fc via the structure userdata.

eval_gj

is a user-supplied function that must have the following signature:

ipc_ eval_gj( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[],
              rpc_ j[], const void *userdata )

The components of the gradient \(g = g(x)\) of the objective and Jacobian \(J = \nabla_x c(x\)) of the constraints must be assigned to g and to j, in the same order as presented to expo_import, and the function return value set to 0. If the evaluation is impossible at x, return should be set to a nonzero value. Data may be passed into eval_gj via the structure userdata.

eval_hl

is a user-supplied function that must have the following signature:

ipc_ eval_hl( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[],
             const rpc_ y[], rpc_ h[], const void *userdata )

The nonzeros of the matrix \(H_L(x,y) = \nabla_{xx}f(x) -\sum_i y_i \nabla_{xx}c_i(x)\) of the Hessian of the Lagrangian function evaluated at x= \(x\) and y= \(y\) must be assigned to h in the same order as presented to expo_import, and the function return value set to 0. If the evaluation is impossible at x, return should be set to a nonzero value. Data may be passed into eval_hl via the structure userdata.

void expo_information(void **data, struct expo_inform_type* inform, ipc_ *status)

Provides output information

Parameters:

data

holds private internal data

inform

is a struct containing output information (see expo_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 expo_terminate(
    void **data,
    struct expo_control_type* control,
    struct expo_inform_type* inform
)

Deallocate all internal private storage

Parameters:

data

holds private internal data

control

is a struct containing control information (see expo_control_type)

inform

is a struct containing output information (see expo_inform_type)

available structures#

expo_control_type structure#

#include <galahad_expo.h>

struct expo_control_type {
    // components

    bool f_indexing;
    ipc_ error;
    ipc_ out;
    ipc_ print_level;
    ipc_ start_print;
    ipc_ stop_print;
    ipc_ print_gap;
    ipc_ max_it;
    ipc_ max_eval;
    ipc_ alive_unit;
    char alive_file[31];
    ipc_ update_multipliers_itmin;
    rpc_ update_multipliers_tol;
    rpc_ infinity;
    rpc_ stop_abs_p;
    rpc_ stop_rel_p;
    rpc_ stop_abs_d;
    rpc_ stop_rel_d;
    rpc_ stop_abs_c;
    rpc_ stop_rel_c;
    rpc_ stop_s;
    rpc_ initial_mu;
    rpc_ mu_reduce;
    rpc_ obj_unbounded;
    rpc_ try_advanced_start;
    rpc_ try_sqp_start;
    rpc_ stop_advanced_start;
    rpc_ cpu_time_limit;
    rpc_ clock_time_limit;
    bool hessian_available;
    bool subproblem_direct;
    bool space_critical;
    bool deallocate_error_fatal;
    char prefix[31];
    struct bsc_control_type bsc_control;
    struct tru_control_type tru_control;
    struct ssls_control_type ssls_control;
};

detailed documentation#

control derived type as a C struct

components#

bool f_indexing

use C or Fortran sparse matrix indexing

ipc_ error

error and warning diagnostics occur on stream error

ipc_ out

general output occurs on stream out

ipc_ print_level

the level of output required.

  • \(\leq\) 0 gives no output,

  • = 1 gives a one-line summary for every iteration,

  • = 2 gives a summary of the inner iteration for each iteration,

  • \(\geq\) 3 gives increasingly verbose (debugging) output

ipc_ start_print

any printing will start on this iteration

ipc_ stop_print

any printing will stop on this iteration

ipc_ print_gap

the number of iterations between printing

ipc_ max_it

the maximum number of iterations permitted

ipc_ max_eval

the maximum number of function evaluations permitted

ipc_ alive_unit

removal of the file alive_file from unit alive_unit terminates execution

char alive_file[31]

see alive_unit

ipc_ update_multipliers_itmin

update the Lagrange multipliers/dual variables from iteration .update_multipliers_itmin (<0 means never) and once the primal infeasibility is below .update_multipliers_tol

rpc_ update_multipliers_tol

see update_multipliers_itmin

rpc_ infinity

any bound larger than infinity in modulus will be regarded as infinite

rpc_ stop_abs_p

the required absolute and relative accuracies for the primal infeasibility

rpc_ stop_rel_p

see stop_abs_p

rpc_ stop_abs_d

the required absolute and relative accuracies for the dual infeasibility

rpc_ stop_rel_d

see stop_abs_d

rpc_ stop_abs_c

the required absolute and relative accuracies for complementary slackness

rpc_ stop_rel_c

see stop_abs_c

rpc_ stop_s

the smallest the norm of the step can be before termination

rpc_ initial_mu

initial value for the penalty parameter (<=0 means set automatically)

rpc_ mu_reduce

the amount by which the penalty parameter is decreased

rpc_ obj_unbounded

the smallest value the objective function may take before the problem is marked as unbounded

rpc_ try_advanced_start

try an advanced start at the end of every iteration when the KKT residuals are smaller than .try_advanced_start (-ve means never)

rpc_ try_sqp_start

try an advanced SQP start at the end of every iteration when the KKT residuals are smaller than .try_sqp_start (-ve means never)

rpc_ stop_advanced_start

stop the advanced start search once the residuals small tham .stop_advanced_start

rpc_ cpu_time_limit

the maximum CPU time allowed (-ve means infinite)

rpc_ clock_time_limit

the maximum elapsed clock time allowed (-ve means infinite)

bool hessian_available

is the Hessian matrix of second derivatives available or is access only via matrix-vector products (coming soon)?

bool subproblem_direct

use a direct (factorization) or (preconditioned) iterative method (coming soon) to find the search direction

bool space_critical

if .space_critical true, every effort will be made to use as little space as possible. This may result in longer computation time

bool deallocate_error_fatal

if .deallocate_error_fatal is true, any array/pointer deallocation error will terminate execution. Otherwise, computation will continue

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 bsc_control_type bsc_control

control parameters for BSC

struct tru_control_type tru_control

control parameters for TRU

struct ssls_control_type ssls_control

control parameters for SSLS

expo_time_type structure#

#include <galahad_expo.h>

struct expo_time_type {
    // components

    spc_ total;
    spc_ preprocess;
    spc_ analyse;
    spc_ factorize;
    spc_ solve;
    rpc_ clock_total;
    rpc_ clock_preprocess;
    rpc_ clock_analyse;
    rpc_ clock_factorize;
    rpc_ clock_solve;
};

detailed documentation#

time derived type as a C struct

components#

spc_ total

the total CPU time spent in the package

spc_ preprocess

the CPU time spent preprocessing the problem

spc_ analyse

the CPU time spent analysing the required matrices prior to factorization

spc_ factorize

the CPU time spent factorizing the required matrices

spc_ solve

the CPU time spent computing the search direction

rpc_ clock_total

the total clock time spent in the package

rpc_ clock_preprocess

the clock time spent preprocessing the problem

rpc_ clock_analyse

the clock time spent analysing the required matrices prior to factorization

rpc_ clock_factorize

the clock time spent factorizing the required matrices

rpc_ clock_solve

the clock time spent computing the search direction

expo_inform_type structure#

#include <galahad_expo.h>

struct expo_inform_type {
    // components

    ipc_ status;
    ipc_ alloc_status;
    char bad_alloc[81];
    char bad_eval[13];
    ipc_ iter;
    ipc_ fc_eval;
    ipc_ gj_eval;
    rpc_ obj;
    rpc_ primal_infeasibility;
    rpc_ dual_infeasibility;
    rpc_ complementary_slackness;
    struct expo_time_type time;
    struct tru_inform_type tru_inform;
    struct bsc_inform_type bsc_inform;
    struct ssls_inform_type ssls_inform;
};

detailed documentation#

inform derived type as a C struct

components#

ipc_ status

return status. See EXPO_solve for details

ipc_ alloc_status

the status of the last attempted allocation/deallocation

char bad_alloc[81]

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

char bad_eval[13]

the name of the user-supplied evaluation routine for which an error occurred

ipc_ iter

the total number of iterations performed

ipc_ fc_eval

the total number of evaluations of the objective f(x) and constraints c(x)

ipc_ gj_eval

the total number of evaluations of the gradient g(x) of f(x) and Jacobian J(x) of c(x)

ipc_ hl_eval

the total number of evaluations of the Hessian H(x,y) of the Lagrangian

rpc_ obj

the value of the objective function \(f(x)\) at the best estimate the solution, x, determined by EXPO_solve

rpc_ primal_infeasibility

the norm of the primal infeasibility (1) at the best estimate of the solution x, determined by EXPO_solve

rpc_ dual_infeasibility

the norm of the dual infeasibility (2) at the best estimate, (x,y,z), of the solution determined by EXPO_solve

rpc_ complementary_slackness

the norm of the complementary slackness (3) at the best estimate, (x,y,z), of the solution determined by EXPO_solve

struct expo_time_type time

timings (see above)

struct bsc_inform_type bsc_inform

inform parameters for BSC

struct tru_inform_type tru_inform

inform parameters for TRU

struct ssls_inform_type ssls_inform

inform parameters for SSLS

example calls#

This is an example of how to use the package to solve a nonlinearly constrained optimization problem; the code is available in $GALAHAD/src/expo/C/expot.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.

/* expot.c */
/* Full test for the EXPO C interface using C sparse matrix indexing */
/* Jari Fowkes & Nick Gould, STFC-Rutherford Appleton Laboratory, 2025 */

#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_expo.h"
#ifdef REAL_128
#include <quadmath.h>
#endif

// Custom userdata struct
struct userdata_type {
   rpc_ p;
};

// Function prototypes
ipc_ fc( ipc_ n, ipc_ m, const rpc_ x[], rpc_ *f, rpc_ c[], const void * );
ipc_ gj( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[],
         rpc_ jval[], const void * );
ipc_ hl( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
         rpc_ hval[], const void * );
ipc_ gj_dense( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[],
               rpc_ jval[],  const void * );
ipc_ hl_dense( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
               rpc_ hval[], const void * );

int main(void) {

    // Derived types
    void *data;
    struct expo_control_type control;
    struct expo_inform_type inform;

    // Set user data
    struct userdata_type userdata;
    userdata.p = 9.0;

    // Set problem data
    ipc_ n = 2; // # variables
    ipc_ m = 5; // # constraints
    ipc_ j_ne = 10; // Jacobian elements
    ipc_ h_ne = 2; // Hesssian elements
    ipc_ j_ne_dense = 10; // Jacobian elements
    ipc_ h_ne_dense = 3; // dense Hesssian elements
    // 0-based indices
    ipc_ J_row[] = {0, 0, 1, 1, 2, 2, 3, 3, 4, 4}; // Jacobian J
    ipc_ J_col[] = {0, 1, 0, 1, 0, 1, 0, 1, 0, 1}; //
    ipc_ J_ptr[] = {0, 2, 4, 6, 8, 10 };           // row pointers
    ipc_ H_row[] = {0, 1};     // Hessian H
    ipc_ H_col[] = {0, 1};     // NB lower triangle
    ipc_ H_ptr[] = {0, 1, 2};  // row pointers

    // Set storage
    rpc_ y[m]; // multipliers
    rpc_ z[n]; // dual variables
    rpc_ c[m]; // constraints
    rpc_ gl[n]; // gradients
    rpc_ xl[] = {-50.0, -50.0}; // lower variable bounds
    rpc_ xu[] = {50.0, 50.0}; // upper variable bounds
    rpc_ cl[] = {0.0, 0.0, 0.0, 0.0, 0.0}; // lower constraint bounds
    rpc_ cu[] = {INFINITY, INFINITY, INFINITY,
                 INFINITY, INFINITY}; // upper constraint bounds
    char st = ' ';
    ipc_ status;

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

    printf(" tests options for all-in-one storage format\n\n");

//  for( ipc_ d=1; d <= 3; d++){
    for( ipc_ d=1; d <= 4; d++){

        // Initialize EXPO
        expo_initialize( &data, &control, &inform );

        // Set user-defined control options
        control.f_indexing = false; // C sparse matrix indexing
        //control.print_level = 1;
        control.max_it = 20;
        control.max_eval = 100;
        control.stop_abs_p = 0.00001;
        control.stop_abs_d = 0.00001;
        control.stop_abs_c = 0.00001;

        rpc_ x[] = {3.0,1.0}; // starting point

        switch(d){
            case 1: // sparse co-ordinate storage
                st = 'C';
                expo_import( &control, &data, &status, n, m,
                             "coordinate", j_ne, J_row, J_col, NULL,
                             "coordinate", h_ne, H_row, H_col, NULL );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne, h_ne,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj, hl );
                break;
            case 2: // sparse by rows
                st = 'R';
                expo_import( &control, &data, &status, n, m,
                             "sparse_by_rows", j_ne, NULL, J_col, J_ptr,
                             "sparse_by_rows", h_ne, NULL, H_col, H_ptr );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne, h_ne,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj, hl );
                break;
            case 3: // dense
                st = 'D';
                expo_import( &control, &data, &status, n, m,
                             "dense", j_ne_dense, NULL, NULL, NULL,
                             "dense", h_ne_dense, NULL, NULL, NULL );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne_dense, h_ne_dense,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj_dense, hl_dense );
                break;
            case 4: // diagonal
                st = 'I';
                expo_import( &control, &data, &status, n, m,
                             "sparse_by_rows", j_ne, NULL, J_col, J_ptr,
                             "diagonal", n, NULL, NULL, NULL );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne, n,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj, hl );
                break;
        }

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

        if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
            printf("%c:%6" i_ipc_ " iterations. Optimal objective "
                   "value = %.2f status = %1" i_ipc_ "\n",
                   st, inform.iter, inform.obj, inform.status);
#endif
        }else{
            printf("%c: EXPO_solve exit status = %1" i_ipc_ "\n",
                   st, inform.status);
        }
        // Delete internal workspace
        expo_terminate( &data, &control, &inform );
    }
}

// compute the function and constraint values
ipc_ fc( ipc_ n, ipc_ m, const rpc_ x[], rpc_ *f, rpc_ c[],
         const void *userdata ){
    struct userdata_type *myuserdata = ( struct userdata_type * ) userdata;
    rpc_ p = myuserdata->p;
    *f = pow(x[0],2.0) + pow(x[1],2.0);
    c[0] = x[0] + x[1] - 1.0;
    c[1] = pow(x[0],2.0) + pow(x[1],2.0) - 1.0;
    c[2] = p * pow(x[0],2.0) + pow(x[1],2.0) - p;
    c[3] = pow(x[0],2.0) - x[1];
    c[4] = pow(x[1],2.0) - x[0];
    return 0;
}

// compute the gradient and Jacobian
ipc_ gj( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[], rpc_ jval[],
         const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    g[0] = 2.0 * x[0];
    g[1] = 2.0 * x[1];
    jval[0] = 1.0;
    jval[1] = 1.0;
    jval[2] = 2.0 * x[0];
    jval[3] = 2.0 * x[1];
    jval[4] = 2.0 * p * x[0];
    jval[5] = 2.0 * x[1];
    jval[6] = 2.0 * x[0];
    jval[7] = - 1.0;
    jval[8] = - 1.0;
    jval[9] = 2.0 * x[1];
    return 0;
}

// compute the Hessian of the Lagrangian
ipc_ hl( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
         rpc_ hval[], const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    hval[0] = 2.0 - 2.0 * (y[1] + p * y[2] + y[3]);
    hval[1] = 2.0 - 2.0 * (y[1] + y[2] + y[4]);
    return 0;
}

// compute the gradient and dense Jacobian
ipc_ gj_dense( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[],
               rpc_ jval[], const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    g[0] = 2.0 * x[0];
    g[1] = 2.0 * x[1];
    jval[0] = 1.0;
    jval[1] = 1.0;
    jval[2] = 2.0 * x[0];
    jval[3] = 2.0 * x[1];
    jval[4] = 2.0 * p * x[0];
    jval[5] = 2.0 * x[1];
    jval[6] = 2.0 * x[0];
    jval[7] = - 1.0;
    jval[8] = - 1.0;
    jval[9] = 2.0 * x[1];
    return 0;
}

// compute the dense Hessian of the Lagrangian
ipc_ hl_dense( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
               rpc_ hval[], const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    hval[0] = 2.0 - 2.0 * (y[1] + p * y[2] +  y[3]);
    hval[1] = 0.0;
    hval[2] = 2.0 - 2.0 * (y[1] + y[2] + y[4]);
    return 0;
}

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

/* expot.c */
/* Full test for the EXPO C interface using Fortran sparse matrix indexing */
/* Jari Fowkes & Nick Gould, STFC-Rutherford Appleton Laboratory, 2025 */

#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_expo.h"
#ifdef REAL_128
#include <quadmath.h>
#endif

// Custom userdata struct
struct userdata_type {
   rpc_ p;
};

// Function prototypes
ipc_ fc( ipc_ n, ipc_ m, const rpc_ x[], rpc_ *f, rpc_ c[], const void * );
ipc_ gj( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[],
         rpc_ jval[], const void * );
ipc_ hl( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
         rpc_ hval[], const void * );
ipc_ gj_dense( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[],
               rpc_ jval[],  const void * );
ipc_ hl_dense( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
               rpc_ hval[], const void * );

int main(void) {

    // Derived types
    void *data;
    struct expo_control_type control;
    struct expo_inform_type inform;

    // Set user data
    struct userdata_type userdata;
    userdata.p = 9.0;

    // Set problem data
    ipc_ n = 2; // # variables
    ipc_ m = 5; // # constraints
    ipc_ j_ne = 10; // Jacobian elements
    ipc_ h_ne = 2; // Hesssian elements
    ipc_ j_ne_dense = 10; // dense Jacobian elements
    ipc_ h_ne_dense = 3; // dense Hesssian elements
    // 1-based indices
    ipc_ J_row[] = {1, 1, 2, 2, 3, 3, 4, 4, 5, 5}; // Jacobian J
    ipc_ J_col[] = {1, 2, 1, 2, 1, 2, 1, 2, 1, 2}; //
    ipc_ J_ptr[] = {1, 3, 5, 7, 9, 11 };           // row pointers
    ipc_ H_row[] = {1, 2};     // Hessian H
    ipc_ H_col[] = {1, 2};     // NB lower triangle
    ipc_ H_ptr[] = {1, 2, 3};  // row pointers

    // Set storage
    rpc_ y[m]; // multipliers
    rpc_ z[n]; // dual variables
    rpc_ c[m]; // constraints
    rpc_ gl[n]; // gradients
    rpc_ xl[] = {-50.0, -50.0}; // lower variable bounds
    rpc_ xu[] = {50.0, 50.0}; // upper variable bounds
    rpc_ cl[] = {0.0, 0.0, 0.0, 0.0, 0.0}; // lower constraint bounds
    rpc_ cu[] = {INFINITY, INFINITY, INFINITY,
                 INFINITY, INFINITY}; // upper constraint bounds
    char st = ' ';
    ipc_ status;

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

    printf(" tests options for all-in-one storage format\n\n");

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

        // Initialize EXPO
        expo_initialize( &data, &control, &inform );

        // Set user-defined control options
        control.f_indexing = true; // Fortran sparse matrix indexing
        //control.print_level = 1;
        control.max_it = 20;
        control.max_eval = 100;
        control.stop_abs_p = 0.00001;
        control.stop_abs_d = 0.00001;
        control.stop_abs_c = 0.00001;

        rpc_ x[] = {3.0,1.0}; // starting point

        switch(d){
            case 1: // sparse co-ordinate storage
                st = 'C';
                expo_import( &control, &data, &status, n, m,
                             "coordinate", j_ne, J_row, J_col, NULL,
                             "coordinate", h_ne, H_row, H_col, NULL );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne, h_ne,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj, hl );
                break;
            case 2: // sparse by rows
                st = 'R';
                expo_import( &control, &data, &status, n, m,
                             "sparse_by_rows", j_ne, NULL, J_col, J_ptr,
                             "sparse_by_rows", h_ne, NULL, H_col, H_ptr );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne, h_ne,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj, hl );
                break;
            case 3: // dense
                st = 'D';
                expo_import( &control, &data, &status, n, m,
                             "dense", j_ne_dense, NULL, NULL, NULL,
                             "dense", h_ne_dense, NULL, NULL, NULL );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne_dense, h_ne_dense,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj_dense, hl_dense );
                break;
            case 4: // diagonal
                st = 'I';
                expo_import( &control, &data, &status, n, m,
                             "sparse_by_rows", j_ne, NULL, J_col, J_ptr,
                             "diagonal", n, NULL, NULL, NULL );
                expo_solve_hessian_direct( &data, &userdata, &status,
                                           n, m, j_ne, n,
                                           cl, cu, xl, xu, x, y, z, c, gl,
                                           fc, gj, hl );
                break;
        }

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

        if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
            printf("%c:%6" i_ipc_ " iterations. Optimal objective "
                   "value = %.2f status = %1" i_ipc_ "\n",
                   st, inform.iter, inform.obj, inform.status);
#endif
        }else{
            printf("%c: EXPO_solve exit status = %1" i_ipc_ "\n",
                   st, inform.status);
        }
        // Delete internal workspace
        expo_terminate( &data, &control, &inform );
    }
}

// compute the function and constraint values
ipc_ fc( ipc_ n, ipc_ m, const rpc_ x[], rpc_ *f, rpc_ c[],
         const void *userdata ){
    struct userdata_type *myuserdata = ( struct userdata_type * ) userdata;
    rpc_ p = myuserdata->p;
    *f = pow(x[0],2.0) + pow(x[1],2.0);
    c[0] = x[0] + x[1] - 1.0;
    c[1] = pow(x[0],2.0) + pow(x[1],2.0) - 1.0;
    c[2] = p * pow(x[0],2.0) + pow(x[1],2.0) - p;
    c[3] = pow(x[0],2.0) - x[1];
    c[4] = pow(x[1],2.0) - x[0];
    return 0;
}

// compute the gradient and Jacobian
ipc_ gj( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[], rpc_ jval[],
         const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    g[0] = 2.0 * x[0];
    g[1] = 2.0 * x[1];
    jval[0] = 1.0;
    jval[1] = 1.0;
    jval[2] = 2.0 * x[0];
    jval[3] = 2.0 * x[1];
    jval[4] = 2.0 * p * x[0];
    jval[5] = 2.0 * x[1];
    jval[6] = 2.0 * x[0];
    jval[7] = - 1.0;
    jval[8] = - 1.0;
    jval[9] = 2.0 * x[1];
    return 0;
}

// compute the Hessian of the Lagrangian
ipc_ hl( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
         rpc_ hval[], const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    hval[0] = 2.0 - 2.0 * (y[1] + p * y[2] + y[3]);
    hval[1] = 2.0 - 2.0 * (y[1] + y[2] + y[4]);
    return 0;
}

// compute the gradient and dense Jacobian
ipc_ gj_dense( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ g[],
               rpc_ jval[], const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    g[0] = 2.0 * x[0];
    g[1] = 2.0 * x[1];
    jval[0] = 1.0;
    jval[1] = 1.0;
    jval[2] = 2.0 * x[0];
    jval[3] = 2.0 * x[1];
    jval[4] = 2.0 * p * x[0];
    jval[5] = 2.0 * x[1];
    jval[6] = 2.0 * x[0];
    jval[7] = - 1.0;
    jval[8] = - 1.0;
    jval[9] = 2.0 * x[1];
    return 0;
}

// compute the dense Hessian of the Lagrangian
ipc_ hl_dense( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[],
               rpc_ hval[], const void *userdata ){
    struct userdata_type *myuserdata = (struct userdata_type *) userdata;
    rpc_ p = myuserdata->p;
    hval[0] = 2.0 - 2.0 * (y[1] + p * y[2] +  y[3]);
    hval[1] = 0.0;
    hval[2] = 2.0 - 2.0 * (y[1] + y[2] + y[4]);
    return 0;
}