GALAHAD TRU package#
purpose#
The tru
package uses a trust-region 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/tru.pdf for additional details.
method#
A trust-region 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)\) within a trust region \(\|s_k\| \leq \Delta_k\) for some specified positive “radius” \(\Delta_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 radius is reduced by powers of a given reduction factor until it is smaller than \(\|s_k\|\). If the ratio is larger than \(\eta_v \geq \eta_d\), the radius will be increased so that it exceeds \(\|s_k\|\) by a given increase factor. 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 trust region will be used, but this may change if preconditioning is employed.
An approximate minimizer of the model within the trust region
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 TRS
or DPS
(depending on the norm chosen). The iterative approach
uses GLTR
, 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.
reference#
The generic trust-region method 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 \(0 <= j <= i <= n-1\)) need be held. In this case the lower triangle should be stored by rows, that is component \(i * i / 2 + j\) of the storage array H_val will hold the value \(H_{ij}\) (and, by symmetry, \(H_{ji}\)) for \(0 <= j <= i <= n-1\).
Sparse co-ordinate storage format: Only the nonzero entries of the matrices are stored. For the \(l\)-th entry, \(0 <= l <= ne-1\), of \(H\), its row index i, column index j and value \(H_{ij}\), \(0 <= j <= i <= n-1\), are stored as the \(l\)-th components of the integer arrays H_row and H_col and real array H_val, respectively, while the number of nonzeros is recorded as H_ne = \(ne\). Note that only the entries in the lower triangle should be stored.
Sparse row-wise storage format: Again only the nonzero entries are stored, but this time they are ordered so that those in row i appear directly before those in row i+1. For the i-th row of \(H\) the i-th component of the integer array H_ptr holds the position of the first entry in this row, while H_ptr(n) holds the total number of entries. The column indices j, \(0 <= j <= 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.
introduction to function calls#
To solve a given problem, functions from the tru package must be called in the following order:
tru_initialize - provide default control parameters and set up initial data structures
tru_read_specfile (optional) - override control values by reading replacement values from a file
tru_import - set up problem data structures and fixed values
tru_reset_control (optional) - possibly change control parameters if a sequence of problems are being solved
solve the problem by calling one of
tru_solve_with_mat - solve using function calls to evaluate function, gradient and Hessian values
tru_solve_without_mat - solve using function calls to evaluate function and gradient values and Hessian-vector products
tru_solve_reverse_with_mat - solve returning to the calling program to obtain function, gradient and Hessian values, or
tru_solve_reverse_without_mat - solve returning to the calling prorgram to obtain function and gradient values and Hessian-vector products
tru_information (optional) - recover information about the solution and solution process
tru_terminate - deallocate data structures
See the examples section for illustrations of use.
callable functions#
overview of functions provided#
// typedefs typedef float spc_; typedef double rpc_; typedef int ipc_; // structs struct tru_control_type; struct tru_inform_type; struct tru_time_type; // function calls void tru_initialize(void **data, struct tru_control_type* control, ipc_ *status); void tru_read_specfile(struct tru_control_type* control, const char specfile[]); void tru_import( struct tru_control_type* control, void **data, ipc_ *status, ipc_ n, const char H_type[], ipc_ ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[] ); void tru_reset_control( struct tru_control_type* control, void **data, ipc_ *status ); void tru_solve_with_mat( void **data, void *userdata, ipc_ *status, ipc_ n, rpc_ x[], rpc_ g[], ipc_ ne, ipc_(*)(ipc_, const rpc_[], rpc_*, const void*) eval_f, ipc_(*)(ipc_, const rpc_[], rpc_[], const void*) eval_g, ipc_(*)(ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_h, ipc_(*)(ipc_, const rpc_[], rpc_[], const rpc_[], const void*) eval_prec ); void tru_solve_without_mat( void **data, void *userdata, ipc_ *status, ipc_ n, rpc_ x[], rpc_ g[], ipc_(*)(ipc_, const rpc_[], rpc_*, const void*) eval_f, ipc_(*)(ipc_, const rpc_[], rpc_[], const void*) eval_g, ipc_(*)(ipc_, const rpc_[], rpc_[], const rpc_[], bool, const void*) eval_hprod, ipc_(*)(ipc_, const rpc_[], rpc_[], const rpc_[], const void*) eval_prec ); void tru_solve_reverse_with_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, rpc_ x[], rpc_ f, rpc_ g[], ipc_ ne, rpc_ H_val[], const rpc_ u[], rpc_ v[] ); void tru_solve_reverse_without_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, rpc_ x[], rpc_ f, rpc_ g[], rpc_ u[], rpc_ v[] ); void tru_information(void **data, struct tru_inform_type* inform, ipc_ *status); void tru_terminate( void **data, struct tru_control_type* control, struct tru_inform_type* inform );
typedefs#
typedef float spc_
spc_
is real single precision
typedef double rpc_
rpc_
is the real working precision used, but may be changed to float
by
defining the preprocessor variable REAL_32
or (if supported) to
__real128
using the variable REAL_128
.
typedef int ipc_
ipc_
is the default integer word length used, but may be changed to
int64_t
by defining the preprocessor variable INTEGER_64
.
function calls#
void tru_initialize(void **data, struct tru_control_type* control, ipc_ *status)
Set default control values and initialize private data
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see tru_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void tru_read_specfile(struct tru_control_type* control, const char specfile[])
Read the content of a specification file, and assign values associated with given keywords to the corresponding control parameters. An in-depth discussion of specification files is available, and a detailed list of keywords with associated default values is provided in $GALAHAD/src/tru/TRU.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/tru.pdf for a list of how these keywords relate to the components of the control structure.
Parameters:
control |
is a struct containing control information (see tru_control_type) |
specfile |
is a character string containing the name of the specification file |
void tru_import( struct tru_control_type* control, void **data, ipc_ *status, ipc_ n, const char H_type[], ipc_ ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[] )
Import problem data into internal storage prior to solution.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see tru_control_type) |
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are:
|
n |
is a scalar variable of type ipc_, that holds the number of variables |
H_type |
is a one-dimensional array of type char 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 ipc_, that holds the number of entries in the lower triangular part of H in the sparse co-ordinate storage scheme. It need not be set for any of the other three schemes. |
H_row |
is a one-dimensional array of size ne and type ipc_, that holds the row indices of the lower triangular part of H in the sparse co-ordinate storage scheme. It need not be set for any of the other three schemes, and in this case can be NULL |
H_col |
is a one-dimensional array of size ne and type ipc_, that holds the column indices of the lower triangular part of H in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the dense or diagonal storage schemes are used, and in this case can be NULL |
H_ptr |
is a one-dimensional array of size n+1 and type ipc_, that holds the starting position of each row of the lower triangular part of H, as well as the total number of entries, in the sparse row-wise storage scheme. It need not be set when the other schemes are used, and in this case can be NULL |
void tru_reset_control( struct tru_control_type* control, void **data, ipc_ *status )
Reset control parameters after import if required.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see tru_control_type) |
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are:
|
void tru_solve_with_mat( void **data, void *userdata, ipc_ *status, ipc_ n, rpc_ x[], rpc_ g[], ipc_ ne, ipc_(*)(ipc_, const rpc_[], rpc_*, const void*) eval_f, ipc_(*)(ipc_, const rpc_[], rpc_[], const void*) eval_g, ipc_(*)(ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_h, ipc_(*)(ipc_, const rpc_[], rpc_[], const rpc_[], const void*) eval_prec )
Find a local minimizer of a given function using a 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 ipc_, 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 ipc_, that holds the number of variables |
x |
is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\). |
g |
is a one-dimensional array of size n and type rpc_, that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 0, … , n-1, contains \(g_j\). |
ne |
is a scalar variable of type ipc_, 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: ipc_ eval_f( ipc_ n, const rpc_ x[], rpc_ *f, const void *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: ipc_ eval_g( ipc_ n, const rpc_ x[], rpc_ g[], const void *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: ipc_ eval_h( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ h[], const void *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 tru_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 NULL. If non-NULL, it must have the following signature: ipc_ eval_prec( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[], const void *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 |
void tru_solve_without_mat( void **data, void *userdata, ipc_ *status, ipc_ n, rpc_ x[], rpc_ g[], ipc_(*)(ipc_, const rpc_[], rpc_*, const void*) eval_f, ipc_(*)(ipc_, const rpc_[], rpc_[], const void*) eval_g, ipc_(*)(ipc_, const rpc_[], rpc_[], const rpc_[], bool, const void*) eval_hprod, ipc_(*)(ipc_, const rpc_[], rpc_[], const rpc_[], const void*) eval_prec )
Find a local minimizer of a given function using a 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 ipc_, 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 ipc_, that holds the number of variables |
x |
is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\). |
g |
is a one-dimensional array of size n and type rpc_, that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 0, … , n-1, contains \(g_j\). |
eval_f |
is a user-supplied function that must have the following signature: ipc_ eval_f( ipc_ n, const rpc_ x[], rpc_ *f, const void *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: ipc_ eval_g( ipc_ n, const rpc_ x[], rpc_ g[], const void *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: ipc_ eval_hprod( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[], bool got_h, const void *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_prec |
is an optional user-supplied function that may be NULL. If non-NULL, it must have the following signature: ipc_ eval_prec( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[], const void *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 |
void tru_solve_reverse_with_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, rpc_ x[], rpc_ f, rpc_ g[], ipc_ ne, rpc_ H_val[], const rpc_ u[], rpc_ v[] )
Find a local minimizer of a given function using a 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 ipc_, 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 ipc_, that is used to indicate if objective function/gradient/Hessian values can be provided (see above) |
n |
is a scalar variable of type ipc_, that holds the number of variables |
x |
is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\). |
f |
is a scalar variable pointer of type rpc_, that holds the value of the objective function. |
g |
is a one-dimensional array of size n and type rpc_, that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 0, … , n-1, contains \(g_j\). |
ne |
is a scalar variable of type ipc_, that holds the number of entries in the lower triangular part of the Hessian matrix \(H\). |
H_val |
is a one-dimensional array of size ne and type rpc_, that holds the values of the entries of the lower triangular part of the Hessian matrix \(H\) in any of the available storage schemes. |
u |
is a one-dimensional array of size n and type rpc_, that is used for reverse communication (see above for details) |
v |
is a one-dimensional array of size n and type rpc_, that is used for reverse communication (see above for details) |
void tru_solve_reverse_without_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, rpc_ x[], rpc_ f, rpc_ g[], rpc_ u[], rpc_ v[] )
Find a local minimizer of a given function using a 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 ipc_, 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 ipc_, that is used to indicate if objective function/gradient/Hessian values can be provided (see above) |
n |
is a scalar variable of type ipc_, that holds the number of variables |
x |
is a one-dimensional array of size n and type rpc_, that holds the values \(x\) of the optimization variables. The j-th component of x, j = 0, … , n-1, contains \(x_j\). |
f |
is a scalar variable pointer of type rpc_, that holds the value of the objective function. |
g |
is a one-dimensional array of size n and type rpc_, that holds the gradient \(g = \nabla_xf(x)\) of the objective function. The j-th component of g, j = 0, … , n-1, contains \(g_j\). |
u |
is a one-dimensional array of size n and type rpc_, that is used for reverse communication (see above for details) |
v |
is a one-dimensional array of size n and type rpc_, that is used for reverse communication (see above for details) |
void tru_information(void **data, struct tru_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see tru_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void tru_terminate( void **data, struct tru_control_type* control, struct tru_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see tru_control_type) |
inform |
is a struct containing output information (see tru_inform_type) |
available structures#
tru_control_type structure#
#include <galahad_tru.h> struct tru_control_type { // components bool f_indexing; ipc_ error; ipc_ out; ipc_ print_level; ipc_ start_print; ipc_ stop_print; ipc_ print_gap; ipc_ maxit; ipc_ alive_unit; char alive_file[31]; ipc_ non_monotone; ipc_ model; ipc_ norm; ipc_ semi_bandwidth; ipc_ lbfgs_vectors; ipc_ max_dxg; ipc_ icfs_vectors; ipc_ mi28_lsize; ipc_ mi28_rsize; rpc_ stop_g_absolute; rpc_ stop_g_relative; rpc_ stop_s; ipc_ advanced_start; rpc_ initial_radius; rpc_ maximum_radius; rpc_ eta_successful; rpc_ eta_very_successful; rpc_ eta_too_successful; rpc_ radius_increase; rpc_ radius_reduce; rpc_ radius_reduce_max; rpc_ obj_unbounded; rpc_ cpu_time_limit; rpc_ clock_time_limit; bool hessian_available; bool subproblem_direct; bool retrospective_trust_region; bool renormalize_radius; bool space_critical; bool deallocate_error_fatal; char prefix[31]; struct trs_control_type trs_control; struct gltr_control_type gltr_control; struct dps_control_type dps_control; struct psls_control_type psls_control; struct lms_control_type lms_control; struct lms_control_type lms_control_prec; struct sec_control_type sec_control; struct sha_control_type sha_control; };
detailed documentation#
control derived type as a C struct
components#
bool f_indexing
use C or Fortran sparse matrix indexing
ipc_ error
error and warning diagnostics occur on stream error
ipc_ out
general output occurs on stream out
ipc_ print_level
the level of output required.
\(\leq\) 0 gives no output,
= 1 gives a one-line summary for every iteration,
= 2 gives a summary of the inner iteration for each iteration,
\(\geq\) 3 gives increasingly verbose (debugging) output
ipc_ start_print
any printing will start on this iteration
ipc_ stop_print
any printing will stop on this iteration
ipc_ print_gap
the number of iterations between printing
ipc_ maxit
the maximum number of iterations allowed
ipc_ alive_unit
removal of the file alive_file from unit alive_unit terminates execution
char alive_file[31]
see alive_unit
ipc_ 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.
ipc_ 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 (sparsity-based)
5 secant second-order (limited-memory BFGS, with .lbfgs_vectors history)
6 secant second-order (limited-memory SR1, with .lbfgs_vectors history)
ipc_ norm
the trust-region 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)
ipc_ semi_bandwidth
specify the semi-bandwidth of the band matrix \(P\) if required
ipc_ lbfgs_vectors
number of vectors used by the L-BFGS matrix \(P\) if required
ipc_ max_dxg
number of vectors used by the sparsity-based secant Hessian if required
ipc_ icfs_vectors
number of vectors used by the Lin-More’ incomplete factorization matrix \(P\) if required
ipc_ 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
ipc_ 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
rpc_ 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
rpc_ stop_g_relative
see stop_g_absolute
rpc_ stop_s
see stop_g_absolute
ipc_ advanced_start
try to pick a good initial trust-region radius using .advanced_start iterates of a variant on the strategy of Sartenaer SISC 18(6) 1990:1788-1803
rpc_ initial_radius
initial value for the trust-region radius
rpc_ maximum_radius
maximum permitted trust-region radius
rpc_ 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 trust-region radius will be increased if this relative decrease is greater than .eta_very_successful but smaller than .eta_too_successful
rpc_ eta_very_successful
see eta_successful
rpc_ eta_too_successful
see eta_successful
rpc_ radius_increase
on very successful iterations, the trust-region radius will be increased by the factor .radius_increase, while if the iteration is unsucceful, the radius will be decreased by a factor .radius_reduce but no more than .radius_reduce_max
rpc_ radius_reduce
see radius_increase;
rpc_ radius_reduce_max
see radius_increase;
rpc_ obj_unbounded
the smallest value the objective function may take before the problem is marked as unbounded
rpc_ cpu_time_limit
the maximum CPU time allowed (-ve means infinite)
rpc_ clock_time_limit
the maximum elapsed clock time allowed (-ve means infinite)
bool hessian_available
is the Hessian matrix of second derivatives available or is access only via matrix-vector products?
bool subproblem_direct
use a direct (factorization) or (preconditioned) iterative method to find the search direction
bool retrospective_trust_region
is a retrospective strategy to be used to update the trust-region radius?
bool renormalize_radius
should the radius be renormalized to account for a change in preconditioner?
bool space_critical
if .space_critical is true, every effort will be made to use as little space as possible. This may result in longer computation time
bool deallocate_error_fatal
if .deallocate_error_fatal is true, any array/pointer deallocation error will terminate execution. Otherwise, computation will continue
char prefix[31]
all output lines will be prefixed by .prefix(2:LEN(TRIM(.prefix))-1) where .prefix contains the required string enclosed in quotes, e.g. “string” or ‘string’
struct trs_control_type trs_control
control parameters for TRS
struct gltr_control_type gltr_control
control parameters for GLTR
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 sec_control_type sec_control
control parameters for SEC
struct sha_control_type sha_control
control parameters for SHA
tru_time_type structure#
#include <galahad_tru.h> struct tru_time_type { // components spc_ total; spc_ preprocess; spc_ analyse; spc_ factorize; spc_ solve; rpc_ clock_total; rpc_ clock_preprocess; rpc_ clock_analyse; rpc_ clock_factorize; rpc_ clock_solve; };
detailed documentation#
time derived type as a C struct
components#
spc_ total
the total CPU time spent in the package
spc_ preprocess
the CPU time spent preprocessing the problem
spc_ analyse
the CPU time spent analysing the required matrices prior to factorization
spc_ factorize
the CPU time spent factorizing the required matrices
spc_ solve
the CPU time spent computing the search direction
rpc_ clock_total
the total clock time spent in the package
rpc_ clock_preprocess
the clock time spent preprocessing the problem
rpc_ clock_analyse
the clock time spent analysing the required matrices prior to factorization
rpc_ clock_factorize
the clock time spent factorizing the required matrices
rpc_ clock_solve
the clock time spent computing the search direction
tru_inform_type structure#
#include <galahad_tru.h> struct tru_inform_type { // components ipc_ status; ipc_ alloc_status; char bad_alloc[81]; ipc_ iter; ipc_ cg_iter; ipc_ f_eval; ipc_ g_eval; ipc_ h_eval; ipc_ factorization_max; ipc_ factorization_status; int64_t max_entries_factors; int64_t factorization_integer; int64_t factorization_real; rpc_ factorization_average; rpc_ obj; rpc_ norm_g; rpc_ radius; struct tru_time_type time; struct trs_inform_type trs_inform; struct gltr_inform_type gltr_inform; struct dps_inform_type dps_inform; struct psls_inform_type psls_inform; struct lms_inform_type lms_inform; struct lms_inform_type lms_inform_prec; struct sec_inform_type sec_inform; struct sha_inform_type sha_inform; };
detailed documentation#
inform derived type as a C struct
components#
ipc_ status
return status. See TRU_solve for details
ipc_ alloc_status
the status of the last attempted allocation/deallocation
char bad_alloc[81]
the name of the array for which an allocation/deallocation error occurred
ipc_ iter
the total number of iterations performed
ipc_ cg_iter
the total number of CG iterations performed
ipc_ f_eval
the total number of evaluations of the objective function
ipc_ g_eval
the total number of evaluations of the gradient of the objective function
ipc_ h_eval
the total number of evaluations of the Hessian of the objective function
ipc_ factorization_max
the maximum number of factorizations in a sub-problem solve
ipc_ factorization_status
the return status from the factorization
int64_t max_entries_factors
the maximum number of entries in the factors
int64_t factorization_integer
the total integer workspace required for the factorization
int64_t factorization_real
the total real workspace required for the factorization
rpc_ factorization_average
the average number of factorizations per sub-problem solve
rpc_ obj
the value of the objective function at the best estimate of the solution determined by the package.
rpc_ norm_g
the norm of the gradient of the objective function at the best estimate of the solution determined by the package.
rpc_ radius
the current value of the trust-region radius
struct tru_time_type time
timings (see above)
struct trs_inform_type trs_inform
inform parameters for TRS
struct gltr_inform_type gltr_inform
inform parameters for GLTR
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 sec_inform_type sec_inform
inform parameters for SEC
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-dimenstional objective function; the code is available in $GALAHAD/src/tru/C/trut.c . A variety of supported Hessian and constraint matrix storage formats are shown.
Notice that C-style indexing is used, and that this is flaggeed by setting
control.f_indexing
to false
. The floating-point type real_wp_
is set in galahad_precision.h
to double
by default, but to float
if the preproccesor variable GALAHAD_SINGLE
is defined.
/* trut.c */
/* Full test for the TRU C interface using C sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_tru.h"
#ifdef REAL_128
#include <quadmath.h>
#endif
// Custom userdata struct
struct userdata_type {
rpc_ p;
};
// Function prototypes
ipc_ fun( ipc_ n, const rpc_ x[], rpc_ *f, const void *);
ipc_ grad( ipc_ n, const rpc_ x[], rpc_ g[], const void *);
ipc_ hess( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[], const void *);
ipc_ hess_dense( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[], const void *);
ipc_ hessprod( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[], bool got_h,
const void *);
ipc_ prec( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[], const void *);
ipc_ fun_diag( ipc_ n, const rpc_ x[], rpc_ *f, const void *);
ipc_ grad_diag( ipc_ n, const rpc_ x[], rpc_ g[], const void *);
ipc_ hess_diag( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[], const void *);
ipc_ hessprod_diag( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
bool got_h, const void *);
int main(void) {
// Derived types
void *data;
struct tru_control_type control;
struct tru_inform_type inform;
// Set user data
struct userdata_type userdata;
userdata.p = 4.0;
// Set problem data
ipc_ n = 3; // dimension
ipc_ ne = 5; // Hesssian elements
ipc_ H_row[] = {0, 1, 2, 2, 2}; // Hessian H
ipc_ H_col[] = {0, 1, 0, 1, 2}; // NB lower triangle
ipc_ H_ptr[] = {0, 1, 2, 5}; // row pointers
// Set storage
rpc_ g[n]; // gradient
char st = ' ';
ipc_ status;
printf(" C sparse matrix indexing\n\n");
printf(" tests options for all-in-one storage format\n\n");
for( ipc_ d=1; d <= 5; d++){
// Initialize TRU
tru_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
//control.print_level = 1;
// Start from 1.5
rpc_ x[] = {1.5,1.5,1.5};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
// control.print_level = 1;
tru_import( &control, &data, &status, n, "coordinate",
ne, H_row, H_col, NULL );
tru_solve_with_mat( &data, &userdata, &status,
n, x, g, ne, fun, grad, hess, prec );
break;
case 2: // sparse by rows
st = 'R';
tru_import( &control, &data, &status, n, "sparse_by_rows",
ne, NULL, H_col, H_ptr);
tru_solve_with_mat( &data, &userdata, &status,
n, x, g, ne, fun, grad, hess, prec );
break;
case 3: // dense
st = 'D';
tru_import( &control, &data, &status, n, "dense",
ne, NULL, NULL, NULL );
tru_solve_with_mat( &data, &userdata, &status,
n, x, g, ne, fun, grad, hess_dense, prec );
break;
case 4: // diagonal
st = 'I';
tru_import( &control, &data, &status, n, "diagonal",
ne, NULL, NULL, NULL );
tru_solve_with_mat( &data, &userdata, &status, n, x, g,
ne, fun_diag, grad_diag, hess_diag, prec );
break;
case 5: // access by products
st = 'P';
tru_import( &control, &data, &status, n, "absent",
ne, NULL, NULL, NULL );
tru_solve_without_mat( &data, &userdata, &status,
n, x, g, fun, grad, hessprod, prec );
break;
}
tru_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: TRU_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
tru_terminate( &data, &control, &inform );
}
printf("\n tests reverse-communication options\n\n");
// reverse-communication input/output
ipc_ eval_status;
rpc_ f = 0.0;
rpc_ u[n], v[n];
rpc_ H_val[ne], H_dense[n*(n+1)/2], H_diag[n];
for( ipc_ d=1; d <= 5; d++){
// Initialize TRU
tru_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = false; // C sparse matrix indexing
// control.print_level = 1;
// control.maxit = 1;
// Start from 1.5
rpc_ x[] = {1.5,1.5,1.5};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
tru_import( &control, &data, &status, n, "coordinate",
ne, H_row, H_col, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, ne, H_val, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess( n, ne, x, H_val, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 2: // sparse by rows
st = 'R';
tru_import( &control, &data, &status, n, "sparse_by_rows", ne,
NULL, H_col, H_ptr );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, ne, H_val, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess( n, ne, x, H_val, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 3: // dense
st = 'D';
tru_import( &control, &data, &status, n, "dense",
ne, NULL, NULL, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, n*(n+1)/2, H_dense, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess_dense( n, n*( n+1)/2, x, H_dense,
&userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 4: // diagonal
st = 'I';
tru_import( &control, &data, &status, n, "diagonal",
ne, NULL, NULL, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, n, H_diag, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun_diag( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad_diag( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess_diag( n, n, x, H_diag, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 5: // access by products
st = 'P';
tru_import( &control, &data, &status, n, "absent",
ne, NULL, NULL, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_without_mat( &data, &status, &eval_status,
n, x, f, g, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 5){ // evaluate H
eval_status = hessprod( n, x, u, v, false, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
}
tru_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: TRU_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
tru_terminate( &data, &control, &inform );
}
}
// Objective function
ipc_ fun( ipc_ n, const rpc_ x[], rpc_ *f, const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
*f = pow(x[0] + x[2] + p, 2) + pow(x[1] + x[2], 2) + cos(x[0]);
return 0;
}
// Gradient of the objective
ipc_ grad( ipc_ n, const rpc_ x[], rpc_ g[], const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
g[0] = 2.0 * ( x[0] + x[2] + p ) - sin(x[0]);
g[1] = 2.0 * ( x[1] + x[2] );
g[2] = 2.0 * ( x[0] + x[2] + p ) + 2.0 * ( x[1] + x[2] );
return 0;
}
// Hessian of the objective
ipc_ hess( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *userdata ){
hval[0] = 2.0 - cos(x[0]);
hval[1] = 2.0;
hval[2] = 2.0;
hval[3] = 2.0;
hval[4] = 4.0;
return 0;
}
// Dense Hessian
ipc_ hess_dense( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *userdata ){
hval[0] = 2.0 - cos(x[0]);
hval[1] = 0.0;
hval[2] = 2.0;
hval[3] = 2.0;
hval[4] = 2.0;
hval[5] = 4.0;
return 0;
}
// Hessian-vector product
ipc_ hessprod( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
bool got_h, const void *userdata ){
u[0] = u[0] + 2.0 * ( v[0] + v[2] ) - cos(x[0]) * v[0];
u[1] = u[1] + 2.0 * ( v[1] + v[2] );
u[2] = u[2] + 2.0 * ( v[0] + v[1] + 2.0 * v[2] );
return 0;
}
// Apply preconditioner
ipc_ prec( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
const void *userdata ){
u[0] = 0.5 * v[0];
u[1] = 0.5 * v[1];
u[2] = 0.25 * v[2];
return 0;
}
// Objective function
ipc_ fun_diag( ipc_ n, const rpc_ x[], rpc_ *f, const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
*f = pow(x[2] + p, 2) + pow(x[1], 2) + cos(x[0]);
return 0;
}
// Gradient of the objective
ipc_ grad_diag( ipc_ n, const rpc_ x[], rpc_ g[], const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
g[0] = -sin(x[0]);
g[1] = 2.0 * x[1];
g[2] = 2.0 * ( x[2] + p );
return 0;
}
// Hessian of the objective
ipc_ hess_diag( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *userdata ){
hval[0] = -cos(x[0]);
hval[1] = 2.0;
hval[2] = 2.0;
return 0;
}
// Hessian-vector product
ipc_ hessprod_diag( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
bool got_h, const void *userdata ){
u[0] = u[0] + - cos(x[0]) * v[0];
u[1] = u[1] + 2.0 * v[1];
u[2] = u[2] + 2.0 * v[2];
return 0;
}
This is the same example, but now fortran-style indexing is used; the code is available in $GALAHAD/src/tru/C/trutf.c .
/* trutf.c */
/* Full test for the TRU C interface using Fortran sparse matrix indexing */
#include <stdio.h>
#include <math.h>
#include "galahad_precision.h"
#include "galahad_cfunctions.h"
#include "galahad_tru.h"
#ifdef REAL_128
#include <quadmath.h>
#endif
// Custom userdata struct
struct userdata_type {
rpc_ p;
};
// Function prototypes
ipc_ fun( ipc_ n, const rpc_ x[], rpc_ *f, const void *);
ipc_ grad( ipc_ n, const rpc_ x[], rpc_ g[], const void *);
ipc_ hess( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[], const void *);
ipc_ hess_dense( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *);
ipc_ hessprod( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
bool got_h, const void *);
ipc_ prec( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
const void *);
ipc_ fun_diag( ipc_ n, const rpc_ x[], rpc_ *f, const void *);
ipc_ grad_diag( ipc_ n, const rpc_ x[], rpc_ g[], const void *);
ipc_ hess_diag( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *);
ipc_ hessprod_diag( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
bool got_h, const void *);
int main(void) {
// Derived types
void *data;
struct tru_control_type control;
struct tru_inform_type inform;
// Set user data
struct userdata_type userdata;
userdata.p = 4.0;
// Set problem data
ipc_ n = 3; // dimension
ipc_ ne = 5; // Hesssian elements
ipc_ H_row[] = {1, 2, 3, 3, 3}; // Hessian H
ipc_ H_col[] = {1, 2, 1, 2, 3}; // NB lower triangle
ipc_ H_ptr[] = {1, 2, 3, 6}; // row pointers
// Set storage
rpc_ g[n]; // gradient
char st = ' ';
ipc_ status;
printf(" Fortran sparse matrix indexing\n\n");
printf(" tests options for all-in-one storage format\n\n");
for( ipc_ d=1; d <= 5; d++){
// Initialize TRU
tru_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = true; // Fortran sparse matrix indexing
//control.print_level = 1;
// Start from 1.5
rpc_ x[] = {1.5,1.5,1.5};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
// control.print_level = 1;
tru_import( &control, &data, &status, n, "coordinate",
ne, H_row, H_col, NULL );
tru_solve_with_mat( &data, &userdata, &status,
n, x, g, ne, fun, grad, hess, prec );
break;
case 2: // sparse by rows
st = 'R';
tru_import( &control, &data, &status, n, "sparse_by_rows",
ne, NULL, H_col, H_ptr);
tru_solve_with_mat( &data, &userdata, &status,
n, x, g, ne, fun, grad, hess, prec );
break;
case 3: // dense
st = 'D';
tru_import( &control, &data, &status, n, "dense",
ne, NULL, NULL, NULL );
tru_solve_with_mat( &data, &userdata, &status,
n, x, g, ne, fun, grad, hess_dense, prec );
break;
case 4: // diagonal
st = 'I';
tru_import( &control, &data, &status, n, "diagonal",
ne, NULL, NULL, NULL );
tru_solve_with_mat( &data, &userdata, &status, n, x, g,
ne, fun_diag, grad_diag, hess_diag, prec );
break;
case 5: // access by products
st = 'P';
tru_import( &control, &data, &status, n, "absent",
ne, NULL, NULL, NULL );
tru_solve_without_mat( &data, &userdata, &status,
n, x, g, fun, grad, hessprod, prec );
break;
}
tru_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: TRU_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
tru_terminate( &data, &control, &inform );
}
printf("\n tests reverse-communication options\n\n");
// reverse-communication input/output
ipc_ eval_status;
rpc_ f = 0.0;
rpc_ u[n], v[n];
rpc_ H_val[ne], H_dense[n*(n+1)/2], H_diag[n];
for( ipc_ d=1; d <= 5; d++){
// Initialize TRU
tru_initialize( &data, &control, &status );
// Set user-defined control options
control.f_indexing = true; // Fortran sparse matrix indexing
// control.print_level = 1;
// control.maxit = 1;
// Start from 1.5
rpc_ x[] = {1.5,1.5,1.5};
switch(d){
case 1: // sparse co-ordinate storage
st = 'C';
tru_import( &control, &data, &status, n, "coordinate",
ne, H_row, H_col, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, ne, H_val, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess( n, ne, x, H_val, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 2: // sparse by rows
st = 'R';
tru_import( &control, &data, &status, n, "sparse_by_rows", ne,
NULL, H_col, H_ptr );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, ne, H_val, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess( n, ne, x, H_val, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 3: // dense
st = 'D';
tru_import( &control, &data, &status, n, "dense",
ne, NULL, NULL, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, n*(n+1)/2, H_dense, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess_dense( n, n*( n+1)/2, x, H_dense,
&userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 4: // diagonal
st = 'I';
tru_import( &control, &data, &status, n, "diagonal",
ne, NULL, NULL, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_with_mat( &data, &status, &eval_status,
n, x, f, g, n, H_diag, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun_diag( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad_diag( n, x, g, &userdata );
}else if(status == 4){ // evaluate H
eval_status = hess_diag( n, n, x, H_diag, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
case 5: // access by products
st = 'P';
tru_import( &control, &data, &status, n, "absent",
ne, NULL, NULL, NULL );
while(true){ // reverse-communication loop
tru_solve_reverse_without_mat( &data, &status, &eval_status,
n, x, f, g, u, v );
if(status == 0){ // successful termination
break;
}else if(status < 0){ // error exit
break;
}else if(status == 2){ // evaluate f
eval_status = fun( n, x, &f, &userdata );
}else if(status == 3){ // evaluate g
eval_status = grad( n, x, g, &userdata );
}else if(status == 5){ // evaluate H
eval_status = hessprod( n, x, u, v, false, &userdata );
}else if(status == 6){ // evaluate the product with P
eval_status = prec( n, x, u, v, &userdata );
}else{
printf(" the value %1" i_ipc_ " of status should not occur\n",
status);
break;
}
}
break;
}
tru_information( &data, &inform, &status );
if(inform.status == 0){
#ifdef REAL_128
// interim replacement for quad output: $GALAHAD/include/galahad_pquad_f.h
#include "galahad_pquad_f.h"
#else
printf("%c:%6" i_ipc_ " iterations. Optimal objective "
"value = %.2f status = %1" i_ipc_ "\n",
st, inform.iter, inform.obj, inform.status);
#endif
}else{
printf("%c: TRU_solve exit status = %1" i_ipc_ "\n",
st, inform.status);
}
//printf("x: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", x[i]);
//printf("\n");
//printf("gradient: ");
//for( ipc_ i = 0; i < n; i++) printf("%f ", g[i]);
//printf("\n");
// Delete internal workspace
tru_terminate( &data, &control, &inform );
}
}
// Objective function
ipc_ fun( ipc_ n, const rpc_ x[], rpc_ *f, const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
*f = pow(x[0] + x[2] + p, 2) + pow(x[1] + x[2], 2) + cos(x[0]);
return 0;
}
// Gradient of the objective
ipc_ grad( ipc_ n, const rpc_ x[], rpc_ g[], const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
g[0] = 2.0 * ( x[0] + x[2] + p ) - sin(x[0]);
g[1] = 2.0 * ( x[1] + x[2] );
g[2] = 2.0 * ( x[0] + x[2] + p ) + 2.0 * ( x[1] + x[2] );
return 0;
}
// Hessian of the objective
ipc_ hess( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *userdata ){
hval[0] = 2.0 - cos(x[0]);
hval[1] = 2.0;
hval[2] = 2.0;
hval[3] = 2.0;
hval[4] = 4.0;
return 0;
}
// Dense Hessian
ipc_ hess_dense( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *userdata ){
hval[0] = 2.0 - cos(x[0]);
hval[1] = 0.0;
hval[2] = 2.0;
hval[3] = 2.0;
hval[4] = 2.0;
hval[5] = 4.0;
return 0;
}
// Hessian-vector product
ipc_ hessprod( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
bool got_h, const void *userdata ){
u[0] = u[0] + 2.0 * ( v[0] + v[2] ) - cos(x[0]) * v[0];
u[1] = u[1] + 2.0 * ( v[1] + v[2] );
u[2] = u[2] + 2.0 * ( v[0] + v[1] + 2.0 * v[2] );
return 0;
}
// Apply preconditioner
ipc_ prec( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
const void *userdata ){
u[0] = 0.5 * v[0];
u[1] = 0.5 * v[1];
u[2] = 0.25 * v[2];
return 0;
}
// Objective function
ipc_ fun_diag( ipc_ n, const rpc_ x[], rpc_ *f, const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
*f = pow(x[2] + p, 2) + pow(x[1], 2) + cos(x[0]);
return 0;
}
// Gradient of the objective
ipc_ grad_diag( ipc_ n, const rpc_ x[], rpc_ g[], const void *userdata ){
struct userdata_type *myuserdata = (struct userdata_type *) userdata;
rpc_ p = myuserdata->p;
g[0] = -sin(x[0]);
g[1] = 2.0 * x[1];
g[2] = 2.0 * ( x[2] + p );
return 0;
}
// Hessian of the objective
ipc_ hess_diag( ipc_ n, ipc_ ne, const rpc_ x[], rpc_ hval[],
const void *userdata ){
hval[0] = -cos(x[0]);
hval[1] = 2.0;
hval[2] = 2.0;
return 0;
}
// Hessian-vector product
ipc_ hessprod_diag( ipc_ n, const rpc_ x[], rpc_ u[], const rpc_ v[],
bool got_h, const void *userdata ){
u[0] = u[0] + - cos(x[0]) * v[0];
u[1] = u[1] + 2.0 * v[1];
u[2] = u[2] + 2.0 * v[2];
return 0;
}