1eb910715SAlp Dener /* 2eb910715SAlp Dener Context for bounded Newton-Krylov type optimization algorithms 3eb910715SAlp Dener */ 4eb910715SAlp Dener 5eb910715SAlp Dener #if !defined(__TAO_BNK_H) 6eb910715SAlp Dener #define __TAO_BNK_H 7eb910715SAlp Dener #include <petsc/private/taoimpl.h> 8c0f10754SAlp Dener #include <../src/tao/bound/impls/bncg/bncg.h> 9eb910715SAlp Dener 10eb910715SAlp Dener typedef struct { 11c0f10754SAlp Dener /* Embedded TAOBNCG */ 12c0f10754SAlp Dener Tao bncg; 13c0f10754SAlp Dener TAO_BNCG *bncg_ctx; 14e031d6f5SAlp Dener PetscInt max_cg_its, tot_cg_its; 15c0f10754SAlp Dener Vec bncg_sol; 16c0f10754SAlp Dener 17c0f10754SAlp Dener /* Allocated vectors */ 18c0f10754SAlp Dener Vec W, Xwork, Gwork, Xold, Gold; 1909164190SAlp Dener Vec unprojected_gradient, unprojected_gradient_old; 20c0f10754SAlp Dener Vec Diag_min, Diag_max; 21c0f10754SAlp Dener 22c0f10754SAlp Dener /* Unallocated matrices and vectors */ 23*cd929ea3SAlp Dener Mat H_inactive, Hpre_inactive, M, M_inactive; 24c0f10754SAlp Dener Vec Diag, Diag_red, X_inactive, G_inactive, inactive_work, active_work; 252f75a4aaSAlp Dener IS inactive_idx, active_idx, active_lower, active_upper, active_fixed; 26eb910715SAlp Dener 27080d2917SAlp Dener /* Scalar values for the solution and step */ 28080d2917SAlp Dener PetscReal fold, f, gnorm, dnorm; 29eb910715SAlp Dener 302f75a4aaSAlp Dener /* Parameters for active set estimation */ 310a4511e9SAlp Dener PetscReal as_tol; 320a4511e9SAlp Dener PetscReal as_step; 332f75a4aaSAlp Dener 34eb910715SAlp Dener /* Parameters when updating the perturbation added to the Hessian matrix 35eb910715SAlp Dener according to the following scheme: 36eb910715SAlp Dener 37eb910715SAlp Dener pert = sval; 38eb910715SAlp Dener 39eb910715SAlp Dener do until convergence 40eb910715SAlp Dener shift Hessian by pert 41eb910715SAlp Dener solve Newton system 42eb910715SAlp Dener 43eb910715SAlp Dener if (linear solver failed or did not compute a descent direction) 44eb910715SAlp Dener use steepest descent direction and increase perturbation 45eb910715SAlp Dener 46eb910715SAlp Dener if (0 == pert) 47eb910715SAlp Dener initialize perturbation 48eb910715SAlp Dener pert = min(imax, max(imin, imfac * norm(G))) 49eb910715SAlp Dener else 50eb910715SAlp Dener increase perturbation 51eb910715SAlp Dener pert = min(pmax, max(pgfac * pert, pmgfac * norm(G))) 52eb910715SAlp Dener fi 53eb910715SAlp Dener else 54eb910715SAlp Dener use linear solver direction and decrease perturbation 55eb910715SAlp Dener 56eb910715SAlp Dener pert = min(psfac * pert, pmsfac * norm(G)) 57eb910715SAlp Dener if (pert < pmin) 58eb910715SAlp Dener pert = 0 59eb910715SAlp Dener fi 60eb910715SAlp Dener fi 61eb910715SAlp Dener 62eb910715SAlp Dener perform line search 63eb910715SAlp Dener function and gradient evaluation 64eb910715SAlp Dener check convergence 65eb910715SAlp Dener od 66eb910715SAlp Dener */ 67eb910715SAlp Dener PetscReal sval; /* Starting perturbation value, default zero */ 68eb910715SAlp Dener 69eb910715SAlp Dener PetscReal imin; /* Minimum perturbation added during initialization */ 70eb910715SAlp Dener PetscReal imax; /* Maximum perturbation added during initialization */ 71eb910715SAlp Dener PetscReal imfac; /* Merit function factor during initialization */ 72eb910715SAlp Dener 73eb910715SAlp Dener PetscReal pert; /* Current perturbation value */ 74eb910715SAlp Dener PetscReal pmin; /* Minimim perturbation value */ 75eb910715SAlp Dener PetscReal pmax; /* Maximum perturbation value */ 76eb910715SAlp Dener PetscReal pgfac; /* Perturbation growth factor */ 77eb910715SAlp Dener PetscReal psfac; /* Perturbation shrink factor */ 78eb910715SAlp Dener PetscReal pmgfac; /* Merit function growth factor */ 79eb910715SAlp Dener PetscReal pmsfac; /* Merit function shrink factor */ 80eb910715SAlp Dener 81eb910715SAlp Dener /* Parameters when updating the trust-region radius based on steplength 82eb910715SAlp Dener if step < nu1 (very bad step) 83eb910715SAlp Dener radius = omega1 * min(norm(d), radius) 84eb910715SAlp Dener elif step < nu2 (bad step) 85eb910715SAlp Dener radius = omega2 * min(norm(d), radius) 86eb910715SAlp Dener elif step < nu3 (okay step) 87eb910715SAlp Dener radius = omega3 * radius; 88eb910715SAlp Dener elif step < nu4 (good step) 89eb910715SAlp Dener radius = max(omega4 * norm(d), radius) 90eb910715SAlp Dener else (very good step) 91eb910715SAlp Dener radius = max(omega5 * norm(d), radius) 92eb910715SAlp Dener fi 93eb910715SAlp Dener */ 94eb910715SAlp Dener PetscReal nu1; /* used to compute trust-region radius */ 95eb910715SAlp Dener PetscReal nu2; /* used to compute trust-region radius */ 96eb910715SAlp Dener PetscReal nu3; /* used to compute trust-region radius */ 97eb910715SAlp Dener PetscReal nu4; /* used to compute trust-region radius */ 98eb910715SAlp Dener 99eb910715SAlp Dener PetscReal omega1; /* factor used for trust-region update */ 100eb910715SAlp Dener PetscReal omega2; /* factor used for trust-region update */ 101eb910715SAlp Dener PetscReal omega3; /* factor used for trust-region update */ 102eb910715SAlp Dener PetscReal omega4; /* factor used for trust-region update */ 103eb910715SAlp Dener PetscReal omega5; /* factor used for trust-region update */ 104eb910715SAlp Dener 105eb910715SAlp Dener /* Parameters when updating the trust-region radius based on reduction 106eb910715SAlp Dener 107eb910715SAlp Dener kappa = ared / pred 108eb910715SAlp Dener if kappa < eta1 (very bad step) 109eb910715SAlp Dener radius = alpha1 * min(norm(d), radius) 110eb910715SAlp Dener elif kappa < eta2 (bad step) 111eb910715SAlp Dener radius = alpha2 * min(norm(d), radius) 112eb910715SAlp Dener elif kappa < eta3 (okay step) 113eb910715SAlp Dener radius = alpha3 * radius; 114eb910715SAlp Dener elif kappa < eta4 (good step) 115eb910715SAlp Dener radius = max(alpha4 * norm(d), radius) 116eb910715SAlp Dener else (very good step) 117eb910715SAlp Dener radius = max(alpha5 * norm(d), radius) 118eb910715SAlp Dener fi 119eb910715SAlp Dener */ 120eb910715SAlp Dener PetscReal eta1; /* used to compute trust-region radius */ 121eb910715SAlp Dener PetscReal eta2; /* used to compute trust-region radius */ 122eb910715SAlp Dener PetscReal eta3; /* used to compute trust-region radius */ 123eb910715SAlp Dener PetscReal eta4; /* used to compute trust-region radius */ 124eb910715SAlp Dener 125eb910715SAlp Dener PetscReal alpha1; /* factor used for trust-region update */ 126eb910715SAlp Dener PetscReal alpha2; /* factor used for trust-region update */ 127eb910715SAlp Dener PetscReal alpha3; /* factor used for trust-region update */ 128eb910715SAlp Dener PetscReal alpha4; /* factor used for trust-region update */ 129eb910715SAlp Dener PetscReal alpha5; /* factor used for trust-region update */ 130eb910715SAlp Dener 131eb910715SAlp Dener /* Parameters when updating the trust-region radius based on interpolation 132eb910715SAlp Dener 133eb910715SAlp Dener kappa = ared / pred 134eb910715SAlp Dener if kappa >= 1.0 - mu1 (very good step) 135eb910715SAlp Dener choose tau in [gamma3, gamma4] 136eb910715SAlp Dener radius = max(tau * norm(d), radius) 137eb910715SAlp Dener elif kappa >= 1.0 - mu2 (good step) 138eb910715SAlp Dener choose tau in [gamma2, gamma3] 139eb910715SAlp Dener if (tau >= 1.0) 140eb910715SAlp Dener radius = max(tau * norm(d), radius) 141eb910715SAlp Dener else 142eb910715SAlp Dener radius = tau * min(norm(d), radius) 143eb910715SAlp Dener fi 144eb910715SAlp Dener else (bad step) 145eb910715SAlp Dener choose tau in [gamma1, 1.0] 146eb910715SAlp Dener radius = tau * min(norm(d), radius) 147eb910715SAlp Dener fi 148eb910715SAlp Dener */ 149eb910715SAlp Dener PetscReal mu1; /* used for model agreement in interpolation */ 150eb910715SAlp Dener PetscReal mu2; /* used for model agreement in interpolation */ 151eb910715SAlp Dener 152eb910715SAlp Dener PetscReal gamma1; /* factor used for interpolation */ 153eb910715SAlp Dener PetscReal gamma2; /* factor used for interpolation */ 154eb910715SAlp Dener PetscReal gamma3; /* factor used for interpolation */ 155eb910715SAlp Dener PetscReal gamma4; /* factor used for interpolation */ 156eb910715SAlp Dener 157eb910715SAlp Dener PetscReal theta; /* factor used for interpolation */ 158eb910715SAlp Dener 159eb910715SAlp Dener /* Parameters when initializing trust-region radius based on interpolation */ 160eb910715SAlp Dener PetscReal mu1_i; /* used for model agreement in interpolation */ 161eb910715SAlp Dener PetscReal mu2_i; /* used for model agreement in interpolation */ 162eb910715SAlp Dener 163eb910715SAlp Dener PetscReal gamma1_i; /* factor used for interpolation */ 164eb910715SAlp Dener PetscReal gamma2_i; /* factor used for interpolation */ 165eb910715SAlp Dener PetscReal gamma3_i; /* factor used for interpolation */ 166eb910715SAlp Dener PetscReal gamma4_i; /* factor used for interpolation */ 167eb910715SAlp Dener 168eb910715SAlp Dener PetscReal theta_i; /* factor used for interpolation */ 169eb910715SAlp Dener 170eb910715SAlp Dener /* Other parameters */ 171eb910715SAlp Dener PetscReal min_radius; /* lower bound on initial radius value */ 172eb910715SAlp Dener PetscReal max_radius; /* upper bound on trust region radius */ 173eb910715SAlp Dener PetscReal epsilon; /* tolerance used when computing ared/pred */ 17462675beeSAlp Dener PetscReal dmin, dmax; /* upper and lower bounds for the Hessian diagonal vector */ 175eb910715SAlp Dener 176eb910715SAlp Dener PetscInt newt; /* Newton directions attempted */ 177eb910715SAlp Dener PetscInt bfgs; /* BFGS directions attempted */ 178eb910715SAlp Dener PetscInt sgrad; /* Scaled gradient directions attempted */ 179eb910715SAlp Dener PetscInt grad; /* Gradient directions attempted */ 180eb910715SAlp Dener 18162675beeSAlp Dener PetscInt as_type; /* Active set estimation method */ 182eb910715SAlp Dener PetscInt pc_type; /* Preconditioner for the code */ 183eb910715SAlp Dener PetscInt bfgs_scale_type; /* Scaling matrix to used for the bfgs preconditioner */ 184eb910715SAlp Dener PetscInt init_type; /* Trust-region initialization method */ 185eb910715SAlp Dener PetscInt update_type; /* Trust-region update method */ 186eb910715SAlp Dener 1872f75a4aaSAlp Dener /* Trackers for KSP solution type and convergence reasons */ 188eb910715SAlp Dener PetscInt ksp_atol; 189eb910715SAlp Dener PetscInt ksp_rtol; 190eb910715SAlp Dener PetscInt ksp_ctol; 191eb910715SAlp Dener PetscInt ksp_negc; 192eb910715SAlp Dener PetscInt ksp_dtol; 193eb910715SAlp Dener PetscInt ksp_iter; 194eb910715SAlp Dener PetscInt ksp_othr; 195eb910715SAlp Dener PetscBool is_nash, is_stcg, is_gltr; 196eb910715SAlp Dener } TAO_BNK; 197eb910715SAlp Dener 198eb910715SAlp Dener #endif /* if !defined(__TAO_BNK_H) */ 199eb910715SAlp Dener 200eb910715SAlp Dener #define BNK_NEWTON 0 201eb910715SAlp Dener #define BNK_BFGS 1 202eb910715SAlp Dener #define BNK_SCALED_GRADIENT 2 203eb910715SAlp Dener #define BNK_GRADIENT 3 204eb910715SAlp Dener 205eb910715SAlp Dener #define BNK_PC_NONE 0 206eb910715SAlp Dener #define BNK_PC_AHESS 1 207eb910715SAlp Dener #define BNK_PC_BFGS 2 208eb910715SAlp Dener #define BNK_PC_PETSC 3 209eb910715SAlp Dener #define BNK_PC_TYPES 4 210eb910715SAlp Dener 211eb910715SAlp Dener #define BFGS_SCALE_AHESS 0 212eb910715SAlp Dener #define BFGS_SCALE_PHESS 1 213eb910715SAlp Dener #define BFGS_SCALE_BFGS 2 214eb910715SAlp Dener #define BFGS_SCALE_TYPES 3 215eb910715SAlp Dener 216eb910715SAlp Dener #define BNK_INIT_CONSTANT 0 217eb910715SAlp Dener #define BNK_INIT_DIRECTION 1 218eb910715SAlp Dener #define BNK_INIT_INTERPOLATION 2 219eb910715SAlp Dener #define BNK_INIT_TYPES 3 220eb910715SAlp Dener 221eb910715SAlp Dener #define BNK_UPDATE_STEP 0 222eb910715SAlp Dener #define BNK_UPDATE_REDUCTION 1 223eb910715SAlp Dener #define BNK_UPDATE_INTERPOLATION 2 224eb910715SAlp Dener #define BNK_UPDATE_TYPES 3 225eb910715SAlp Dener 2262f75a4aaSAlp Dener #define BNK_AS_NONE 0 2272f75a4aaSAlp Dener #define BNK_AS_BERTSEKAS 1 2282f75a4aaSAlp Dener #define BNK_AS_TYPES 2 2292f75a4aaSAlp Dener 230eb910715SAlp Dener PETSC_INTERN PetscErrorCode TaoCreate_BNK(Tao); 2319b6ef848SAlp Dener PETSC_INTERN PetscErrorCode TaoSetUp_BNK(Tao); 232eb910715SAlp Dener 233*cd929ea3SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKPreconBFGS(PC, Vec, Vec); 234c0f10754SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKInitialize(Tao, PetscInt, PetscBool*); 23508752603SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKEstimateActiveSet(Tao, PetscInt); 23662675beeSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKComputeHessian(Tao); 237a1318120SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKBoundStep(Tao, PetscInt, Vec); 238c0f10754SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKTakeCGSteps(Tao, PetscBool*); 23962675beeSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKComputeStep(Tao, PetscBool, KSPConvergedReason*); 2405e9b73cbSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKRecomputePred(Tao, Vec, PetscReal*); 241e465cd6fSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKSafeguardStep(Tao, KSPConvergedReason, PetscInt*); 242937a31a1SAlp Dener PETSC_INTERN PetscErrorCode TaoBNKPerformLineSearch(Tao, PetscInt*, PetscReal*, TaoLineSearchConvergedReason*); 24328017e9fSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKUpdateTrustRadius(Tao, PetscReal, PetscReal, PetscInt, PetscInt, PetscBool*); 24462675beeSAlp Dener PETSC_INTERN PetscErrorCode TaoBNKAddStepCounts(Tao, PetscInt); 245