nls_control_type structure#

#include <galahad_nls.h>

struct nls_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_ jacobian_available;
    ipc_ hessian_available;
    ipc_ model;
    ipc_ norm;
    ipc_ non_monotone;
    ipc_ weight_update_strategy;
    rpc_ stop_c_absolute;
    rpc_ stop_c_relative;
    rpc_ stop_g_absolute;
    rpc_ stop_g_relative;
    rpc_ stop_s;
    rpc_ power;
    rpc_ initial_weight;
    rpc_ minimum_weight;
    rpc_ initial_inner_weight;
    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_ reduce_gap;
    rpc_ tiny_gap;
    rpc_ large_root;
    rpc_ switch_to_newton;
    rpc_ cpu_time_limit;
    rpc_ clock_time_limit;
    bool subproblem_direct;
    bool renormalize_weight;
    bool magic_step;
    bool print_obj;
    bool space_critical;
    bool deallocate_error_fatal;
    char prefix[31];
    struct rqs_control_type rqs_control;
    struct glrt_control_type glrt_control;
    struct psls_control_type psls_control;
    struct bsc_control_type bsc_control;
    struct roots_control_type roots_control;
    struct nls_subproblem_control_type subproblem_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_ jacobian_available

is the Jacobian matrix of first derivatives available (\(\geq\) 2), is access only via matrix-vector products (=1) or is it not available (\(\leq\) 0) ?

ipc_ hessian_available

is the Hessian matrix of second derivatives available (\(\geq\) 2), is access only via matrix-vector products (=1) or is it not available (\(\leq\) 0) ?

ipc_ model

the model used.

Possible values are

  • 0 dynamic (not yet implemented)

  • 1 first-order (no Hessian)

  • 2 barely second-order (identity Hessian)

  • 3 Gauss-Newton (\(J^T J\) Hessian)

  • 4 second-order (exact Hessian)

  • 5 Gauss-Newton to Newton transition

  • 6 tensor Gauss-Newton treated as a least-squares model

  • 7 tensor Gauss-Newton treated as a general model

  • 8 tensor Gauss-Newton transition from a least-squares to a general mode

ipc_ norm

the regularization norm used.

The norm is defined via \(\|v\|^2 = v^T S v\), and will define the preconditioner used for iterative methods. Possible values for \(S\) are

  • -3 user’s own regularization norm

  • -2 \(S\) = limited-memory BFGS matrix (with .PSLS_control.lbfgs_vectors history) (not yet implemented)

  • -1 identity (= Euclidan two-norm)

  • 0 automatic (not yet implemented)

  • 1 diagonal, \(S\) = diag( max(\(J^TJ\) Hessian, .PSLS_control.min_diagonal ) )

  • 2 diagonal, \(S\) = diag( max( Hessian, .PSLS_control.min_diagonal ) )

  • 3 banded, \(S\) = band( Hessian ) with semi-bandwidth .PSLS_control.semi_bandwidth

  • 4 re-ordered band, P=band(order(A)) with semi-bandwidth .PSLS_control.semi_bandwidth

  • 5 full factorization, \(S\) = Hessian, Schnabel-Eskow modification

  • 6 full factorization, \(S\) = Hessian, GMPS modification (not yet implemented)

  • 7 incomplete factorization of Hessian, Lin-More’

  • 8 incomplete factorization of Hessian, HSL_MI28

  • 9 incomplete factorization of Hessian, Munskgaard (not yet implemented)

  • 10 expanding band of Hessian (not yet implemented)

ipc_ non_monotone

non-monotone \(\leq\) 0 monotone strategy used, anything else non-monotone strategy with this history length used

ipc_ weight_update_strategy

define the weight-update strategy: 1 (basic), 2 (reset to zero when very successful), 3 (imitate TR), 4 (increase lower bound), 5 (GPT)

rpc_ stop_c_absolute

overall convergence tolerances. The iteration will terminate when \(||c(x)||_2 \leq\) MAX( .stop_c_absolute, .stop_c_relative \(* \|c(x_{\mbox{initial}})\|_2\) or when the norm of the gradient, \(g = J^T(x) c(x) / \|c(x)\|_2\), of \|\|c(x)\|\|_2 satisfies \(\|g\|_2 \leq\) MAX( .stop_g_absolute, .stop_g_relative \(* \|g_{\mbox{initial}}\|_2\), or if the step is less than .stop_s

rpc_ stop_c_relative

see stop_c_absolute

rpc_ stop_g_absolute

see stop_c_absolute

rpc_ stop_g_relative

see stop_c_absolute

rpc_ stop_s

see stop_c_absolute

rpc_ power

the regularization power (<2 => chosen according to the model)

rpc_ initial_weight

initial value for the regularization weight (-ve => \(1/\|g_0\|)\))

rpc_ minimum_weight

minimum permitted regularization weight

rpc_ initial_inner_weight

initial value for the inner regularization weight for tensor GN (-ve => 0)

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 decreaed 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_ 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 unsucceful, 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_ reduce_gap
expert parameters as suggested in Gould, Porcelli and 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_ switch_to_newton

if the Gauss-Newto to Newton model is specified, switch to Newton as soon as the norm of the gradient g is smaller than switch_to_newton

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 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 scaling?

bool magic_step

allow the user to perform a “magic” step to improve the objective

bool print_obj

print values of the objective/gradient rather than ||c|| and its gradient

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 psls_control_type psls_control

control parameters for PSLS

struct bsc_control_type bsc_control

control parameters for BSC

struct roots_control_type roots_control

control parameters for ROOTS

struct nls_subproblem_control_type subproblem_control

control parameters for the step-finding subproblem