overview of functions provided#
// typedefs typedef float spc_; typedef double rpc_; typedef int ipc_; // structs struct rqs_control_type; struct rqs_history_type; struct rqs_inform_type; struct rqs_time_type; // global functions void rqs_initialize(void **data, struct rqs_control_type* control, ipc_ *status); void rqs_read_specfile(struct rqs_control_type* control, const char specfile[]); void rqs_import( struct rqs_control_type* control, void **data, ipc_ *status, ipc_ n, const char H_type[], ipc_ H_ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[] ); void rqs_import_m( void **data, ipc_ *status, ipc_ n, const char M_type[], ipc_ M_ne, const ipc_ M_row[], const ipc_ M_col[], const ipc_ M_ptr[] ); void rqs_import_a( void **data, ipc_ *status, ipc_ m, const char A_type[], ipc_ A_ne, const ipc_ A_row[], const ipc_ A_col[], const ipc_ A_ptr[] ); void rqs_reset_control( struct rqs_control_type* control, void **data, ipc_ *status ); void rqs_solve_problem( void **data, ipc_ *status, ipc_ n, const rpc_ power, const rpc_ weight, const rpc_ f, const rpc_ c[], ipc_ H_ne, const rpc_ H_val[], rpc_ x[], ipc_ M_ne, const rpc_ M_val[], ipc_ m, ipc_ A_ne, const rpc_ A_val[], rpc_ y[] ); void rqs_information(void **data, struct rqs_inform_type* inform, ipc_ *status); void rqs_terminate( void **data, struct rqs_control_type* control, struct rqs_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 rqs_initialize(void **data, struct rqs_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 rqs_control_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void rqs_read_specfile(struct rqs_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/rqs/RQS.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/rqs.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 rqs_control_type) |
specfile |
is a character string containing the name of the specification file |
void rqs_import( struct rqs_control_type* control, void **data, ipc_ *status, ipc_ n, const char H_type[], ipc_ H_ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[] )
Import problem data into internal storage prior to solution.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see rqs_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 rows (and columns) of H. |
H_type |
is a one-dimensional array of type char that specifies the symmetric storage scheme used for the Hessian, \(H\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, or ‘diagonal’; lower or upper case variants are allowed. |
H_ne |
is a scalar variable of type ipc_, that holds the number of entries in the lower triangular part of \(H\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. |
H_row |
is a one-dimensional array of size H_ne and type ipc_, that holds the row indices of the lower triangular part of \(H\) in the sparse co-ordinate storage scheme. It need not be set for any of the other three schemes, and in this case can be NULL. |
H_col |
is a one-dimensional array of size H_ne and type ipc_, that holds the column indices of the lower triangular part of \(H\) 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 rqs_import_m( void **data, ipc_ *status, ipc_ n, const char M_type[], ipc_ M_ne, const ipc_ M_row[], const ipc_ M_col[], const ipc_ M_ptr[] )
Import data for the scaling matrix M into internal storage prior to solution.
Parameters:
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 rows (and columns) of M. |
M_type |
is a one-dimensional array of type char that specifies the symmetric storage scheme used for the scaling matrix, \(M\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’, or ‘diagonal’; lower or upper case variants are allowed. |
M_ne |
is a scalar variable of type ipc_, that holds the number of entries in the lower triangular part of \(M\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. |
M_row |
is a one-dimensional array of size M_ne and type ipc_, that holds the row indices of the lower triangular part of \(M\) 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. |
M_col |
is a one-dimensional array of size M_ne and type ipc_, that holds the column indices of the lower triangular part of \(M\) in either the sparse co-ordinate, or the sparse row-wise storage scheme. It need not be set when the dense, diagonal or identity storage schemes are used, and in this case can be NULL. |
M_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 \(M\), 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 rqs_import_a( void **data, ipc_ *status, ipc_ m, const char A_type[], ipc_ A_ne, const ipc_ A_row[], const ipc_ A_col[], const ipc_ A_ptr[] )
Import data for the constraint matrix A into internal storage prior to solution.
Parameters:
data |
holds private internal data |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are:
|
m |
is a scalar variable of type ipc_, that holds the number of general linear constraints, i.e., the number of rows of A, if any. m must be non-negative. |
A_type |
is a one-dimensional array of type char that specifies the unsymmetric storage scheme used for the constraint Jacobian, \(A\) if any. It should be one of ‘coordinate’, ‘sparse_by_rows’ or ‘dense’; lower or upper case variants are allowed. |
A_ne |
is a scalar variable of type ipc_, that holds the number of entries in \(A\), if used, in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. |
A_row |
is a one-dimensional array of size A_ne and type ipc_, that holds the row indices of \(A\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes, and in this case can be NULL. |
A_col |
is a one-dimensional array of size A_ne and type ipc_, that holds the column indices of \(A\) 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. |
A_ptr |
is a one-dimensional array of size n+1 and type ipc_, that holds the starting position of each row of \(A\), 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 rqs_reset_control( struct rqs_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 rqs_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 rqs_solve_problem( void **data, ipc_ *status, ipc_ n, const rpc_ power, const rpc_ weight, const rpc_ f, const rpc_ c[], ipc_ H_ne, const rpc_ H_val[], rpc_ x[], ipc_ M_ne, const rpc_ M_val[], ipc_ m, ipc_ A_ne, const rpc_ A_val[], rpc_ y[] )
Solve the regularised quadratic problem.
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:
|
n |
is a scalar variable of type ipc_, that holds the number of variables. |
power |
is a scalar of type rpc_, that holds the order of regularisation, \(p\), used. power must be no smaller than 2. |
weight |
is a scalar of type rpc_, that holds the regularisation weight, \(\sigma\), used. weight must be strictly positive. |
c |
is a one-dimensional array of size n and type rpc_, that holds the linear term \(c\) of the objective function. The j-th component of c, j = 0, … , n-1, contains \(c_j\). |
f |
is a scalar of type rpc_, that holds the constant term \(f\) of the objective function. |
H_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 h_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. |
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\). |
M_ne |
is a scalar variable of type ipc_, that holds the number of entries in the scaling matrix \(M\) if it not the identity matrix. |
M_val |
is a one-dimensional array of size M_ne and type rpc_, that holds the values of the entries of the scaling matrix \(M\), if it is not the identity matrix, in any of the available storage schemes. If M_val is NULL, M will be taken to be the identity matrix. |
m |
is a scalar variable of type ipc_, that holds the number of general linear constraints, if any. m must be non-negative. |
A_ne |
is a scalar variable of type ipc_, that holds the number of entries in the constraint Jacobian matrix \(A\) if used. A_ne must be non-negative. |
A_val |
is a one-dimensional array of size A_ne and type rpc_, that holds the values of the entries of the constraint Jacobian matrix \(A\), if used, in any of the available storage schemes. If A_val is NULL, no constraints will be enforced. |
y |
is a one-dimensional array of size m and type rpc_, that holds the values \(y\) of the Lagrange multipliers for the equality constraints \(A x = 0\) if used. The i-th component of y, i = 0, … , m-1, contains \(y_i\). |
void rqs_information(void **data, struct rqs_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see rqs_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void rqs_terminate( void **data, struct rqs_control_type* control, struct rqs_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see rqs_control_type) |
inform |
is a struct containing output information (see rqs_inform_type) |