GALAHAD ARC package#

purpose#

The arc package uses a regularization method to find a (local) minimizer of a differentiable objective function \(f(x)\) of many variables \(x\). The method offers the choice of direct and iterative solution of the key subproblems, and is most suitable for large problems. First derivatives are required, and if second derivatives can be calculated, they will be exploited.

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

method#

An adaptive cubic regularization method is used. In this, an improvement to a current estimate of the required minimizer, \(x_k\) is sought by computing a step \(s_k\). The step is chosen to approximately minimize a model \(m_k(s)\) of \(f(x_k + s)\) that includes a weighted term \(\sigma_k \|s_k\|^3\) for some specified positive weight \(\sigma_k\). The quality of the resulting step \(s_k\) is assessed by computing the “ratio” \((f(x_k) - f(x_k + s_k))/ (m_k(0) - m_k(s_k))\). The step is deemed to have succeeded if the ratio exceeds a given \(\eta_s > 0\), and in this case \(x_{k+1} = x_k + s_k\). Otherwise \(x_{k+1} = x_k\), and the weight is increased by powers of a given increase factor up to a given limit. If the ratio is larger than \(\eta_v \geq \eta_d\), the weight will be decreased by powers of a given decrease factor again up to a given limit. The method will terminate as soon as \(\|\nabla_x f(x_k)\|\) is smaller than a specified value.

Either linear or quadratic models \(m_k(s)\) may be used. The former will be taken as the first two terms \(f(x_k) + s^T \nabla_x f(x_k)\) of a Taylor series about \(x_k\), while the latter uses an approximation to the first three terms \(f(x_k) + s^T \nabla_x f(x_k) + \frac{1}{2} s^T B_k s\), for which \(B_k\) is a symmetric approximation to the Hessian \(\nabla_{xx} f(x_k)\); possible approximations include the true Hessian, limited-memory secant and sparsity approximations and a scaled identity matrix. Normally a two-norm regularization will be used, but this may change if preconditioning is employed.

An approximate minimizer of the cubic model is found using either a direct approach involving factorization or an iterative (conjugate-gradient/Lanczos) approach based on approximations to the required solution from a so-called Krlov subspace. The direct approach is based on the knowledge that the required solution satisfies the linear system of equations \((B_k + \lambda_k I) s_k = - \nabla_x f(x_k)\) involving a scalar Lagrange multiplier \(\lambda_k\). This multiplier is found by uni-variate root finding, using a safeguarded Newton-like process, by RQS or DPS (depending on the norm chosen). The iterative approach uses GLRT, and is best accelerated by preconditioning with good approximations to \(B_k\) using PSLS. The iterative approach has the advantage that only matrix-vector products \(B_k v\) are required, and thus \(B_k\) is not required explicitly. However when factorizations of \(B_k\) are possible, the direct approach is often more efficient.

references#

The generic adaptive cubic regularization method is described in detail in

C. Cartis, N. I. M. Gould and Ph. L. Toint, ``Adaptive cubic regularisation methods for unconstrained optimization. Part I: motivation, convergence and numerical results’’ Mathematical Programming 127(2) (2011) 245–295,

and uses ``tricks’’ as suggested in

N. I. M. Gould, M. Porcelli and Ph. L. Toint, ``Updating the regularization parameter in the adaptive cubic regularization algorithm’’. Computational Optimization and Applications 53(1) (2012) 1-22.

matrix storage#

