overview of functions provided#
// typedefs typedef float spc_; typedef double rpc_; typedef int ipc_; // structs struct arc_control_type; struct arc_inform_type; struct arc_time_type; // function calls void arc_initialize(void **data, struct arc_control_type* control, ipc_ *status); void arc_read_specfile(struct arc_control_type* control, const char specfile[]); void arc_import( struct arc_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 arc_reset_control( struct arc_control_type* control, void **data, ipc_ *status ); void arc_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 arc_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 arc_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 arc_solve_reverse_without_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, rpc_ x[], rpc_ f, rpc_ g[], rpc_ u[], rpc_ v[] ); void arc_information(void **data, struct arc_inform_type* inform, ipc_ *status); void arc_terminate( void **data, struct arc_control_type* control, struct arc_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 SINGLE
.
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 arc_initialize(void **data, struct arc_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 arc_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void arc_read_specfile(struct arc_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/arc/ARC.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/arc.pdf for a list of how these keywords relate to the components of the control structure.
Parameters:
control |
is a struct containing control information (see arc_control_type) |
specfile |
is a character string containing the name of the specification file |
void arc_import( struct arc_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 arc_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 arc_reset_control( struct arc_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 arc_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 arc_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 regularization method.
This call is for the case where \(H = \nabla_{xx}f(x)\) is provided specifically, and all function/derivative information is available by function calls.
Parameters:
data |
holds private internal data |
userdata |
is a structure that allows data to be passed into the function and derivative evaluation programs. |
status |
is a scalar variable of type 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 arc_import, and the function return value set to 0. If the evaluation is impossible at x, return should be set to a nonzero value. Data may be passed into |
eval_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 arc_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 regularization method.
This call is for the case where access to \(H = \nabla_{xx}f(x)\) is provided by Hessian-vector products, and all function/derivative information is available by function calls.
Parameters:
data |
holds private internal data |
userdata |
is a structure that allows data to be passed into the function and derivative evaluation programs. |
status |
is a scalar variable of type 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 arc_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 regularization method.
This call is for the case where \(H = \nabla_{xx}f(x)\) is provided specifically, but function/derivative information is only available by returning to the calling procedure
Parameters:
data |
holds private internal data |
status |
is a scalar variable of type 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 arc_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 regularization method.
This call is for the case where access to \(H = \nabla_{xx}f(x)\) is provided by Hessian-vector products, but function/derivative information is only available by returning to the calling procedure.
Parameters:
data |
holds private internal data |
status |
is a scalar variable of type 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 arc_information(void **data, struct arc_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see arc_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void arc_terminate( void **data, struct arc_control_type* control, struct arc_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see arc_control_type) |
inform |
is a struct containing output information (see arc_inform_type) |