1 #pragma once 2 3 /* 4 Context for bound-constrained nonlinear conjugate gradient method 5 */ 6 7 #include <petsc/private/taoimpl.h> 8 9 typedef struct { 10 Mat B; 11 Mat pc; 12 Vec G_old, X_old, W, work; 13 Vec g_work, y_work, d_work; 14 Vec sk, yk; 15 Vec unprojected_gradient, unprojected_gradient_old; 16 Vec inactive_grad, inactive_step; 17 18 IS active_lower, active_upper, active_fixed, active_idx, inactive_idx, inactive_old, new_inactives; 19 20 PetscReal alpha; /* convex ratio in the scalar scaling */ 21 PetscReal hz_theta; 22 PetscReal xi; /* Parameter for KD, DK, and HZ methods. */ 23 PetscReal theta; /* The convex combination parameter in the SSML_Broyden method. */ 24 PetscReal zeta; /* Another parameter, exclusive to Kou-Dai, modifying the y_k scalar contribution */ 25 PetscReal hz_eta, dk_eta; 26 PetscReal bfgs_scale, dfp_scale; /* Scaling of the bfgs/dfp tau parameter in the bfgs and broyden methods. Default 1. */ 27 PetscReal tau_bfgs, tau_dfp; 28 PetscReal as_step, as_tol, yts, yty, sts; 29 PetscReal f, delta_min, delta_max; 30 PetscReal epsilon; /* Machine precision unless changed by user */ 31 PetscReal eps_23; /* Two-thirds power of machine precision */ 32 33 TaoBNCGType cg_type; /* Formula to use */ 34 PetscInt min_restart_num; /* Restarts every x*n iterations, where n is the dimension */ 35 PetscInt ls_fails, resets, descent_error, skipped_updates, pure_gd_steps; 36 PetscInt iter_quad, min_quad; /* Dynamic restart variables in Dai-Kou, SIAM J. Optim. Vol 23, pp. 296-320, Algorithm 4.1 */ 37 PetscInt as_type; 38 39 PetscBool inv_sig; 40 PetscReal tol_quad; /* tolerance for Dai-Kou dynamic restart */ 41 PetscBool dynamic_restart; /* Keeps track of whether or not to do a dynamic (KD) restart */ 42 PetscBool spaced_restart; /* If true, restarts the CG method every x*n iterations */ 43 PetscBool use_dynamic_restart; 44 PetscBool neg_xi; 45 PetscBool unscaled_restart; /* Gradient descent restarts are done without rescaling*/ 46 PetscBool diag_scaling; 47 PetscBool no_scaling; 48 49 } TAO_BNCG; 50 51 PETSC_INTERN PetscErrorCode TaoBNCGEstimateActiveSet(Tao, PetscInt); 52 PETSC_INTERN PetscErrorCode TaoBNCGBoundStep(Tao, PetscInt, Vec); 53 PETSC_INTERN PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal, PetscReal, PetscReal, PetscReal *, PetscReal); 54 PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao, PetscReal); 55 PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscBool); 56 PETSC_INTERN PetscErrorCode TaoBNCGComputeDiagScaling(Tao, PetscReal, PetscReal); 57 PETSC_INTERN PetscErrorCode TaoBNCGResetUpdate(Tao, PetscReal); 58 PETSC_INTERN PetscErrorCode TaoBNCGCheckDynamicRestart(Tao, PetscReal, PetscReal, PetscReal, PetscBool *, PetscReal); 59