GALAHAD DGO package#
purpose#
The dgo
package uses a deterministic partition-and-bound trust-region
method to find an approximation to the global minimizer of a
differentiable objective function \(f(x)\) of n variables \(x\),
subject to simple bounds \(x^l <= x <= x^u\) on the variables.
Here, any of the components of the vectors of bounds \(x^l\) and \(x^u\)
may be infinite. The method offers the choice of direct and
iterative solution of the key trust-region subproblems, and
is 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.
Although there are theoretical guarantees, these may require a large
number of evaluations as the dimension and nonconvexity increase.
The alternative package bgo
may sometimes be preferred.
See Section 4 of $GALAHAD/doc/dgo.pdf for additional details.
method#
Starting with the initial box \(x^l \leq x \leq x^u\), a sequence of
boxes is generated by considering the current set, and partitioning
a promising candidate into three equally-sized sub-boxes by splitting
along one of the box dimensions. Each partition requires only a pair of
new function and derivative evaluations, and these values, together with
estimates of Lipschitz constants, makes it possible to remove other boxes
from further consideration as soon as they cannot contain a global minimizer.
Efficient control of the dictionary of vertices of the sub-boxes
is handled using a suitable hashing procedure provided by
HASH
; each sub-box is indexed by the concatenated
coordinates of a pair of opposite vertices. At various
stages, local minimization in a promising sub-box, using
TRB
, may be used to improve the best-known upper bound
on the global minimizer.
If \(n=1\), the specialised univariate global minimization package
UGO
is called directly.
We reiterate that although there are theoretical guarantees, these may require a large number of evaluations as the dimension and nonconvexity increase. Thus the method should best be viewed as a heuristic to try to find a reasonable approximation of the global minimum.
references#
The global minimization method employed is an extension of that due to
Ya. D. Sergeyev and D. E. Kasov, ``A deterministic global optimization using smooth diagonal auxiliary functions’’, Communications in Nonlinear Science and Numerical Simulation, 21(1-3) (2015) 99-111.
but adapted to use 2nd derivatives, while in the special case when \(n=1\), a simplification based on the ideas in
D. Lera and Ya. D. Sergeyev (2013), ``Acceleration of univariate global optimization algorithms working with Lipschitz functions and Lipschitz first derivatives’’ SIAM J. Optimization 23(1) (2013) 508–529.
is used instead. The generic bound-constrained trust-region method used for local minimization is described in detail in
A. R. Conn, N. I. M. Gould and Ph. L. Toint, Trust-region methods. SIAM/MPS Series on Optimization (2000).
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 dgo package must be called in the following order:
dgo_initialize - provide default control parameters and set up initial data structures
dgo_read_specfile (optional) - override control values by reading replacement values from a file
dgo_import - set up problem data structures and fixed values
dgo_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved
solve the problem by calling one of
dgo_solve_with_mat - solve using function calls to evaluate function, gradient and Hessian values
dgo_solve_without_mat - solve using function calls to evaluate function and gradient values and Hessian-vector products
dgo_solve_reverse_with_mat - solve returning to the calling program to obtain function, gradient and Hessian values, or
dgo_solve_reverse_without_mat - solve returning to the calling prorgram to obtain function and gradient values and Hessian-vector products
dgo_information (optional) - recover information about the solution and solution process
dgo_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 dgo_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 dgo_control_type) |
status |
is a scalar variable of type Int32 that gives the exit status from the package. Possible values are (currently):
|
function dgo_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/dgo/DGO.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/dgo.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 dgo_control_type) |
specfile |
is a one-dimensional array of type Vararg{Cchar} that must give the name of the specification file |
function dgo_import(T, control, data, status, n, x_l, x_u, 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 dgo_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:
|
n |
is a scalar variable of type Int32 that holds the number of variables** |
x_l |
is a one-dimensional array of size n and type T that holds the values \(x^l\) of the lower bounds on the optimization variables \(x\). The j-th component of x_l, \(j = 1, \ldots, n\), contains \(x^l_j\). |
x_u |
is a one-dimensional array of size n and type T that holds the values \(x^u\) of the upper bounds on the optimization variables \(x\). The j-th component of x_u, \(j = 1, \ldots, n\), contains \(x^u_j\). |
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 dgo_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 dgo_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:
|
function dgo_solve_with_mat(T, data, userdata, status, n, x, g, ne, eval_f, eval_g, eval_h, eval_hprod, eval_prec)
Find an approximation to the global minimizer of a given function subject to simple bounds on the variables using a partition-and-bound trust-region 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:
|
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 |
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 |
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_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_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
dgo_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_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 |
function dgo_solve_without_mat(T, data, userdata, status, n, x, g, eval_f, eval_g, eval_hprod, eval_shprod, eval_prec)
Find an approximation to the global minimizer of a given function subject to simple bounds on the variables using a partition-and-bound trust-region 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:
|
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 |
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 |
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_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_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_shprod |
is a user-supplied function that must have the following signature: function eval_shprod(n, x, nnz_v, index_nz_v, v, nnz_u, index_nz_u, u, got_h, userdata) The product \(u = \nabla_{xx}f(x) v\) of the Hessian
\(\nabla_{xx}f(x)\) of the objective function evaluated
at \(x\) with the sparse vector v=\(v\) must be returned
in u, and the function return value set to 0. Only the
components index_nz_v[0:nnz_v-1] of v are nonzero, and
the remaining components may not have been be set. On
exit, the user must indicate the nnz_u indices of u
that are nonzero in index_nz_u[0:nnz_u-1], and only
these components of u need be set. 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_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 |
function dgo_solve_reverse_with_mat(T, data, status, eval_status, n, x, f, g, ne, H_val, u, v)
Find an approximation to the global minimizer of a given function subject to simple bounds on the variables using a partition-and-bound trust-region 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:
|
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 |
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 |
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 dgo_solve_reverse_without_mat(T, data, status, eval_status, n, x, f, g, u, v, index_nz_v, nnz_v, index_nz_u, nnz_u)
Find an approximation to the global minimizer of a given function subject to simple bounds on the variables using a partition-and-bound trust-region 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:
|
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 |
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 |
u |
is a one-dimensional array of size n and type T that is used for reverse communication (see status=5,6,7 above for details) |
v |
is a one-dimensional array of size n and type T that is used for reverse communication (see status=5,6,7 above for details) |
index_nz_v |
is a one-dimensional array of size n and type Int32 that is used for reverse communication (see status=7 above for details) |
nnz_v |
is a scalar variable of type Int32 that is used for reverse communication (see status=7 above for details) |
index_nz_u |
is a one-dimensional array of size n and type Int32 that is used for reverse communication (see status=7 above for details) |
nnz_u |
is a scalar variable of type Int32 that is used for reverse communication (see status=7 above for details). On initial (status=1) entry, nnz_u should be set to an (arbitrary) nonzero value, and nnz_u=0 is recommended |
function dgo_information(T, data, inform, status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a structure containing output information (see dgo_inform_type) |
status |
is a scalar variable of type Int32 that gives the exit status from the package. Possible values are (currently):
|
function dgo_terminate(T, data, control, inform)
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a structure containing control information (see dgo_control_type) |
inform |
is a structure containing output information (see dgo_inform_type) |
available structures#
dgo_control_type structure#
struct dgo_control_type{T} f_indexing::Bool error::Int32 out::Int32 print_level::Int32 start_print::Int32 stop_print::Int32 print_gap::Int32 maxit::Int32 max_evals::Int32 dictionary_size::Int32 alive_unit::Int32 alive_file::NTuple{31,Cchar} infinity::T lipschitz_lower_bound::T lipschitz_reliability::T lipschitz_control::T stop_length::T stop_f::T obj_unbounded::T cpu_time_limit::T clock_time_limit::T hessian_available::Bool prune::Bool perform_local_optimization::Bool space_critical::Bool deallocate_error_fatal::Bool prefix::NTuple{31,Cchar} hash_control::hash_control_type ugo_control::ugo_control_type{T} trb_control::trb_control_type{T}
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. Possible values are:
\(\leq\) 0 no output,
1 a one-line summary for every improvement
2 a summary of each iteration
\(\geq\) 3 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 max_evals
the maximum number of function evaluations made
Int32 dictionary_size
the size of the initial hash dictionary
Int32 alive_unit
removal of the file alive_file from unit alive_unit terminates execution
char alive_file[31]
see alive_unit
T infinity
any bound larger than infinity in modulus will be regarded as infinite
T lipschitz_lower_bound
a small positive constant (<= 1e-6) that ensure that the estimted gradient Lipschitz constant is not too small
T lipschitz_reliability
the Lipschitz reliability parameter, the Lipschiz constant used will be a factor lipschitz_reliability times the largest value observed
T lipschitz_control
the reliablity control parameter, the actual reliability parameter used will be .lipschitz_reliability
MAX( 1, n - 1 ) * .lipschitz_control / iteration
T stop_length
the iteration will stop if the length, delta, of the diagonal in the box with the smallest-found objective function is smaller than .stop_length times that of the original bound box, delta_0
T stop_f
the iteration will stop if the gap between the best objective value found and the smallest lower bound is smaller than .stop_f
T obj_unbounded
the smallest value the objective 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 prune
should boxes that cannot contain the global minimizer be pruned (i.e., removed from further consideration)?
Bool perform_local_optimization
should approximate minimizers be impoved by judicious local minimization?
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 hash_control_type hash_control
control parameters for HASH
struct ugo_control_type ugo_control
control parameters for UGO
struct trb_control_type trb_control
control parameters for TRB
dgo_time_type structure#
struct dgo_time_type{T} total::Float32 univariate_global::Float32 multivariate_local::Float32 clock_total::T clock_univariate_global::T clock_multivariate_local::T
detailed documentation#
time derived type as a Julia structure
components#
Float32 total
the total CPU time spent in the package
Float32 univariate_global
the CPU time spent performing univariate global optimization
Float32 multivariate_local
the CPU time spent performing multivariate local optimization
T clock_total
the total clock time spent in the package
T clock_univariate_global
the clock time spent performing univariate global optimization
T clock_multivariate_local
the clock time spent performing multivariate local optimization
dgo_inform_type structure#
struct dgo_inform_type{T} status::Int32 alloc_status::Int32 bad_alloc::NTuple{81,Cchar} iter::Int32 f_eval::Int32 g_eval::Int32 h_eval::Int32 obj::T norm_pg::T length_ratio::T f_gap::T why_stop::NTuple{2,Cchar} time::dgo_time_type{T} hash_inform::hash_inform_type ugo_inform::ugo_inform_type{T} trb_inform::trb_inform_type{T}
detailed documentation#
inform derived type as a Julia structure
components#
Int32 status
return status. See DGO_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 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 function
Int32 h_eval
the total number of evaluations of the Hessian of the objective function
T obj
the value of the objective function at the best estimate of the solution determined by DGO_solve
T norm_pg
the norm of the projected gradient of the objective function at the best estimate of the solution determined by DGO_solve
T length_ratio
the ratio of the final to the initial box lengths
T f_gap
the gap between the best objective value found and the lowest bound
char why_stop[2]
why did the iteration stop? This wil be ‘D’ if the box length is small enough, ‘F’ if the objective gap is small enough, and ‘ ‘ otherwise
struct dgo_time_type time
timings (see above)
struct hash_inform_type hash_inform
inform parameters for HASH
struct ugo_inform_type ugo_inform
inform parameters for UGO
struct trb_inform_type trb_inform
inform parameters for TRB
example calls#
This is an example of how to use the package to minimize a multi-dimensional objective; the code is available in $GALAHAD/src/dgo/Julia/test_dgo.jl . A variety of supported Hessian and constraint matrix storage formats are shown.
# test_dgo.jl
# Simple code to test the Julia interface to DGO
using GALAHAD
using Test
using Printf
using Accessors
# Custom userdata struct
struct userdata_dgo{T}
p::T
freq::T
mag::T
end
function test_dgo(::Type{T}) where T
# Objective function
function fun(n::Int, x::Vector{T}, f::Ref{T}, userdata::userdata_dgo)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
f[] = (x[1] + x[3] + p)^2 + (x[2] + x[3])^2 + mag * cos(freq * x[1]) + x[1] + x[2] + x[3]
return 0
end
# Gradient of the objective
function grad(n::Int, x::Vector{T}, g::Vector{T}, userdata::userdata_dgo)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
g[1] = 2.0 * (x[1] + x[3] + p) - mag * freq * sin(freq * x[1]) + 1.0
g[2] = 2.0 * (x[2] + x[3]) + 1.0
g[3] = 2.0 * (x[1] + x[3] + p) + 2.0 * (x[2] + x[3]) + 1.0
return 0
end
# Hessian of the objective
function hess(n::Int, ne::Int, x::Vector{T}, hval::Vector{T},
userdata::userdata_dgo)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
hval[1] = 2.0 - mag * freq^2 * cos(freq * x[1])
hval[2] = 2.0
hval[3] = 2.0
hval[4] = 4.0
return 0
end
# Dense Hessian
function hess_dense(n::Int, ne::Int, x::Vector{T}, hval::Vector{T},
userdata::userdata_dgo)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
hval[1] = 2.0 - mag * freq^2 * cos(freq * x[1])
hval[2] = 0.0
hval[3] = 2.0
hval[4] = 2.0
hval[5] = 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_dgo)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
u[1] = u[1] + 2.0 * (v[1] + v[3]) - mag * freq^2 * cos(freq * 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
# Sparse Hessian-vector product
function shessprod(n::Int, x::Vector{T}, nnz_v::Cint, index_nz_v::Vector{Cint},
v::Vector{T}, nnz_u::Ref{Cint}, index_nz_u::Vector{Cint},
u::Vector{T}, got_h::Bool, userdata::userdata_dgo)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
p = zeros(T, 3)
used = falses(3)
for i in 1:nnz_v
j = index_nz_v[i]
if j == 1
p[1] = p[1] + 2.0 * v[1] - mag * freq^2 * cos(freq * x[1]) * v[1]
used[1] = true
p[3] = p[3] + 2.0 * v[1]
used[3] = true
elseif j == 2
p[2] = p[2] + 2.0 * v[2]
used[2] = true
p[3] = p[3] + 2.0 * v[2]
used[3] = true
elseif j == 3
p[1] = p[1] + 2.0 * v[3]
used[1] = true
p[2] = p[2] + 2.0 * v[3]
used[2] = true
p[3] = p[3] + 4.0 * v[3]
used[3] = true
end
end
nnz_u[] = 0
for j in 1:3
if used[j]
u[j] = p[j]
nnz_u[] += 1
index_nz_u[nnz_u[]] = j
end
end
return 0
end
# Apply preconditioner
function prec(n::Int, x::Vector{T}, u::Vector{T}, v::Vector{T},
userdata::userdata_dgo)
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)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
f[] = (x[3] + p)^2 + x[2]^2 + mag * cos(freq * x[1]) + sum(x)
return 0
end
# Gradient of the objective
function grad_diag(n::Int, x::Vector{T}, g::Vector{T}, userdata)
p = userdata.p
freq = userdata.freq
mag = userdata.mag
g[1] = -mag * freq * sin(freq * x[1]) + 1
g[2] = 2.0 * x[2] + 1
g[3] = 2.0 * (x[3] + p) + 1
return 0
end
# Hessian of the objective
function hess_diag(n::Int, ne::Int, x::Vector{T}, hval::Vector{T}, userdata)
freq = userdata.freq
mag = userdata.mag
hval[1] = -mag * freq^2 * cos(freq * 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)
freq = userdata.freq
mag = userdata.mag
u[1] += -mag * freq^2 * cos(freq * x[1]) * v[1]
u[2] += 2.0 * v[2]
u[3] += 2.0 * v[3]
return 0
end
# Sparse Hessian-vector product
function shessprod_diag(n::Int, x::Vector{T}, nnz_v::Cint, index_nz_v::Vector{Cint},
v::Vector{T}, nnz_u::Ref{Cint}, index_nz_u::Vector{Cint},
u::Vector{T}, got_h::Bool, userdata)
freq = userdata.freq
mag = userdata.mag
p = zeros(3)
used = falses(3)
for i in 1:nnz_v
j = index_nz_v[i]
if j == 1
p[1] -= mag * freq^2 * cos(freq * x[1]) * v[1]
used[1] = true
elseif j == 2
p[2] += 2.0 * v[2]
used[2] = true
elseif j == 3
p[3] += 2.0 * v[3]
used[3] = true
end
end
nnz_u[] = 0
for j in 1:3
if used[j]
u[j] = p[j]
nnz_u[] += 1
index_nz_u[nnz_u[]] = j
end
end
return 0
end
# Derived types
data = Ref{Ptr{Cvoid}}()
control = Ref{dgo_control_type{T}}()
inform = Ref{dgo_inform_type{T}}()
# Set user data
userdata = userdata_dgo(4.0, 10.0, 1000.0)
# Set problem data
n = 3 # dimension
ne = 5 # Hesssian elements
x_l = T[-10, -10, -10]
x_u = T[0.5, 0.5, 0.5]
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}()
nnz_u = Ref{Cint}()
nnz_v = Ref{Cint}()
f = Ref{T}(0.0)
u = zeros(T, n)
v = zeros(T, n)
index_nz_u = zeros(Cint, n)
index_nz_v = zeros(Cint, 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 DGO
dgo_initialize(T, data, control, status)
# Set user-defined control options
@reset control[].f_indexing = true # Fortran sparse matrix indexing
@reset control[].maxit = Cint(2500)
# @reset control[].trb_control[].maxit = Cint(100)
# @reset control[].print_level = Cint(1)
# Start from 0
x = T[0.0, 0.0, 0.0]
# sparse co-ordinate storage
if d == 1
st = 'C'
dgo_import(T, control, data, status, n, x_l, x_u, "coordinate", ne, H_row, H_col, C_NULL)
terminated = false
while !terminated # reverse-communication loop
dgo_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[] == 5 # evaluate Hv product
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)
elseif status[] == 23 # evaluate f and g
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
elseif status[] == 25 # evaluate f and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 35 # evaluate g and Hv product
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 235 # evaluate f, g and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
else
@printf(" the value %1i of status should not occur\n", status)
end
end
end
# sparse by rows
if d == 2
st = 'R'
dgo_import(T, control, data, status, n, x_l, x_u,
"sparse_by_rows", ne, C_NULL, H_col, H_ptr)
terminated = false
while !terminated # reverse-communication loop
dgo_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[] == 5 # evaluate Hv product
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)
elseif status[] == 23 # evaluate f and g
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
elseif status[] == 25 # evaluate f and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 35 # evaluate g and Hv product
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 235 # evaluate f, g and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
else
@printf(" the value %1i of status should not occur\n", status)
end
end
end
# dense
if d == 3
st = 'D'
dgo_import(T, control, data, status, n, x_l, x_u,
"dense", ne, C_NULL, C_NULL, C_NULL)
terminated = false
while !terminated # reverse-communication loop
dgo_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[] == 5 # evaluate Hv product
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)
elseif status[] == 23 # evaluate f and g
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
elseif status[] == 25 # evaluate f and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 35 # evaluate g and Hv product
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 235 # evaluate f, g and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
else
@printf(" the value %1i of status should not occur\n", status)
end
end
end
# diagonal
if d == 4
st = 'I'
dgo_import(T, control, data, status, n, x_l, x_u,
"diagonal", ne, C_NULL, C_NULL, C_NULL)
terminated = false
while !terminated # reverse-communication loop
dgo_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[] == 5 # evaluate Hv product
eval_status[] = hessprod_diag(n, x, u, v, false,
userdata)
elseif status[] == 6 # evaluate the product with P
eval_status[] = prec(n, x, u, v, userdata)
elseif status[] == 23 # evaluate f and g
eval_status[] = fun_diag(n, x, f, userdata)
eval_status[] = grad_diag(n, x, g, userdata)
elseif status[] == 25 # evaluate f and Hv product
eval_status[] = fun_diag(n, x, f, userdata)
eval_status[] = hessprod_diag(n, x, u, v, false,
userdata)
elseif status[] == 35 # evaluate g and Hv product
eval_status[] = grad_diag(n, x, g, userdata)
eval_status[] = hessprod_diag(n, x, u, v, false,
userdata)
elseif status[] == 235 # evaluate f, g and Hv product
eval_status[] = fun_diag(n, x, f, userdata)
eval_status[] = grad_diag(n, x, g, userdata)
eval_status[] = hessprod_diag(n, x, u, v, false,
userdata)
else
@printf(" the value %1i of status should not occur\n", status)
end
end
end
# access by products
if d == 5
st = 'P'
dgo_import(T, control, data, status, n, x_l, x_u,
"absent", ne, C_NULL, C_NULL, C_NULL)
nnz_u = Ref{Cint}(0)
terminated = false
while !terminated # reverse-communication loop
dgo_solve_reverse_without_mat(T, data, status, eval_status,
n, x, f[], g, u, v, index_nz_v,
nnz_v, index_nz_u, nnz_u[])
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 Hv product
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)
elseif status[] == 7 # evaluate sparse Hess-vect product
eval_status[] = shessprod(n, x, nnz_v[], index_nz_v, v,
nnz_u, index_nz_u, u,
false, userdata)
elseif status[] == 23 # evaluate f and g
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
elseif status[] == 25 # evaluate f and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 35 # evaluate g and Hv product
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
elseif status[] == 235 # evaluate f, g and Hv product
eval_status[] = fun(n, x, f, userdata)
eval_status[] = grad(n, x, g, userdata)
eval_status[] = hessprod(n, x, u, v, false, userdata)
else
@printf(" the value %1i of status should not occur\n", status)
end
end
end
# Record solution information
dgo_information(T, data, inform, status)
if inform[].status[] == 0
@printf("%c:%6i evaluations. Optimal objective value = %5.2f status = %1i\n", st,
inform[].f_eval, inform[].obj, inform[].status)
elseif inform[].status[] == -18
@printf("%c:%6i evaluations. Best objective value = %5.2f status = %1i\n", st,
inform[].f_eval, inform[].obj, inform[].status)
else
@printf("%c: DGO_solve exit status = %1i\n", st, inform[].status)
end
# @printf("x: ")
# for i in 1:n
# @printf("%f ", x[i])
# end
# @printf("\n")
# @printf("gradient: ")
# for i in 1:n
# @printf("%f ", g[i])
# end
# @printf("\n")
# Delete internal workspace
dgo_terminate(T, data, control, inform)
end
return 0
end
@testset "DGO" begin
@test test_dgo(Float32) == 0
@test test_dgo(Float64) == 0
end