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