arc_control_type structure#

#include <galahad_arc.h>

struct arc_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;
    ipc_ advanced_start;
    rpc_ stop_g_absolute;
    rpc_ stop_g_relative;
    rpc_ stop_s;
    rpc_ initial_weight;
    rpc_ minimum_weight;
    rpc_ reduce_gap;
    rpc_ tiny_gap;
    rpc_ large_root;
    rpc_ eta_successful;
    rpc_ eta_very_successful;
    rpc_ eta_too_successful;
    rpc_ weight_decrease_min;
    rpc_ weight_decrease;
    rpc_ weight_increase;
    rpc_ weight_increase_max;
    rpc_ obj_unbounded;
    rpc_ cpu_time_limit;
    rpc_ clock_time_limit;
    bool hessian_available;
    bool subproblem_direct;
    bool renormalize_weight;
    bool quadratic_ratio_test;
    bool space_critical;
    bool deallocate_error_fatal;
    char prefix[31];
    struct rqs_control_type rqs_control;
    struct glrt_control_type glrt_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 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 performed

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 (limited-memory BFGS, with .lbfgs_vectors history) (not yet implemented)

  • 5 secant second-order (limited-memory SR1, with .lbfgs_vectors history) (not yet implemented)

ipc_ norm

the regularization 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

ipc_ advanced_start

try to pick a good initial regularization weight using .advanced_start iterates of a variant on the strategy of Sartenaer SISC 18(6) 1990:1788-1803

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

rpc_ initial_weight

Initial value for the regularisation weight (-ve => 1/||g_0||)

rpc_ minimum_weight

minimum permitted regularisation weight

rpc_ reduce_gap

expert parameters as suggested in Gould, Porcelli & Toint, “Updating the regularization parameter in the adaptive cubic regularization algorithm” RAL-TR-2011-007, Rutherford Appleton Laboratory, England (2011), http://epubs.stfc.ac.uk/bitstream/6181/RAL-TR-2011-007.pdf (these are denoted beta, epsilon_chi and alpha_max in the paper)

rpc_ tiny_gap

see reduce_gap

rpc_ large_root

see reduce_gap

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 regularization weight will be decreased if this relative decrease is greater than .eta_very_successful but smaller than .eta_too_successful (the first is eta in Gould, Porcell and Toint, 2011)

rpc_ eta_very_successful

see eta_successful

rpc_ eta_too_successful

see eta_successful

rpc_ weight_decrease_min

on very successful iterations, the regularization weight will be reduced by the factor .weight_decrease but no more than .weight_decrease_min while if the iteration is unsuccessful, the weight will be increased by a factor .weight_increase but no more than .weight_increase_max (these are delta_1, delta_2, delta3 and delta_max in Gould, Porcelli and Toint, 2011)

rpc_ weight_decrease

see weight_decrease_min

rpc_ weight_increase

see weight_decrease_min

rpc_ weight_increase_max

see weight_decrease_min

rpc_ obj_unbounded

the smallest value the onjective 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 renormalize_weight

should the weight be renormalized to account for a change in preconditioner?

bool quadratic_ratio_test

should the test for acceptance involve the quadratic model or the cubic?

bool space_critical

if .space_critical true, every effort will be made to use as little space as possible. This may result in longer computation time

bool deallocate_error_fatal

if .deallocate_error_fatal is true, any array/pointer deallocation error will terminate execution. Otherwise, computation will continue

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 rqs_control_type rqs_control

control parameters for RQS

struct glrt_control_type glrt_control

control parameters for GLRT

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 sha_control_type sha_control

control parameters for SHA