1ac9112b8SAlp Dener /* 2ac9112b8SAlp Dener Context for bound-constrained nonlinear conjugate gradient method 3ac9112b8SAlp Dener */ 4ac9112b8SAlp Dener 5ac9112b8SAlp Dener 6ac9112b8SAlp Dener #ifndef __TAO_BNCG_H 7ac9112b8SAlp Dener #define __TAO_BNCG_H 8ac9112b8SAlp Dener 9ac9112b8SAlp Dener #include <petsc/private/taoimpl.h> 10ac9112b8SAlp Dener 11ac9112b8SAlp Dener typedef struct { 1250b47da0SAdam Denchfield Mat B; 13484c7b14SAdam Denchfield Mat pc; 1450b47da0SAdam Denchfield Vec G_old, X_old, W, work; 15c8bcdf1eSAdam Denchfield Vec g_work, y_work, d_work; 16c8bcdf1eSAdam Denchfield Vec sk, yk; 1789da521bSAlp Dener Vec unprojected_gradient, unprojected_gradient_old; 18ac9112b8SAlp Dener Vec inactive_grad, inactive_step; 1961be54a6SAlp Dener 20484c7b14SAdam Denchfield IS active_lower, active_upper, active_fixed, active_idx, inactive_idx, inactive_old, new_inactives; 21484c7b14SAdam Denchfield 2250b47da0SAdam Denchfield PetscReal alpha; /* convex ratio in the scalar scaling */ 2350b47da0SAdam Denchfield PetscReal hz_theta; 2450b47da0SAdam Denchfield PetscReal xi; /* Parameter for KD, DK, and HZ methods. */ 2550b47da0SAdam Denchfield PetscReal theta; /* The convex combination parameter in the SSML_Broyden method. */ 2650b47da0SAdam Denchfield PetscReal zeta; /* Another parameter, exclusive to Kou-Dai, modifying the y_k scalar contribution */ 2750b47da0SAdam Denchfield PetscReal hz_eta, dk_eta; 2850b47da0SAdam Denchfield PetscReal bfgs_scale, dfp_scale; /* Scaling of the bfgs/dfp tau parameter in the bfgs and broyden methods. Default 1. */ 2950b47da0SAdam Denchfield PetscReal tau_bfgs, tau_dfp; 3050b47da0SAdam Denchfield PetscReal as_step, as_tol, yts, yty, sts; 31484c7b14SAdam Denchfield PetscReal f, delta_min, delta_max; 3250b47da0SAdam Denchfield PetscReal epsilon; /* Machine precision unless changed by user */ 33c8bcdf1eSAdam Denchfield PetscReal eps_23; /* Two-thirds power of machine precision */ 3450b47da0SAdam Denchfield 3550b47da0SAdam Denchfield PetscInt cg_type; /* Formula to use */ 3650b47da0SAdam Denchfield PetscInt min_restart_num; /* Restarts every x*n iterations, where n is the dimension */ 37484c7b14SAdam Denchfield PetscInt ls_fails, resets, descent_error, skipped_updates, pure_gd_steps; 3850b47da0SAdam Denchfield PetscInt iter_quad, min_quad; /* Dynamic restart variables in Dai-Kou, SIAM J. Optim. Vol 23, pp. 296-320, Algorithm 4.1 */ 3950b47da0SAdam Denchfield PetscInt as_type; 4050b47da0SAdam Denchfield 41c0f10754SAlp Dener PetscBool recycle; 42c8bcdf1eSAdam Denchfield PetscBool inv_sig; 43c8bcdf1eSAdam Denchfield PetscReal tol_quad; /* tolerance for Dai-Kou dynamic restart */ 44c8bcdf1eSAdam Denchfield PetscBool dynamic_restart; /* Keeps track of whether or not to do a dynamic (KD) restart */ 45c8bcdf1eSAdam Denchfield PetscBool spaced_restart; /* If true, restarts the CG method every x*n iterations */ 46c8bcdf1eSAdam Denchfield PetscBool use_dynamic_restart; 47c8bcdf1eSAdam Denchfield PetscBool neg_xi; 4850b47da0SAdam Denchfield PetscBool unscaled_restart; /* Gradient descent restarts are done without rescaling*/ 49c8bcdf1eSAdam Denchfield PetscBool diag_scaling; 50484c7b14SAdam Denchfield PetscBool no_scaling; 5150b47da0SAdam Denchfield 52ac9112b8SAlp Dener } TAO_BNCG; 53ac9112b8SAlp Dener 54ac9112b8SAlp Dener #endif /* ifndef __TAO_BNCG_H */ 55ac9112b8SAlp Dener 5661be54a6SAlp Dener PETSC_INTERN PetscErrorCode TaoBNCGEstimateActiveSet(Tao, PetscInt); 57a1318120SAlp Dener PETSC_INTERN PetscErrorCode TaoBNCGBoundStep(Tao, PetscInt, Vec); 58c0f10754SAlp Dener PETSC_EXTERN PetscErrorCode TaoBNCGSetRecycleFlag(Tao, PetscBool); 59c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal, PetscReal, PetscReal, PetscReal*, PetscReal); 60c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao, PetscReal); 61*8ca2df50S PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao, PetscReal, PetscReal, PetscReal, PetscReal, PetscReal, PetscBool); 62c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGComputeDiagScaling(Tao, PetscReal, PetscReal); 63c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGResetUpdate(Tao, PetscReal); 64c8bcdf1eSAdam Denchfield PETSC_INTERN PetscErrorCode TaoBNCGCheckDynamicRestart(Tao, PetscReal, PetscReal, PetscReal, PetscBool*, PetscReal); 65