overview of functions provided#
// typedefs typedef float spc_; typedef double rpc_; typedef int ipc_; // structs struct nls_subproblem_control_type; struct nls_control_type; struct nls_subproblem_inform_type; struct nls_inform_type; struct nls_time_type; // function calls void nls_initialize( void **data, struct nls_control_type* control, struct nls_inform_type* inform ); void nls_read_specfile(struct nls_control_type* control, const char specfile[]); void nls_import( struct nls_control_type* control, void **data, ipc_ *status, ipc_ n, ipc_ m, const char J_type[], ipc_ J_ne, const ipc_ J_row[], const ipc_ J_col[], const ipc_ J_ptr[], const char H_type[], ipc_ H_ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[], const char P_type[], ipc_ P_ne, const ipc_ P_row[], const ipc_ P_col[], const ipc_ P_ptr[], const rpc_ w[] ); void nls_reset_control( struct nls_control_type* control, void **data, ipc_ *status ); void nls_solve_with_mat( void **data, void *userdata, ipc_ *status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], ipc_(*)(ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_c, ipc_ j_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_j, ipc_ h_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], const void*) eval_h, ipc_ p_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], bool, const void*) eval_hprods ); void nls_solve_without_mat( void **data, void *userdata, ipc_ *status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], ipc_(*)(ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_c, ipc_(*)(ipc_, ipc_, const rpc_[], const bool, rpc_[], const rpc_[], bool, const void*) eval_jprod, ipc_(*)(ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], const rpc_[], bool, const void*) eval_hprod, ipc_ p_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], bool, const void*) eval_hprods ); void nls_solve_reverse_with_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], ipc_ j_ne, rpc_ J_val[], const rpc_ y[], ipc_ h_ne, rpc_ H_val[], rpc_ v[], ipc_ p_ne, rpc_ P_val[] ); void nls_solve_reverse_without_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], bool* transpose, rpc_ u[], rpc_ v[], rpc_ y[], ipc_ p_ne, rpc_ P_val[] ); void nls_information(void **data, struct nls_inform_type* inform, ipc_ *status); void nls_terminate( void **data, struct nls_control_type* control, struct nls_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 nls_initialize( void **data, struct nls_control_type* control, struct nls_inform_type* inform )
Set default control values and initialize private data
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see nls_control_type) |
inform |
is a struct containing output information (see nls_inform_type) |
void nls_read_specfile(struct nls_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/nls/NLS.template. See also Table 2.1 in the Fortran documentation provided in $GALAHAD/doc/nls.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 nls_control_type) |
specfile |
is a character string containing the name of the specification file |
void nls_import( struct nls_control_type* control, void **data, ipc_ *status, ipc_ n, ipc_ m, const char J_type[], ipc_ J_ne, const ipc_ J_row[], const ipc_ J_col[], const ipc_ J_ptr[], const char H_type[], ipc_ H_ne, const ipc_ H_row[], const ipc_ H_col[], const ipc_ H_ptr[], const char P_type[], ipc_ P_ne, const ipc_ P_row[], const ipc_ P_col[], const ipc_ P_ptr[], const rpc_ w[] )
Import problem data into internal storage prior to solution.
Parameters:
control |
is a struct whose members provide control paramters for the remaining prcedures (see nls_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. |
m |
is a scalar variable of type ipc_, that holds the number of residuals. |
J_type |
is a one-dimensional array of type char that specifies the unsymmetric storage scheme used for the Jacobian, \(J\). It should be one of ‘coordinate’, ‘sparse_by_rows’, ‘dense’ or ‘absent’, the latter if access to the Jacobian is via matrix-vector products; lower or upper case variants are allowed. |
J_ne |
is a scalar variable of type ipc_, that holds the number of entries in \(J\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. |
J_row |
is a one-dimensional array of size J_ne and type ipc_, that holds the row indices of \(J\) 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. |
J_col |
is a one-dimensional array of size J_ne and type ipc_, that holds the column indices of \(J\) 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. |
J_ptr |
is a one-dimensional array of size m+1 and type ipc_, that holds the starting position of each row of \(J\), 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. |
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’, ‘diagonal’ or ‘absent’, the latter if access to \(H\) is via matrix-vector products; 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 three 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. |
P_type |
is a one-dimensional array of type char that specifies the unsymmetric storage scheme used for the residual-Hessians-vector product matrix, \(P\). It should be one of ‘coordinate’, ‘sparse_by_columns’, ‘dense_by_columns’ or ‘absent’, the latter if access to \(P\) is via matrix-vector products; lower or upper case variants are allowed. |
P_ne |
is a scalar variable of type ipc_, that holds the number of entries in \(P\) in the sparse co-ordinate storage scheme. It need not be set for any of the other schemes. |
P_row |
is a one-dimensional array of size P_ne and type ipc_, that holds the row indices of \(P\) in either the sparse co-ordinate, or the sparse column-wise storage scheme. It need not be set when the dense storage scheme is used, and in this case can be NULL. |
P_col |
is a one-dimensional array of size P_ne and type ipc_, that holds the row indices of \(P\) 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. |
P_ptr |
is a one-dimensional array of size n+1 and type ipc_, that holds the starting position of each row of \(P\), 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. |
w |
is a one-dimensional array of size m and type rpc_, that holds the values \(w\) of the weights on the residuals in the least-squares objective function. It need not be set if the weights are all ones, and in this case can be NULL. |
void nls_reset_control( struct nls_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 nls_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 nls_solve_with_mat( void **data, void *userdata, ipc_ *status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], ipc_(*)(ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_c, ipc_ j_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_j, ipc_ h_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], const void*) eval_h, ipc_ p_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], bool, const void*) eval_hprods )
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. |
m |
is a scalar variable of type ipc_, that holds the number of residuals. |
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\). |
c |
is a one-dimensional array of size m and type rpc_, that holds the residual \(c(x)\). The i-th component of c, j = 0, … , n-1, contains \(c_j(x)\). |
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_c |
is a user-supplied function that must have the following signature: ipc_ eval_c( ipc_ n, const rpc_ x[], rpc_ c[], const void *userdata ) The componnts of the residual function \(c(x)\) evaluated at x= \(x\) must be assigned to c, 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 |
j_ne |
is a scalar variable of type ipc_, that holds the number of entries in the Jacobian matrix \(J\). |
eval_j |
is a user-supplied function that must have the following signature: ipc_ eval_j( ipc_ n, ipc_ m, ipc_ jne, const rpc_ x[], rpc_ j[], const void *userdata ) The components of the Jacobian \(J = \nabla_x c(x\)) of the residuals must be assigned to j in the same order as presented to nls_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 |
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\) if it is used. |
eval_h |
is a user-supplied function that must have the following signature: ipc_ eval_h( ipc_ n, ipc_ m, ipc_ hne, const rpc_ x[], const rpc_ y[], rpc_ h[], const void *userdata ) The nonzeros of the matrix \(H = \sum_{i=1}^m y_i \nabla_{xx}c_i(x)\) of the weighted residual Hessian evaluated at x= \(x\) and y= \(y\) must be assigned to h in the same order as presented to nls_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 |
p_ne |
is a scalar variable of type ipc_, that holds the number of entries in the residual-Hessians-vector product matrix \(P\) if it is used. |
eval_hprods |
is an optional user-supplied function that may be NULL. If non-NULL, it must have the following signature: ipc_ eval_hprods( ipc_ n, ipc_ m, ipc_ pne, const rpc_ x[], const rpc_ v[], rpc_ p[], bool got_h, const void *userdata ) ); The entries of the matrix \(P\), whose i-th column is the product \(\nabla_{xx}c_i(x) v\) between \(\nabla_{xx}c_i(x)\), the Hessian of the i-th component of the residual \(c(x)\) at x= \(x\), and v= \(v\) must be returned in p 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 nls_solve_without_mat( void **data, void *userdata, ipc_ *status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], ipc_(*)(ipc_, ipc_, const rpc_[], rpc_[], const void*) eval_c, ipc_(*)(ipc_, ipc_, const rpc_[], const bool, rpc_[], const rpc_[], bool, const void*) eval_jprod, ipc_(*)(ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], const rpc_[], bool, const void*) eval_hprod, ipc_ p_ne, ipc_(*)(ipc_, ipc_, ipc_, const rpc_[], const rpc_[], rpc_[], bool, const void*) eval_hprods )
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 |
m |
is a scalar variable of type ipc_, that holds the number of residuals. |
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\). |
c |
is a one-dimensional array of size m and type rpc_, that holds the residual \(c(x)\). The i-th component of c, j = 0, … , n-1, contains \(c_j(x)\). |
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_c |
is a user-supplied function that must have the following signature: ipc_ eval_c( ipc_ n, const rpc_ x[], rpc_ c[], const void *userdata ) The componnts of the residual function \(c(x)\) evaluated at x= \(x\) must be assigned to c, 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_jprod |
is a user-supplied function that must have the following signature: ipc_ eval_jprod( ipc_ n, ipc_ m, const rpc_ x[], bool transpose, rpc_ u[], const rpc_ v[], bool got_j, const void *userdata ) The sum \(u + \nabla_{x}c_(x) v\) (if tranpose is false) or The sum \(u + (\nabla_{x}c_(x))^T v\) (if tranpose is true) bewteen the product of the Jacobian \(\nabla_{x}c_(x)\) or its tranpose 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. Data may be passed into |
eval_hprod |
is a user-supplied function that must have the following signature: ipc_ eval_hprod( ipc_ n, ipc_ m, const rpc_ x[], const rpc_ y[], rpc_ u[], const rpc_ v[], bool got_h, const void *userdata ) The sum \(u + \sum_{i=1}^m y_i \nabla_{xx}c_i(x) v\) of the product of the weighted residual Hessian \(H = \sum_{i=1}^m y_i \nabla_{xx}c_i(x)\) evaluated at x= \(x\) and y= \(y\) 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 Hessians have already been evaluated or used at x if got_h is true. Data may be passed into |
p_ne |
is a scalar variable of type ipc_, that holds the number of entries in the residual-Hessians-vector product matrix \(P\) if it is used. |
eval_hprods |
is an optional user-supplied function that may be NULL. If non-NULL, it must have the following signature: ipc_ eval_hprods( ipc_ n, ipc_ m, ipc_ p_ne, const rpc_ x[], const rpc_ v[], rpc_ pval[], bool got_h, const void *userdata ) The entries of the matrix \(P\), whose i-th column is the product \(\nabla_{xx}c_i(x) v\) between \(\nabla_{xx}c_i(x)\), the Hessian of the i-th component of the residual \(c(x)\) at x= \(x\), and v= \(v\) must be returned in pval 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 nls_solve_reverse_with_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], ipc_ j_ne, rpc_ J_val[], const rpc_ y[], ipc_ h_ne, rpc_ H_val[], rpc_ v[], ipc_ p_ne, rpc_ P_val[] )
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 |
m |
is a scalar variable of type ipc_, that holds the number of residuals. |
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\). |
c |
is a one-dimensional array of size m and type rpc_, that holds the residual \(c(x)\). The i-th component of c, j = 0, … , n-1, contains \(c_j(x)\). See status = 2, above, for more details. |
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\). |
j_ne |
is a scalar variable of type ipc_, that holds the number of entries in the Jacobian matrix \(J\). |
J_val |
is a one-dimensional array of size j_ne and type rpc_, that holds the values of the entries of the Jacobian matrix \(J\) in any of the available storage schemes. See status = 3, above, for more details. |
y |
is a one-dimensional array of size m and type rpc_, that is used for reverse communication. See status = 4 above for more details. |
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. See status = 4, above, for more details. |
v |
is a one-dimensional array of size n and type rpc_, that is used for reverse communication. See status = 7, above, for more details. |
p_ne |
is a scalar variable of type ipc_, that holds the number of entries in the residual-Hessians-vector product matrix, \(P\). |
P_val |
is a one-dimensional array of size p_ne and type rpc_, that holds the values of the entries of the residual-Hessians-vector product matrix, \(P\). See status = 7, above, for more details. |
void nls_solve_reverse_without_mat( void **data, ipc_ *status, ipc_ *eval_status, ipc_ n, ipc_ m, rpc_ x[], rpc_ c[], rpc_ g[], bool* transpose, rpc_ u[], rpc_ v[], rpc_ y[], ipc_ p_ne, rpc_ P_val[] )
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 |
m |
is a scalar variable of type ipc_, that holds the number of residuals. |
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\). |
c |
is a one-dimensional array of size m and type rpc_, that holds the residual \(c(x)\). The i-th component of c, j = 0, … , n-1, contains \(c_j(x)\). See status = 2, above, for more details. |
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\). |
transpose |
is a scalar variable of type bool, that indicates whether the product with Jacobian or its transpose should be obtained when status=5. |
u |
is a one-dimensional array of size max(n,m) and type rpc_, that is used for reverse communication. See status = 5,6 above for more details. |
v |
is a one-dimensional array of size max(n,m) and type rpc_, that is used for reverse communication. See status = 5,6,7 above for more details. |
y |
is a one-dimensional array of size m and type rpc_, that is used for reverse communication. See status = 6 above for more details. |
p_ne |
is a scalar variable of type ipc_, that holds the number of entries in the residual-Hessians-vector product matrix, \(P\). |
P_val |
is a one-dimensional array of size P_ne and type rpc_, that holds the values of the entries of the residual-Hessians-vector product matrix, \(P\). See status = 7, above, for more details. |
void nls_information(void **data, struct nls_inform_type* inform, ipc_ *status)
Provides output information
Parameters:
data |
holds private internal data |
inform |
is a struct containing output information (see nls_inform_type) |
status |
is a scalar variable of type ipc_, that gives the exit status from the package. Possible values are (currently):
|
void nls_terminate( void **data, struct nls_control_type* control, struct nls_inform_type* inform )
Deallocate all internal private storage
Parameters:
data |
holds private internal data |
control |
is a struct containing control information (see nls_control_type) |
inform |
is a struct containing output information (see nls_inform_type) |