The symmetric \(n\) by \(n\) matrix \(H = \nabla^2_{xx}f\) 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 \(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.

introduction to function calls#

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

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

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

  • arc_import - set up problem data structures and fixed values

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

  • solve the problem by calling one of

    • arc_solve_with_mat - solve using function calls to evaluate function, gradient and Hessian values

    • arc_solve_without_mat - solve using function calls to evaluate function and gradient values and Hessian-vector products

    • arc_solve_reverse_with_mat - solve returning to the calling program to obtain function, gradient and Hessian values, or

    • arc_solve_reverse_without_mat - solve returning to the calling prorgram to obtain function and gradient values and Hessian-vector products

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

  • arc_terminate - deallocate data structures

See the examples section for illustrations of use.

parametric real type T#

Below, the symbol T refers to a parametric real type that may be Float32 (single precision) or Float64 (double precision).

callable functions#

    function arc_initialize(T, 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 arc_control_type)

status

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

  • 0

    The initialization was successful.

    function arc_read_specfile(T, 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/arc/ARC.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/arc.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 arc_control_type)

specfile

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

    function arc_import(T, control, data, status, n, H_type, 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 arc_control_type)

data

holds private internal data

status

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

n

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

H_type

is a one-dimensional array of type Vararg{Cchar} that specifies the symmetric storage scheme used for the Hessian. It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, ‘diagonal’ or ‘absent’, the latter if access to the Hessian is via matrix-vector products; lower or upper case variants are allowed

ne

is a scalar variable of type Int32 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 three schemes.

H_row

is a one-dimensional array of size ne and type Int32 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 ne and type Int32 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 Int32 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 arc_reset_control(T, 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 arc_control_type)

data

holds private internal data

status

is a scalar variable of type Int32 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

    function arc_solve_with_mat(T, data, userdata, status, n, x, g, ne,
                                eval_f, eval_g, eval_h, eval_prec)

Find a local minimizer of a given function using a regularization method.

This call is for the case where \(H = \nabla_{xx}f(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 Int32 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 or requirement that type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, ‘diagonal’ or ‘absent’ has been violated.

  • -7

    The objective function appears to be unbounded from below

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

  • -18

    Too many iterations have been performed. This may happen if control.maxit 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 Int32 that holds the number of variables

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

g

is a one-dimensional array of size n and type T that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 1, … , n, contains \(g_j\).

ne

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

eval_f

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

function eval_f(n, x, f, userdata)

The value of the objective function \(f(x)\) evaluated at x=\(x\) must be assigned to f, 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_f via the structure userdata.

eval_g

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

function eval_g(n, x, g, userdata)

The components of the gradient \(g = \nabla_x f(x\)) of the objective function evaluated at x=\(x\) must be assigned to g, 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_g via the structure userdata.

eval_h

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

function eval_h(n, ne, x, h, userdata)

The nonzeros of the Hessian \(H = \nabla_{xx}f(x)\) of the objective function evaluated at x=\(x\) must be assigned to h in the same order as presented to arc_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_h via the structure userdata.

eval_prec

is an optional user-supplied function that may be C_NULL. If non-NULL, it must have the following signature:

function eval_prec(n, x, u, v, userdata)

The product \(u = P(x) v\) of the user’s preconditioner \(P(x)\) evaluated at \(x\) with the vector v=\(v\), the result \(u\) must be retured in u, 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_prec via the structure userdata.

    function arc_solve_without_mat(T, data, userdata, status, n, x, g,
                                   eval_f, eval_g, eval_hprod, eval_prec)

Find a local minimizer of a given function using a regularization method.

This call is for the case where access to \(H = \nabla_{xx}f(x)\) is provided by Hessian-vector products, 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 Int32 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 or requirement that type contains its relevant string ‘dense’, ‘coordinate’, ‘sparse_by_rows’, ‘diagonal’ or ‘absent’ has been violated.

  • -7

    The objective function appears to be unbounded from below

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

  • -18

    Too many iterations have been performed. This may happen if control.maxit 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 Int32 that holds the number of variables

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

g

is a one-dimensional array of size n and type T that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 1, … , n, contains \(g_j\).

eval_f

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

function eval_f(n, x, f, userdata)

The value of the objective function \(f(x)\) evaluated at x=\(x\) must be assigned to f, 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_f via the structure userdata.

eval_g

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

function eval_g(n, x, g, userdata)

The components of the gradient \(g = \nabla_x f(x\)) of the objective function evaluated at x=\(x\) must be assigned to g, 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_g via the structure userdata.

eval_hprod

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

function eval_hprod(n, x, u, v, got_h, userdata)

The sum \(u + \nabla_{xx}f(x) v\) of the product of the Hessian \(\nabla_{xx}f(x)\) of the objective function evaluated at x=\(x\) with the vector v=\(v\) and the vector $ \(u\) must be returned in u, and the function return value set to 0. If the evaluation is impossible at x, return should be set to a nonzero value. The Hessian has already been evaluated or used at x if got_h is true. Data may be passed into eval_hprod via the structure userdata.

eval_prec

is an optional user-supplied function that may be C_NULL. If non-NULL, it must have the following signature:

function eval_prec(n, x, u, v, userdata)

The product \(u = P(x) v\) of the user’s preconditioner \(P(x)\) evaluated at \(x\) with the vector v=\(v\), the result \(u\) must be retured in u, 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_prec via the structure userdata.

    function arc_solve_reverse_with_mat(T, data, status, eval_status,
                                        n, x, f, g, ne, H_val, u, v)

Find a local minimizer of a given function using a regularization method.

This call is for the case where \(H = \nabla_{xx}f(x)\) is provided specifically, but function/derivative information is only available by returning to the calling procedure

Parameters:

data

holds private internal data

status

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

  • -7

    The objective function appears to be unbounded from below.

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

  • -18

    Too many iterations have been performed. This may happen if control.maxit 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.

  • 2

    The user should compute the objective function value \(f(x)\) at the point \(x\) indicated in x and then re-enter the function. The required value should be set in f, and eval_status should be set to 0. If the user is unable to evaluate \(f(x)\) for instance, if the function is undefined at \(x\) the user need not set f, but should then set eval_status to a non-zero value.

  • 3

    The user should compute the gradient of the objective function \(\nabla_x f(x)\) at the point \(x\) indicated in x and then re-enter the function. The value of the i-th component of the g radient should be set in g[i], for i = 1, …, n and eval_status should be set to 0. If the user is unable to evaluate a component of \(\nabla_x f(x)\) for instance if a component of the gradient is undefined at \(x\) -the user need not set g, but should then set eval_status to a non-zero value.

  • 4

    The user should compute the Hessian of the objective function \(\nabla_{xx} f(x)\) at the point x indicated in \(x\) and then re-enter the function. The value l-th component of the Hessian stored according to the scheme input in the remainder of \(H\) should be set in H_val[l], for l = 0, …, ne-1 and eval_status should be set to 0. If the user is unable to evaluate a component of \(\nabla_{xx}f(x)\) for instance, if a component of the Hessian is undefined at \(x\) the user need not set H_val, but should then set eval_status to a non-zero value.

  • 6

    The user should compute the product \(u = P(x)v\) of their preconditioner \(P(x)\) at the point x indicated in \(x\) with the vector \(v\) and then re-enter the function. The vector \(v\) is given in v, the resulting vector \(u = P(x)v\) should be set in u and eval_status should be set to 0. If the user is unable to evaluate the product for instance, if a component of the preconditioner is undefined at \(x\) the user need not set u, but should then set eval_status to a non-zero value.

eval_status

is a scalar variable of type Int32 that is used to indicate if objective function/gradient/Hessian values can be provided (see above)

n

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

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

f

is a scalar variable pointer of type T that holds the value of the objective function.

g

is a one-dimensional array of size n and type T that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 1, … , n, contains \(g_j\).

ne

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

u

is a one-dimensional array of size n and type T that is used for reverse communication (see above for details)

v

is a one-dimensional array of size n and type T that is used for reverse communication (see above for details)

    function arc_solve_reverse_without_mat(T, data, status, eval_status,
                                           n, x, f, g, u, v)

Find a local minimizer of a given function using a regularization method.

This call is for the case where access to \(H = \nabla_{xx}f(x)\) is provided by Hessian-vector products, but function/derivative information is only available by returning to the calling procedure.

Parameters:

data

holds private internal data

status

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

  • -7

    The objective function appears to be unbounded from below

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

  • -18

    Too many iterations have been performed. This may happen if control.maxit 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.

  • 2

    The user should compute the objective function value \(f(x)\) at the point \(x\) indicated in x and then re-enter the function. The required value should be set in f, and eval_status should be set to 0. If the user is unable to evaluate \(f(x)\) for instance, if the function is undefined at \(x\) the user need not set f, but should then set eval_status to a non-zero value.

  • 3

    The user should compute the gradient of the objective function \(\nabla_x f(x)\) at the point \(x\) indicated in x and then re-enter the function. The value of the i-th component of the g radient should be set in g[i], for i = 1, …, n and eval_status should be set to 0. If the user is unable to evaluate a component of \(\nabla_x f(x)\) for instance if a component of the gradient is undefined at \(x\) -the user need not set g, but should then set eval_status to a non-zero value.

  • 5

    The user should compute the product \(\nabla_{xx} f(x)v\) of the Hessian of the objective function \(\nabla_{xx} f(x)\) at the point \(x\) indicated in x with the vector \(v\), add the result to the vector \(u\) and then re-enter the function. The vectors \(u\) and \(v\) are given in u and v respectively, the resulting vector \(u + \nabla_{xx} f(x)v\) should be set in u and eval_status should be set to 0. If the user is unable to evaluate the product for instance, if a component of the Hessian is undefined at \(x\) the user need not alter u, but should then set eval_status to a non-zero value.

  • 6

    The user should compute the product \(u = P(x)v\) of their preconditioner \(P(x)\) at the point x indicated in \(x\) with the vector \(v\) and then re-enter the function. The vector \(v\) is given in v, the resulting vector \(u = P(x)v\) should be set in u and eval_status should be set to 0. If the user is unable to evaluate the product for instance, if a component of the preconditioner is undefined at \(x\) the user need not set u, but should then set eval_status to a non-zero value.

eval_status

is a scalar variable of type Int32 that is used to indicate if objective function/gradient/Hessian values can be provided (see above)

n

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

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

f

is a scalar variable pointer of type T that holds the value of the objective function.

g

is a one-dimensional array of size n and type T that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 1, … , n, contains \(g_j\).

u

is a one-dimensional array of size n and type T that is used for reverse communication (see above for details)

v

is a one-dimensional array of size n and type T that is used for reverse communication (see above for details)

    function arc_information(T, data, inform, status)

Provides output information

Parameters:

data

holds private internal data

inform

is a structure containing output information (see arc_inform_type)

status

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

  • 0

    The values were recorded successfully

    function arc_terminate(T, data, control, inform)

Deallocate all internal private storage

Parameters:

data

holds private internal data

control

is a structure containing control information (see arc_control_type)

inform

is a structure containing output information (see arc_inform_type)

available structures#

arc_control_type structure#

    struct arc_control_type{T}
      f_indexing::Bool
      error::Int32
      out::Int32
      print_level::Int32
      start_print::Int32
      stop_print::Int32
      print_gap::Int32
      maxit::Int32
      alive_unit::Int32
      alive_file::NTuple{31,Cchar}
      non_monotone::Int32
      model::Int32
      norm::Int32
      semi_bandwidth::Int32
      lbfgs_vectors::Int32
      max_dxg::Int32
      icfs_vectors::Int32
      mi28_lsize::Int32
      mi28_rsize::Int32
      advanced_start::Int32
      stop_g_absolute::T
      stop_g_relative::T
      stop_s::T
      initial_weight::T
      minimum_weight::T
      reduce_gap::T
      tiny_gap::T
      large_root::T
      eta_successful::T
      eta_very_successful::T
      eta_too_successful::T
      weight_decrease_min::T
      weight_decrease::T
      weight_increase::T
      weight_increase_max::T
      obj_unbounded::T
      cpu_time_limit::T
      clock_time_limit::T
      hessian_available::Bool
      subproblem_direct::Bool
      renormalize_weight::Bool
      quadratic_ratio_test::Bool
      space_critical::Bool
      deallocate_error_fatal::Bool
      prefix::NTuple{31,Cchar}
      rqs_control::rqs_control_type{T}
      glrt_control::glrt_control_type{T}
      dps_control::dps_control_type{T}
      psls_control::psls_control_type{T}
      lms_control::lms_control_type{T}
      lms_control_prec::lms_control_type{T}
      sha_control::sha_control_type

detailed documentation#

control derived type as a Julia structure

components#

Bool f_indexing

use C or Fortran sparse matrix indexing

Int32 error

error and warning diagnostics occur on stream error

Int32 out

general output occurs on stream out

Int32 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

Int32 start_print

any printing will start on this iteration

Int32 stop_print

any printing will stop on this iteration

Int32 print_gap

the number of iterations between printing

Int32 maxit

the maximum number of iterations performed

Int32 alive_unit

removal of the file alive_file from unit alive_unit terminates execution

char alive_file[31]

see alive_unit

Int32 non_monotone

the descent strategy used.

Possible values are

  • <= 0 a monotone strategy is used.

  • anything else, a non-monotone strategy with history length .non_monotine is used.

Int32 model

the model used.

Possible values are

  • 0 dynamic (not yet implemented)

  • 1 first-order (no Hessian)

  • 2 second-order (exact Hessian)

  • 3 barely second-order (identity Hessian)

  • 4 secant second-order (limited-memory BFGS, with .lbfgs_vectors history) (not yet implemented)

  • 5 secant second-order (limited-memory SR1, with .lbfgs_vectors history) (not yet implemented)

Int32 norm

the regularization norm used.

The norm is defined via \(\|v\|^2 = v^T P v\), and will define the preconditioner used for iterative methods. Possible values for \(P\) are

  • -3 users own preconditioner

  • -2 \(P =\) limited-memory BFGS matrix (with .lbfgs_vectors history)

  • -1 identity (= Euclidan two-norm)

  • 0 automatic (not yet implemented)

  • 1 diagonal, \(P =\) diag( max( Hessian, .min_diagonal ) )

  • 2 banded, \(P =\) band( Hessian ) with semi-bandwidth .semi_bandwidth

  • 3 re-ordered band, P=band(order(A)) with semi-bandwidth .semi_bandwidth

  • 4 full factorization, \(P =\) Hessian, Schnabel-Eskow modification

  • 5 full factorization, \(P =\) Hessian, GMPS modification (not yet implemented)

  • 6 incomplete factorization of Hessian, Lin-More’

  • 7 incomplete factorization of Hessian, HSL_MI28

  • 8 incomplete factorization of Hessian, Munskgaard (not yet implemented)

  • 9 expanding band of Hessian (not yet implemented)

  • 10 diagonalizing norm from GALAHAD_DPS (subproblem_direct only)

Int32 semi_bandwidth

specify the semi-bandwidth of the band matrix P if required

Int32 lbfgs_vectors

number of vectors used by the L-BFGS matrix P if required

Int32 max_dxg

number of vectors used by the sparsity-based secant Hessian if required

Int32 icfs_vectors

number of vectors used by the Lin-More’ incomplete factorization matrix P if required

Int32 mi28_lsize

the maximum number of fill entries within each column of the incomplete factor L computed by HSL_MI28. In general, increasing .mi28_lsize improve the quality of the preconditioner but increases the time to compute and then apply the preconditioner. Values less than 0 are treated as 0

Int32 mi28_rsize

the maximum number of entries within each column of the strictly lower triangular matrix \(R\) used in the computation of the preconditioner by HSL_MI28. Rank-1 arrays of size .mi28_rsize \* n are allocated internally to hold \(R\). Thus the amount of memory used, as well as the amount of work involved in computing the preconditioner, depends on .mi28_rsize. Setting .mi28_rsize > 0 generally leads to a higher quality preconditioner than using .mi28_rsize = 0, and choosing .mi28_rsize >= .mi28_lsize is generally recommended

Int32 advanced_start

try to pick a good initial regularization weight using .advanced_start iterates of a variant on the strategy of Sartenaer SISC 18(6) 1990:1788-1803

T stop_g_absolute

overall convergence tolerances. The iteration will terminate when the norm of the gradient of the objective function is smaller than MAX( .stop_g_absolute, .stop_g_relative * norm of the initial gradient ) or if the step is less than .stop_s

T stop_g_relative

see stop_g_absolute

T stop_s

see stop_g_absolute

T initial_weight

Initial value for the regularisation weight (-ve => 1/||g_0||)

T minimum_weight

minimum permitted regularisation weight

T reduce_gap

expert parameters as suggested in Gould, Porcelli & Toint, “Updating the regularization parameter in the adaptive cubic regularization algorithm” RAL-TR-2011-007, Rutherford Appleton Laboratory, England (2011), http://epubs.stfc.ac.uk/bitstream/6181/RAL-TR-2011-007.pdf (these are denoted beta, epsilon_chi and alpha_max in the paper)

T tiny_gap

see reduce_gap

T large_root

see reduce_gap

T eta_successful

a potential iterate will only be accepted if the actual decrease f - f(x_new) is larger than .eta_successful times that predicted by a quadratic model of the decrease. The regularization weight will be decreased if this relative decrease is greater than .eta_very_successful but smaller than .eta_too_successful (the first is eta in Gould, Porcell and Toint, 2011)

T eta_very_successful

see eta_successful

T eta_too_successful

see eta_successful

T weight_decrease_min

on very successful iterations, the regularization weight will be reduced by the factor .weight_decrease but no more than .weight_decrease_min while if the iteration is unsuccessful, the weight will be increased by a factor .weight_increase but no more than .weight_increase_max (these are delta_1, delta_2, delta3 and delta_max in Gould, Porcelli and Toint, 2011)

T weight_decrease

see weight_decrease_min

T weight_increase

see weight_decrease_min

T weight_increase_max

see weight_decrease_min

T obj_unbounded

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

T cpu_time_limit

the maximum CPU time allowed (-ve means infinite)

T 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?

Bool subproblem_direct

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

Bool renormalize_weight

should the weight be renormalized to account for a change in preconditioner?

Bool quadratic_ratio_test

should the test for acceptance involve the quadratic model or the cubic?

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

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 rqs_control_type rqs_control

control parameters for RQS

struct glrt_control_type glrt_control

control parameters for GLRT

struct dps_control_type dps_control

control parameters for DPS

struct psls_control_type psls_control

control parameters for PSLS

struct lms_control_type lms_control

control parameters for LMS

struct lms_control_type lms_control_prec

control parameters for LMS used for preconditioning

struct sha_control_type sha_control

control parameters for SHA

arc_time_type structure#

    struct arc_time_type{T}
      total::Float32
      preprocess::Float32
      analyse::Float32
      factorize::Float32
      solve::Float32
      clock_total::T
      clock_preprocess::T
      clock_analyse::T
      clock_factorize::T
      clock_solve::T

detailed documentation#

time derived type as a Julia structure

components#

Float32 total

the total CPU time spent in the package

Float32 preprocess

the CPU time spent preprocessing the problem

Float32 analyse

the CPU time spent analysing the required matrices prior to factorizatio

Float32 factorize

the CPU time spent factorizing the required matrices

Float32 solve

the CPU time spent computing the search direction

T clock_total

the total clock time spent in the package

T clock_preprocess

the clock time spent preprocessing the problem

T clock_analyse

the clock time spent analysing the required matrices prior to factorizat

T clock_factorize

the clock time spent factorizing the required matrices

T clock_solve

the clock time spent computing the search direction

arc_inform_type structure#

    struct arc_inform_type{T}
      status::Int32
      alloc_status::Int32
      bad_alloc::NTuple{81,Cchar}
      iter::Int32
      cg_iter::Int32
      f_eval::Int32
      g_eval::Int32
      h_eval::Int32
      factorization_status::Int32
      factorization_max::Int32
      max_entries_factors::Int64
      factorization_integer::Int64
      factorization_real::Int64
      factorization_average::T
      obj::T
      norm_g::T
      weight::T
      time::arc_time_type{T}
      rqs_inform::rqs_inform_type{T}
      glrt_inform::glrt_inform_type{T}
      dps_inform::dps_inform_type{T}
      psls_inform::psls_inform_type{T}
      lms_inform::lms_inform_type{T}
      lms_inform_prec::lms_inform_type{T}
      sha_inform::sha_inform_type

detailed documentation#

inform derived type as a Julia structure

components#

Int32 status

return status. See ARC_solve for details

Int32 alloc_status

the status of the last attempted allocation/deallocation

NTuple{81,Cchar} bad_alloc

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

Int32 iter

the total number of iterations performed

Int32 cg_iter

the total number of CG iterations performed

Int32 f_eval

the total number of evaluations of the objective function

Int32 g_eval

the total number of evaluations of the gradient of the objective functio

Int32 h_eval

the total number of evaluations of the Hessian of the objective function

Int32 factorization_status

the return status from the factorization

Int32 factorization_max

the maximum number of factorizations in a sub-problem solve

Int64 max_entries_factors

the maximum number of entries in the factors

Int64 factorization_integer

the total integer workspace required for the factorization

Int64 factorization_real

the total real workspace required for the factorization

T factorization_average

the average number of factorizations per sub-problem solve

T obj

the value of the objective function at the best estimate of the solution determined by the package.

T norm_g

the norm of the gradient of the objective function at the best estimate of the solution determined by the package.

T weight

the current value of the regularization weight

struct arc_time_type time

timings (see above)

struct rqs_inform_type rqs_inform

inform parameters for RQS

struct glrt_inform_type glrt_inform

inform parameters for GLRT

struct dps_inform_type dps_inform

inform parameters for DPS

struct psls_inform_type psls_inform

inform parameters for PSLS

struct lms_inform_type lms_inform

inform parameters for LMS

struct lms_inform_type lms_inform_prec

inform parameters for LMS used for preconditioning

struct sha_inform_type sha_inform

inform parameters for SHA

example calls#

This is an example of how to use the package to minimize a multi-dimensional objective function; the code is available in $GALAHAD/src/arc/Julia/test_arc.jl . A variety of supported Hessian and constraint matrix storage formats are shown.

# test_arc.jl
# Simple code to test the Julia interface to ARC

using GALAHAD
using Test
using Printf
using Accessors

# Custom userdata struct
struct userdata_arc{T}
  p::T
end

function test_arc(::Type{T}) where T

  # Objective function
  function fun(n::Int, x::Vector{T}, f::Ref{T}, userdata::userdata_arc)
    p = userdata.p
    f[] = (x[1] + x[3] + p)^2 + (x[2] + x[3])^2 + cos(x[1])
    return 0
  end

  # Gradient of the objective
  function grad(n::Int, x::Vector{T}, g::Vector{T}, userdata::userdata_arc)
    p = userdata.p
    g[1] = 2.0 * (x[1] + x[3] + p) - sin(x[1])
    g[2] = 2.0 * (x[2] + x[3])
    g[3] = 2.0 * (x[1] + x[3] + p) + 2.0 * (x[2] + x[3])
    return 0
  end

  # Hessian of the objective
  function hess(n::Int, ne::Int, x::Vector{T}, hval::Vector{T},
                userdata::userdata_arc)
    hval[1] = 2.0 - cos(x[1])
    hval[2] = 2.0
    hval[3] = 2.0
    hval[4] = 2.0
    hval[5] = 4.0
    return 0
  end

  # Dense Hessian
  function hess_dense(n::Int, ne::Int, x::Vector{T}, hval::Vector{T},
                      userdata::userdata_arc)
    hval[1] = 2.0 - cos(x[1])
    hval[2] = 0.0
    hval[3] = 2.0
    hval[4] = 2.0
    hval[5] = 2.0
    hval[6] = 4.0
    return 0
  end

  # Hessian-vector product
  function hessprod(n::Int, x::Vector{T}, u::Vector{T}, v::Vector{T},
                    got_h::Bool, userdata::userdata_arc)
    u[1] = u[1] + 2.0 * (v[1] + v[3]) - cos(x[1]) * v[1]
    u[2] = u[2] + 2.0 * (v[2] + v[3])
    u[3] = u[3] + 2.0 * (v[1] + v[2] + 2.0 * v[3])
    return 0
  end

  # Apply preconditioner
  function prec(n::Int, x::Vector{T}, u::Vector{T}, v::Vector{T},
                userdata::userdata_arc)
    u[1] = 0.5 * v[1]
    u[2] = 0.5 * v[2]
    u[3] = 0.25 * v[3]
    return 0
  end

  # Objective function
  function fun_diag(n::Int, x::Vector{T}, f::Ref{T}, userdata::userdata_arc)
    p = userdata.p
    f[] = (x[3] + p)^2 + x[2]^2 + cos(x[1])
    return 0
  end

  # Gradient of the objective
  function grad_diag(n::Int, x::Vector{T}, g::Vector{T}, userdata::userdata_arc)
    p = userdata.p
    g[1] = -sin(x[1])
    g[2] = 2.0 * x[2]
    g[3] = 2.0 * (x[3] + p)
    return 0
  end

  # Hessian of the objective
  function hess_diag(n::Int, ne::Int, x::Vector{T}, hval::Vector{T},
                     userdata::userdata_arc)
    hval[1] = -cos(x[1])
    hval[2] = 2.0
    hval[3] = 2.0
    return 0
  end

  # Hessian-vector product
  function hessprod_diag(n::Int, x::Vector{T}, u::Vector{T}, v::Vector{T},
                         got_h::Bool, userdata::userdata_arc)
    u[1] = u[1] + -cos(x[1]) * v[1]
    u[2] = u[2] + 2.0 * v[2]
    u[3] = u[3] + 2.0 * v[3]
    return 0
  end

  # Derived types
  data = Ref{Ptr{Cvoid}}()
  control = Ref{arc_control_type{T}}()
  inform = Ref{arc_inform_type{T}}()

  # Set user data
  userdata = userdata_arc(4.0)

  # Set problem data
  n = 3 # dimension
  ne = 5 # Hesssian elements
  H_row = Cint[1, 2, 3, 3, 3]  # Hessian H
  H_col = Cint[1, 2, 1, 2, 3]  # NB lower triangle
  H_ptr = Cint[1, 2, 3, 6]  # row pointers

  # Set storage
  g = zeros(T, n) # gradient
  st = ' '
  status = Ref{Cint}()

  @printf(" Fortran sparse matrix indexing\n\n")
  @printf(" tests reverse-communication options\n\n")

  # reverse-communication input/output
  eval_status = Ref{Cint}()
  f = Ref{T}(0.0)
  u = zeros(T, n)
  v = zeros(T, n)
  H_val = zeros(T, ne)
  H_dense = zeros(T, div(n * (n + 1), 2))
  H_diag = zeros(T, n)

  for d in 1:5

    # Initialize ARC
    arc_initialize(T, data, control, status)

    # Set user-defined control options
    @reset control[].f_indexing = true # Fortran sparse matrix indexing
    # @reset control[].print_level = Cint(1)

    # Start from 1.5
    x = T[1.5, 1.5, 1.5]

    # sparse co-ordinate storage
    if d == 1
      st = 'C'
      arc_import(T, control, data, status, n, "coordinate",
                 ne, H_row, H_col, C_NULL)

      terminated = false
      while !terminated # reverse-communication loop
        arc_solve_reverse_with_mat(T, data, status, eval_status, n, x, f[], g, ne, H_val, u, v)
        if status[] == 0 # successful termination
          terminated = true
        elseif status[] < 0 # error exit
          terminated = true
        elseif status[] == 2 # evaluate f
          eval_status[] = fun(n, x, f, userdata)
        elseif status[] == 3 # evaluate g
          eval_status[] = grad(n, x, g, userdata)
        elseif status[] == 4 # evaluate H
          eval_status[] = hess(n, ne, x, H_val, userdata)
        elseif status[] == 6 # evaluate the product with P
          eval_status[] = prec(n, x, u, v, userdata)
        else
          @printf(" the value %1i of status should not occur\n", status)
        end
      end
    end

    # sparse by rows
    if d == 2
      st = 'R'
      arc_import(T, control, data, status, n, "sparse_by_rows", ne,
                 C_NULL, H_col, H_ptr)

      terminated = false
      while !terminated # reverse-communication loop
        arc_solve_reverse_with_mat(T, data, status, eval_status,
                                   n, x, f[], g, ne, H_val, u, v)
        if status[] == 0 # successful termination
          terminated = true
        elseif status[] < 0 # error exit
          terminated = true
        elseif status[] == 2 # evaluate f
          eval_status[] = fun(n, x, f, userdata)
        elseif status[] == 3 # evaluate g
          eval_status[] = grad(n, x, g, userdata)
        elseif status[] == 4 # evaluate H
          eval_status[] = hess(n, ne, x, H_val, userdata)
        elseif status[] == 6 # evaluate the product with P
          eval_status[] = prec(n, x, u, v, userdata)
        else
          @printf(" the value %1i of status should not occur\n", status)
        end
      end
    end

    # dense
    if d == 3
      st = 'D'
      arc_import(T, control, data, status, n, "dense",
                 ne, C_NULL, C_NULL, C_NULL)

      terminated = false
      while !terminated # reverse-communication loop
        arc_solve_reverse_with_mat(T, data, status, eval_status,
                                   n, x, f[], g, div(n * (n + 1), 2), H_dense, u, v)
        if status[] == 0 # successful termination
          terminated = true
        elseif status[] < 0 # error exit
          terminated = true
        elseif status[] == 2 # evaluate f
          eval_status[] = fun(n, x, f, userdata)
        elseif status[] == 3 # evaluate g
          eval_status[] = grad(n, x, g, userdata)
        elseif status[] == 4 # evaluate H
          eval_status[] = hess_dense(n, div(n * (n + 1), 2), x, H_dense, userdata)
        elseif status[] == 6 # evaluate the product with P
          eval_status[] = prec(n, x, u, v, userdata)
        else
          @printf(" the value %1i of status should not occur\n", status)
        end
      end
    end

    # diagonal
    if d == 4
      st = 'I'
      arc_import(T, control, data, status, n, "diagonal",
                 ne, C_NULL, C_NULL, C_NULL)

      terminated = false
      while !terminated # reverse-communication loop
        arc_solve_reverse_with_mat(T, data, status, eval_status,
                                   n, x, f[], g, n, H_diag, u, v)
        if status[] == 0 # successful termination
          terminated = true
        elseif status[] < 0 # error exit
          terminated = true
        elseif status[] == 2 # evaluate f
          eval_status[] = fun_diag(n, x, f, userdata)
        elseif status[] == 3 # evaluate g
          eval_status[] = grad_diag(n, x, g, userdata)
        elseif status[] == 4 # evaluate H
          eval_status[] = hess_diag(n, n, x, H_diag, userdata)
        elseif status[] == 6 # evaluate the product with P
          eval_status[] = prec(n, x, u, v, userdata)
        else
          @printf(" the value %1i of status should not occur\n", status)
        end
      end
    end

    # access by products
    if d == 5
      st = 'P'
      arc_import(T, control, data, status, n, "absent",
                 ne, C_NULL, C_NULL, C_NULL)

      terminated = false
      while !terminated # reverse-communication loop
        arc_solve_reverse_without_mat(T, data, status, eval_status,
                                      n, x, f[], g, u, v)
        if status[] == 0 # successful termination
          terminated = true
        elseif status[] < 0 # error exit
          terminated = true
        elseif status[] == 2 # evaluate f
          eval_status[] = fun(n, x, f, userdata)
        elseif status[] == 3 # evaluate g
          eval_status[] = grad(n, x, g, userdata)
        elseif status[] == 5 # evaluate H
          eval_status[] = hessprod(n, x, u, v, false, userdata)
        elseif status[] == 6 # evaluate the product with P
          eval_status[] = prec(n, x, u, v, userdata)
        else
          @printf(" the value %1i of status should not occur\n", status)
        end
      end
    end

    arc_information(T, data, inform, status)

    if inform[].status[] == 0
      @printf("%c:%6i iterations. Optimal objective value = %5.2f status = %1i\n", st,
              inform[].iter, inform[].obj, inform[].status)
    else
      @printf("%c: ARC_solve exit status = %1i\n", st, inform[].status)
    end

    # @printf("x: ")
    # for i = 1:n
    #   @printf("%f ", x[i])
    # end
    # @printf("\n")
    # @printf("gradient: ")
    # for i = 1:n
    #   @printf("%f ", g[i])
    # end
    # @printf("\n")

    # Delete internal workspace
    arc_terminate(T, data, control, inform)
  end
  return 0
end

@testset "ARC" begin
  @test test_arc(Float32) == 0
  @test test_arc(Float64) == 0
end