1 #include <petsctaolinesearch.h> 2 #include <../src/tao/bound/impls/bncg/bncg.h> 3 #include <petscksp.h> 4 5 #define CG_GradientDescent 0 6 #define CG_HestenesStiefel 1 7 #define CG_FletcherReeves 2 8 #define CG_PolakRibierePolyak 3 9 #define CG_PolakRibierePlus 4 10 #define CG_DaiYuan 5 11 #define CG_HagerZhang 6 12 #define CG_DaiKou 7 13 #define CG_KouDai 8 14 #define CG_SSML_BFGS 9 15 #define CG_SSML_DFP 10 16 #define CG_SSML_BROYDEN 11 17 #define CG_PCGradientDescent 12 18 #define CGTypes 13 19 20 static const char *CG_Table[64] = {"gd", "hs", "fr", "pr", "prp", "dy", "hz", "dk", "kd", "ssml_bfgs", "ssml_dfp", "ssml_brdn", "pcgd"}; 21 22 #define CG_AS_NONE 0 23 #define CG_AS_BERTSEKAS 1 24 #define CG_AS_SIZE 2 25 26 static const char *CG_AS_TYPE[64] = {"none", "bertsekas"}; 27 28 PetscErrorCode TaoBNCGSetRecycleFlag(Tao tao, PetscBool recycle) 29 { 30 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 31 32 PetscFunctionBegin; 33 cg->recycle = recycle; 34 PetscFunctionReturn(0); 35 } 36 37 PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType) 38 { 39 PetscErrorCode ierr; 40 TAO_BNCG *cg = (TAO_BNCG *)tao->data; 41 42 PetscFunctionBegin; 43 ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 44 if (cg->inactive_idx) { 45 ierr = ISDuplicate(cg->inactive_idx, &cg->inactive_old);CHKERRQ(ierr); 46 ierr = ISCopy(cg->inactive_idx, cg->inactive_old);CHKERRQ(ierr); 47 } 48 switch (asType) { 49 case CG_AS_NONE: 50 ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 51 ierr = VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx);CHKERRQ(ierr); 52 ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 53 ierr = ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx);CHKERRQ(ierr); 54 break; 55 56 case CG_AS_BERTSEKAS: 57 /* Use gradient descent to estimate the active set */ 58 ierr = VecCopy(cg->unprojected_gradient, cg->W);CHKERRQ(ierr); 59 ierr = VecScale(cg->W, -1.0);CHKERRQ(ierr); 60 ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->work, cg->as_step, &cg->as_tol, 61 &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx);CHKERRQ(ierr); 62 break; 63 64 default: 65 break; 66 } 67 PetscFunctionReturn(0); 68 } 69 70 PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step) 71 { 72 PetscErrorCode ierr; 73 TAO_BNCG *cg = (TAO_BNCG *)tao->data; 74 75 PetscFunctionBegin; 76 switch (asType) { 77 case CG_AS_NONE: 78 ierr = VecISSet(step, cg->active_idx, 0.0);CHKERRQ(ierr); 79 break; 80 81 case CG_AS_BERTSEKAS: 82 ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step);CHKERRQ(ierr); 83 break; 84 85 default: 86 break; 87 } 88 PetscFunctionReturn(0); 89 } 90 91 static PetscErrorCode TaoSolve_BNCG(Tao tao) 92 { 93 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 94 PetscErrorCode ierr; 95 PetscReal step=1.0,gnorm,gnorm2, resnorm; 96 PetscInt nDiff; 97 98 PetscFunctionBegin; 99 /* Project the current point onto the feasible set */ 100 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 101 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 102 103 /* Project the initial point onto the feasible region */ 104 ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr); 105 106 if (nDiff > 0 || !cg->recycle){ 107 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient);CHKERRQ(ierr); 108 } 109 ierr = VecNorm(cg->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr); 110 if (PetscIsInfOrNanReal(cg->f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 111 112 /* Estimate the active set and compute the projected gradient */ 113 ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 114 115 /* Project the gradient and calculate the norm */ 116 ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 117 ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 118 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 119 gnorm2 = gnorm*gnorm; 120 121 /* Initialize counters */ 122 tao->niter = 0; 123 cg->ls_fails = cg->descent_error = 0; 124 cg->resets = -1; 125 cg->skipped_updates = cg->pure_gd_steps = 0; 126 cg->iter_quad = 0; 127 128 /* Convergence test at the starting point. */ 129 tao->reason = TAO_CONTINUE_ITERATING; 130 131 ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 132 ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 133 if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 134 ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 135 ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 136 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 137 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 138 /* Calculate initial direction. */ 139 if (!cg->recycle) { 140 /* We are not recycling a solution/history from a past TaoSolve */ 141 ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 142 } 143 /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */ 144 while(1) { 145 ierr = TaoBNCGConductIteration(tao, gnorm);CHKERRQ(ierr); 146 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 147 } 148 PetscFunctionReturn(0); 149 } 150 151 static PetscErrorCode TaoSetUp_BNCG(Tao tao) 152 { 153 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 if (!tao->gradient) { 158 ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 159 } 160 if (!tao->stepdirection) { 161 ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 162 } 163 if (!cg->W) { 164 ierr = VecDuplicate(tao->solution,&cg->W);CHKERRQ(ierr); 165 } 166 if (!cg->work) { 167 ierr = VecDuplicate(tao->solution,&cg->work);CHKERRQ(ierr); 168 } 169 if (!cg->sk) { 170 ierr = VecDuplicate(tao->solution,&cg->sk);CHKERRQ(ierr); 171 } 172 if (!cg->yk) { 173 ierr = VecDuplicate(tao->gradient,&cg->yk);CHKERRQ(ierr); 174 } 175 if (!cg->X_old) { 176 ierr = VecDuplicate(tao->solution,&cg->X_old);CHKERRQ(ierr); 177 } 178 if (!cg->G_old) { 179 ierr = VecDuplicate(tao->gradient,&cg->G_old);CHKERRQ(ierr); 180 } 181 if (cg->diag_scaling){ 182 ierr = VecDuplicate(tao->solution,&cg->d_work);CHKERRQ(ierr); 183 ierr = VecDuplicate(tao->solution,&cg->y_work);CHKERRQ(ierr); 184 ierr = VecDuplicate(tao->solution,&cg->g_work);CHKERRQ(ierr); 185 } 186 if (!cg->unprojected_gradient) { 187 ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient);CHKERRQ(ierr); 188 } 189 if (!cg->unprojected_gradient_old) { 190 ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient_old);CHKERRQ(ierr); 191 } 192 ierr = MatLMVMAllocate(cg->B, cg->sk, cg->yk);CHKERRQ(ierr); 193 if (cg->pc){ 194 ierr = MatLMVMSetJ0(cg->B, cg->pc);CHKERRQ(ierr); 195 } 196 PetscFunctionReturn(0); 197 } 198 199 static PetscErrorCode TaoDestroy_BNCG(Tao tao) 200 { 201 TAO_BNCG *cg = (TAO_BNCG*) tao->data; 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 if (tao->setupcalled) { 206 ierr = VecDestroy(&cg->W);CHKERRQ(ierr); 207 ierr = VecDestroy(&cg->work);CHKERRQ(ierr); 208 ierr = VecDestroy(&cg->X_old);CHKERRQ(ierr); 209 ierr = VecDestroy(&cg->G_old);CHKERRQ(ierr); 210 ierr = VecDestroy(&cg->unprojected_gradient);CHKERRQ(ierr); 211 ierr = VecDestroy(&cg->unprojected_gradient_old);CHKERRQ(ierr); 212 ierr = VecDestroy(&cg->g_work);CHKERRQ(ierr); 213 ierr = VecDestroy(&cg->d_work);CHKERRQ(ierr); 214 ierr = VecDestroy(&cg->y_work);CHKERRQ(ierr); 215 ierr = VecDestroy(&cg->sk);CHKERRQ(ierr); 216 ierr = VecDestroy(&cg->yk);CHKERRQ(ierr); 217 } 218 ierr = ISDestroy(&cg->active_lower);CHKERRQ(ierr); 219 ierr = ISDestroy(&cg->active_upper);CHKERRQ(ierr); 220 ierr = ISDestroy(&cg->active_fixed);CHKERRQ(ierr); 221 ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 222 ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 223 ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 224 ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 225 ierr = MatDestroy(&cg->B);CHKERRQ(ierr); 226 if (cg->pc) { 227 ierr = MatDestroy(&cg->pc);CHKERRQ(ierr); 228 } 229 ierr = PetscFree(tao->data);CHKERRQ(ierr); 230 PetscFunctionReturn(0); 231 } 232 233 static PetscErrorCode TaoSetFromOptions_BNCG(PetscOptionItems *PetscOptionsObject,Tao tao) 234 { 235 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 236 PetscErrorCode ierr; 237 238 PetscFunctionBegin; 239 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 240 ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");CHKERRQ(ierr); 241 ierr = PetscOptionsEList("-tao_bncg_type","cg formula", "", CG_Table, CGTypes, CG_Table[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr); 242 if (cg->cg_type != CG_SSML_BFGS){ 243 cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */ 244 } 245 if (CG_GradientDescent == cg->cg_type){ 246 cg->cg_type = CG_PCGradientDescent; 247 /* Set scaling equal to none or, at best, scalar scaling. */ 248 cg->unscaled_restart = PETSC_TRUE; 249 cg->diag_scaling = PETSC_FALSE; 250 } 251 ierr = PetscOptionsEList("-tao_bncg_as_type","active set estimation method", "", CG_AS_TYPE, CG_AS_SIZE, CG_AS_TYPE[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr); 252 253 ierr = PetscOptionsReal("-tao_bncg_hz_eta","(developer) cutoff tolerance for HZ", "", cg->hz_eta,&cg->hz_eta,NULL);CHKERRQ(ierr); 254 ierr = PetscOptionsReal("-tao_bncg_eps","(developer) cutoff value for restarts", "", cg->epsilon,&cg->epsilon,NULL);CHKERRQ(ierr); 255 ierr = PetscOptionsReal("-tao_bncg_dk_eta","(developer) cutoff tolerance for DK", "", cg->dk_eta,&cg->dk_eta,NULL);CHKERRQ(ierr); 256 ierr = PetscOptionsReal("-tao_bncg_xi","(developer) Parameter in the KD method", "", cg->xi,&cg->xi,NULL);CHKERRQ(ierr); 257 ierr = PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL);CHKERRQ(ierr); 258 ierr = PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL);CHKERRQ(ierr); 259 ierr = PetscOptionsReal("-tao_bncg_alpha","(developer) parameter for the scalar scaling","",cg->alpha,&cg->alpha,NULL);CHKERRQ(ierr); 260 ierr = PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL);CHKERRQ(ierr); 261 ierr = PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL);CHKERRQ(ierr); 262 ierr = PetscOptionsBool("-tao_bncg_diag_scaling","Enable diagonal Broyden-like preconditioning","",cg->diag_scaling,&cg->diag_scaling,NULL);CHKERRQ(ierr); 263 ierr = PetscOptionsBool("-tao_bncg_dynamic_restart","(developer) use dynamic restarts as in HZ, DK, KD","",cg->use_dynamic_restart,&cg->use_dynamic_restart,NULL);CHKERRQ(ierr); 264 ierr = PetscOptionsBool("-tao_bncg_unscaled_restart","(developer) use unscaled gradient restarts","",cg->unscaled_restart,&cg->unscaled_restart,NULL);CHKERRQ(ierr); 265 ierr = PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL);CHKERRQ(ierr); 266 ierr = PetscOptionsInt("-tao_bncg_min_quad", "(developer) Number of iterations with approximate quadratic behavior needed for restart", "", cg->min_quad, &cg->min_quad, NULL);CHKERRQ(ierr); 267 ierr = PetscOptionsInt("-tao_bncg_min_restart_num", "(developer) Number of iterations between restarts (times dimension)", "", cg->min_restart_num, &cg->min_restart_num, NULL);CHKERRQ(ierr); 268 ierr = PetscOptionsBool("-tao_bncg_recycle","enable recycling the existing solution, gradient, and diagonal scaling vector at the start of a new TaoSolve() call","",cg->recycle,&cg->recycle,NULL);CHKERRQ(ierr); 269 ierr = PetscOptionsBool("-tao_bncg_spaced_restart","(developer) Enable regular steepest descent restarting every fixed number of iterations","",cg->spaced_restart,&cg->spaced_restart,NULL);CHKERRQ(ierr); 270 ierr = PetscOptionsBool("-tao_bncg_no_scaling","Disable all scaling except in restarts","",cg->no_scaling,&cg->no_scaling,NULL);CHKERRQ(ierr); 271 if (cg->no_scaling){ 272 cg->diag_scaling = PETSC_FALSE; 273 cg->alpha = -1.0; 274 } 275 if (cg->alpha == -1.0 && cg->cg_type == CG_KouDai && !cg->diag_scaling){ /* Some more default options that appear to be good. */ 276 cg->neg_xi = PETSC_TRUE; 277 } 278 ierr = PetscOptionsBool("-tao_bncg_neg_xi","(developer) Use negative xi when it might be a smaller descent direction than necessary","",cg->neg_xi,&cg->neg_xi,NULL);CHKERRQ(ierr); 279 ierr = PetscOptionsReal("-tao_bncg_as_tol", "(developer) initial tolerance used when estimating actively bounded variables","",cg->as_tol,&cg->as_tol,NULL);CHKERRQ(ierr); 280 ierr = PetscOptionsReal("-tao_bncg_as_step", "(developer) step length used when estimating actively bounded variables","",cg->as_step,&cg->as_step,NULL);CHKERRQ(ierr); 281 ierr = PetscOptionsReal("-tao_bncg_delta_min", "(developer) minimum scaling factor used for scaled gradient restarts","",cg->delta_min,&cg->delta_min,NULL);CHKERRQ(ierr); 282 ierr = PetscOptionsReal("-tao_bncg_delta_max", "(developer) maximum scaling factor used for scaled gradient restarts","",cg->delta_max,&cg->delta_max,NULL);CHKERRQ(ierr); 283 284 ierr = PetscOptionsTail();CHKERRQ(ierr); 285 ierr = MatSetFromOptions(cg->B);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 289 static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer) 290 { 291 PetscBool isascii; 292 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 293 PetscErrorCode ierr; 294 295 PetscFunctionBegin; 296 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 297 if (isascii) { 298 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 299 ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]);CHKERRQ(ierr); 300 ierr = PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %i\n", cg->skipped_updates);CHKERRQ(ierr); 301 ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %i\n", cg->resets);CHKERRQ(ierr); 302 ierr = PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %i\n", cg->pure_gd_steps);CHKERRQ(ierr); 303 ierr = PetscViewerASCIIPrintf(viewer, "Not a descent direction: %i\n", cg->descent_error);CHKERRQ(ierr); 304 ierr = PetscViewerASCIIPrintf(viewer, "Line search fails: %i\n", cg->ls_fails);CHKERRQ(ierr); 305 if (cg->diag_scaling){ 306 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 307 if (isascii) { 308 ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 309 ierr = MatView(cg->B, viewer);CHKERRQ(ierr); 310 ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 311 } 312 } 313 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 314 } 315 PetscFunctionReturn(0); 316 } 317 318 PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha) 319 { 320 PetscReal a, b, c, sig1, sig2; 321 322 PetscFunctionBegin; 323 *scale = 0.0; 324 325 if (1.0 == alpha){ 326 *scale = yts/yty; 327 } else if (0.0 == alpha){ 328 *scale = sts/yts; 329 } 330 else if (-1.0 == alpha) *scale = 1.0; 331 else { 332 a = yty; 333 b = yts; 334 c = sts; 335 a *= alpha; 336 b *= -(2.0*alpha - 1.0); 337 c *= alpha - 1.0; 338 sig1 = (-b + PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 339 sig2 = (-b - PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a); 340 /* accept the positive root as the scalar */ 341 if (sig1 > 0.0) { 342 *scale = sig1; 343 } else if (sig2 > 0.0) { 344 *scale = sig2; 345 } else { 346 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar"); 347 } 348 } 349 PetscFunctionReturn(0); 350 } 351 352 /*MC 353 TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method. 354 355 Options Database Keys: 356 + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled) 357 . -tao_bncg_eta <r> - restart tolerance 358 . -tao_bncg_type <taocg_type> - cg formula 359 . -tao_bncg_as_type <none,bertsekas> - active set estimation method 360 . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation 361 . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation 362 . -tao_bncg_eps <r> - cutoff used for determining whether or not we restart based on steplength each iteration, as well as determining whether or not we continue using the last stepdirection. Defaults to machine precision. 363 . -tao_bncg_theta <r> - convex combination parameter for the Broyden method 364 . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 365 . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods 366 . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method 367 . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method 368 . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method 369 . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method 370 . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True. 371 . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods 372 . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts 373 . -tao_bncg_zeta <r> - Scaling parameter in the KD method 374 . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps 375 . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps 376 . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart 377 . -tao_bncg_min_restart_num <i> - This number, x, makes sure there is a gradient descent step every x*n iterations, where n is the dimension of the problem 378 . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations 379 . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults. 380 . -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions 381 382 Notes: 383 CG formulas are: 384 "gd" - Gradient Descent 385 "fr" - Fletcher-Reeves 386 "pr" - Polak-Ribiere-Polyak 387 "prp" - Polak-Ribiere-Plus 388 "hs" - Hestenes-Steifel 389 "dy" - Dai-Yuan 390 "ssml_bfgs" - Self-Scaling Memoryless BFGS 391 "ssml_dfp" - Self-Scaling Memoryless DFP 392 "ssml_brdn" - Self-Scaling Memoryless Broyden 393 "hz" - Hager-Zhang (CG_DESCENT 5.3) 394 "dk" - Dai-Kou (2013) 395 "kd" - Kou-Dai (2015) 396 Level: beginner 397 M*/ 398 399 PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao) 400 { 401 TAO_BNCG *cg; 402 const char *morethuente_type = TAOLINESEARCHMT; 403 PetscErrorCode ierr; 404 405 PetscFunctionBegin; 406 tao->ops->setup = TaoSetUp_BNCG; 407 tao->ops->solve = TaoSolve_BNCG; 408 tao->ops->view = TaoView_BNCG; 409 tao->ops->setfromoptions = TaoSetFromOptions_BNCG; 410 tao->ops->destroy = TaoDestroy_BNCG; 411 412 /* Override default settings (unless already changed) */ 413 if (!tao->max_it_changed) tao->max_it = 2000; 414 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 415 416 /* Note: nondefault values should be used for nonlinear conjugate gradient */ 417 /* method. In particular, gtol should be less that 0.5; the value used in */ 418 /* Nocedal and Wright is 0.10. We use the default values for the */ 419 /* linesearch because it seems to work better. */ 420 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 421 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 422 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 423 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 424 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 425 426 ierr = PetscNewLog(tao,&cg);CHKERRQ(ierr); 427 tao->data = (void*)cg; 428 ierr = KSPInitializePackage();CHKERRQ(ierr); 429 ierr = MatCreate(PetscObjectComm((PetscObject)tao), &cg->B);CHKERRQ(ierr); 430 ierr = PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1);CHKERRQ(ierr); 431 ierr = MatSetOptionsPrefix(cg->B, "tao_bncg_");CHKERRQ(ierr); 432 ierr = MatSetType(cg->B, MATLMVMDIAGBRDN);CHKERRQ(ierr); 433 434 cg->pc = NULL; 435 436 cg->dk_eta = 0.5; 437 cg->hz_eta = 0.4; 438 cg->dynamic_restart = PETSC_FALSE; 439 cg->unscaled_restart = PETSC_FALSE; 440 cg->no_scaling = PETSC_FALSE; 441 cg->delta_min = 1e-7; 442 cg->delta_max = 100; 443 cg->theta = 1.0; 444 cg->hz_theta = 1.0; 445 cg->dfp_scale = 1.0; 446 cg->bfgs_scale = 1.0; 447 cg->zeta = 0.1; 448 cg->min_quad = 6; 449 cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/ 450 cg->xi = 1.0; 451 cg->neg_xi = PETSC_TRUE; 452 cg->spaced_restart = PETSC_FALSE; 453 cg->tol_quad = 1e-8; 454 cg->as_step = 0.001; 455 cg->as_tol = 0.001; 456 cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); /* Just a little tighter*/ 457 cg->as_type = CG_AS_BERTSEKAS; 458 cg->cg_type = CG_SSML_BFGS; 459 cg->recycle = PETSC_FALSE; 460 cg->alpha = 1.0; 461 cg->diag_scaling = PETSC_TRUE; 462 PetscFunctionReturn(0); 463 } 464 465 466 PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq) 467 { 468 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 469 PetscErrorCode ierr; 470 PetscReal scaling; 471 472 PetscFunctionBegin; 473 ++cg->resets; 474 scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23); 475 scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling)); 476 if (cg->unscaled_restart) { 477 scaling = 1.0; 478 ++cg->pure_gd_steps; 479 } 480 ierr = VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient);CHKERRQ(ierr); 481 /* Also want to reset our diagonal scaling with each restart */ 482 if (cg->diag_scaling) { 483 ierr = MatLMVMReset(cg->B, PETSC_FALSE);CHKERRQ(ierr); 484 } 485 PetscFunctionReturn(0); 486 } 487 488 PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold) 489 { 490 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 491 PetscReal quadinterp; 492 493 PetscFunctionBegin; 494 if (cg->f < cg->min_quad/10) { 495 *dynrestart = PETSC_FALSE; 496 PetscFunctionReturn(0); 497 } /* just skip this since this strategy doesn't work well for functions near zero */ 498 quadinterp = 2.0*(cg->f - fold)/(stepsize*(gd + gd_old)); 499 if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) ++cg->iter_quad; 500 else { 501 cg->iter_quad = 0; 502 *dynrestart = PETSC_FALSE; 503 } 504 if (cg->iter_quad >= cg->min_quad) { 505 cg->iter_quad = 0; 506 *dynrestart = PETSC_TRUE; 507 } 508 PetscFunctionReturn(0); 509 } 510 511 PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscBool pcgd_fallback) 512 { 513 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 514 PetscErrorCode ierr; 515 PetscReal gamma = 1.0, tau_k, beta; 516 PetscReal tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk=1.0, gd; 517 PetscReal gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg; 518 PetscInt dim; 519 PetscBool cg_restart = PETSC_FALSE; 520 PetscFunctionBegin; 521 522 /* Local curvature check to see if we need to restart */ 523 if (tao->niter >= 1 || cg->recycle){ 524 ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 525 ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 526 ynorm2 = ynorm*ynorm; 527 ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 528 if (step*dnorm < PETSC_MACHINE_EPSILON || step*dk_yk < PETSC_MACHINE_EPSILON){ 529 cg_restart = PETSC_TRUE; 530 ++cg->skipped_updates; 531 } 532 if (cg->spaced_restart) { 533 ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr); 534 if (tao->niter % (dim*cg->min_restart_num)) cg_restart = PETSC_TRUE; 535 } 536 } 537 /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */ 538 if (cg->spaced_restart){ 539 ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr); 540 if (0 == tao->niter % (6*dim)) cg_restart = PETSC_TRUE; 541 } 542 /* Compute the diagonal scaling vector if applicable */ 543 if (cg->diag_scaling) { 544 ierr = MatLMVMUpdate(cg->B, tao->solution, tao->gradient);CHKERRQ(ierr); 545 } 546 547 /* A note on diagonal scaling (to be added to paper): 548 For the FR, PR, PRP, and DY methods, the diagonally scaled versions 549 must be derived as a preconditioned CG method rather than as 550 a Hessian initialization like in the Broyden methods. */ 551 552 /* In that case, one writes the objective function as 553 f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k). 554 Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to 555 HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the 556 same under preconditioning. Note that A is diagonal, such that A^T = A. */ 557 558 /* This yields questions like what the dot product d_k^T y_k 559 should look like. HZ mistakenly treats that as the same under 560 preconditioning, but that is not necessarily true. */ 561 562 /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation, 563 we get d_k^T y_k = (d_k^T A_k^{-T} A_k g_k - d_k^T A_k^{-T} A_{k-1} g_{k-1}), 564 yielding d_k^T y_k = d_k^T g_k - d_k^T (A_k^{-T} A_{k-1} g_{k-1}), which is 565 NOT the same if our preconditioning matrix is updated between iterations. 566 This same issue is found when considering dot products of the form g_{k+1}^T y_k. */ 567 568 /* Compute CG step direction */ 569 if (cg_restart) { 570 ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 571 } else if (pcgd_fallback) { 572 /* Just like preconditioned CG */ 573 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 574 ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr); 575 } else if (ynorm2 > PETSC_MACHINE_EPSILON) { 576 switch (cg->cg_type) { 577 case CG_PCGradientDescent: 578 if (!cg->diag_scaling){ 579 if (!cg->no_scaling){ 580 cg->sts = step*step*dnorm*dnorm; 581 ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 582 } else { 583 tau_k = 1.0; 584 ++cg->pure_gd_steps; 585 } 586 ierr = VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient);CHKERRQ(ierr); 587 } else { 588 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 589 ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr); 590 } 591 break; 592 593 case CG_HestenesStiefel: 594 /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */ 595 if (!cg->diag_scaling){ 596 cg->sts = step*step*dnorm*dnorm; 597 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 598 ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 599 beta = tau_k*gkp1_yk/dk_yk; 600 ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 601 } else { 602 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 603 ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 604 beta = gkp1_yk/dk_yk; 605 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 606 } 607 break; 608 609 case CG_FletcherReeves: 610 ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 611 ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 612 ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 613 ynorm2 = ynorm*ynorm; 614 ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 615 if (!cg->diag_scaling){ 616 ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, step*step*dnorm*dnorm, &tau_k, cg->alpha);CHKERRQ(ierr); 617 beta = tau_k*gnorm2/gnorm2_old; 618 ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 619 } else { 620 ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); /* Before it's updated */ 621 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 622 ierr = VecDot(tao->gradient, cg->g_work, &tmp);CHKERRQ(ierr); 623 beta = tmp/gnorm2_old; 624 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 625 } 626 break; 627 628 case CG_PolakRibierePolyak: 629 snorm = step*dnorm; 630 if (!cg->diag_scaling){ 631 ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 632 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 633 ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 634 beta = tau_k*gkp1_yk/gnorm2_old; 635 ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 636 } else { 637 ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); 638 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 639 ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 640 beta = gkp1_yk/gnorm2_old; 641 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 642 } 643 break; 644 645 case CG_PolakRibierePlus: 646 ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr); 647 ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr); 648 ynorm2 = ynorm*ynorm; 649 if (!cg->diag_scaling){ 650 ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 651 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 652 ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 653 beta = tau_k*gkp1_yk/gnorm2_old; 654 beta = PetscMax(beta, 0.0); 655 ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 656 } else { 657 ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); /* Old gtDg */ 658 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 659 ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 660 beta = gkp1_yk/gnorm2_old; 661 beta = PetscMax(beta, 0.0); 662 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 663 } 664 break; 665 666 case CG_DaiYuan: 667 /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property." 668 SIAM Journal on optimization 10, no. 1 (1999): 177-182. */ 669 if (!cg->diag_scaling){ 670 ierr = VecDot(tao->stepdirection, tao->gradient, &gd);CHKERRQ(ierr); 671 ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 672 ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->yts, &tau_k, cg->alpha);CHKERRQ(ierr); 673 beta = tau_k*gnorm2/(gd - gd_old); 674 ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 675 } else { 676 ierr = MatMult(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 677 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 678 ierr = VecDot(cg->g_work, tao->gradient, >Dg);CHKERRQ(ierr); 679 ierr = VecDot(tao->stepdirection, cg->G_old, &gd_old);CHKERRQ(ierr); 680 ierr = VecDot(cg->d_work, cg->g_work, &dk_yk);CHKERRQ(ierr); 681 dk_yk = dk_yk - gd_old; 682 beta = gtDg/dk_yk; 683 ierr = VecScale(cg->d_work, beta);CHKERRQ(ierr); 684 ierr = VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work);CHKERRQ(ierr); 685 } 686 break; 687 688 case CG_HagerZhang: 689 /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent." 690 ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */ 691 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 692 ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 693 ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 694 snorm = dnorm*step; 695 cg->yts = step*dk_yk; 696 if (cg->use_dynamic_restart){ 697 ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr); 698 } 699 if (cg->dynamic_restart){ 700 ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 701 } else { 702 if (!cg->diag_scaling){ 703 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 704 ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 705 /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */ 706 tmp = gd/dk_yk; 707 beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)); 708 /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */ 709 beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm)); 710 /* d <- -t*g + beta*t*d */ 711 ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr); 712 } else { 713 /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 714 cg->yty = ynorm2; 715 cg->sts = snorm*snorm; 716 /* Apply the diagonal scaling to all my vectors */ 717 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 718 ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 719 ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 720 /* Construct the constant ytDgkp1 */ 721 ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 722 /* Construct the constant for scaling Dkyk in the update */ 723 tmp = gd/dk_yk; 724 ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 725 tau_k = -tau_k*gd/(dk_yk*dk_yk); 726 /* beta is the constant which adds the dk contribution */ 727 beta = gkp1_yk/dk_yk + cg->hz_theta*tau_k; /* HZ; (1.15) from DK 2013 */ 728 /* From HZ2013, modified to account for diagonal scaling*/ 729 ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr); 730 ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr); 731 beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm)); 732 /* Do the update */ 733 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 734 } 735 } 736 break; 737 738 case CG_DaiKou: 739 /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search." 740 SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */ 741 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 742 ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 743 ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 744 snorm = step*dnorm; 745 cg->yts = dk_yk*step; 746 if (!cg->diag_scaling){ 747 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 748 ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 749 /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */ 750 tmp = gd/dk_yk; 751 beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) + gd/(dnorm*dnorm)) - step*gd/dk_yk; 752 beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm)); 753 /* d <- -t*g + beta*t*d */ 754 ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 755 } else { 756 /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 757 cg->yty = ynorm2; 758 cg->sts = snorm*snorm; 759 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 760 ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 761 ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr); 762 /* Construct the constant ytDgkp1 */ 763 ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr); 764 ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 765 tau_k = tau_k*gd/(dk_yk*dk_yk); 766 tmp = gd/dk_yk; 767 /* beta is the constant which adds the dk contribution */ 768 beta = gkp1_yk/dk_yk - step*tmp - tau_k; 769 /* Update this for the last term in beta */ 770 ierr = VecDot(cg->y_work, tao->stepdirection, &dk_yk);CHKERRQ(ierr); 771 beta += tmp*dk_yk/(dnorm*dnorm); /* projection of y_work onto dk */ 772 ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr); 773 ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr); 774 beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm)); 775 /* Do the update */ 776 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 777 } 778 break; 779 780 case CG_KouDai: 781 /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden–Fletcher–Goldfarb–Shanno method for unconstrained optimization." 782 Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */ 783 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 784 ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 785 ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 786 snorm = step*dnorm; 787 cg->yts = dk_yk*step; 788 if (cg->use_dynamic_restart){ 789 ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr); 790 } 791 if (cg->dynamic_restart){ 792 ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 793 } else { 794 if (!cg->diag_scaling){ 795 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 796 ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr); 797 beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 798 if (beta < cg->zeta*tau_k*gd/(dnorm*dnorm)) /* 0.1 is KD's zeta parameter */ 799 { 800 beta = cg->zeta*tau_k*gd/(dnorm*dnorm); 801 gamma = 0.0; 802 } else { 803 if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0*gd/dk_yk; 804 /* This seems to be very effective when there's no tau_k scaling. 805 This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */ 806 else { 807 gamma = cg->xi*gd/dk_yk; 808 } 809 } 810 /* d <- -t*g + beta*t*d + t*tmp*yk */ 811 ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 812 } else { 813 /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 814 cg->yty = ynorm2; 815 cg->sts = snorm*snorm; 816 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 817 ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 818 /* Construct the constant ytDgkp1 */ 819 ierr = VecDot(cg->yk, cg->g_work, &gkp1D_yk);CHKERRQ(ierr); 820 /* Construct the constant for scaling Dkyk in the update */ 821 gamma = gd/dk_yk; 822 /* tau_k = -ytDy/(ytd)^2 * gd */ 823 ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr); 824 tau_k = tau_k*gd/(dk_yk*dk_yk); 825 /* beta is the constant which adds the d_k contribution */ 826 beta = gkp1D_yk/dk_yk - step*gamma - tau_k; 827 /* Here is the requisite check */ 828 ierr = VecDot(tao->stepdirection, cg->g_work, &tmp);CHKERRQ(ierr); 829 if (cg->neg_xi){ 830 /* modified KD implementation */ 831 if (gkp1D_yk/dk_yk < 0) gamma = -1.0*gd/dk_yk; 832 else { 833 gamma = cg->xi*gd/dk_yk; 834 } 835 if (beta < cg->zeta*tmp/(dnorm*dnorm)){ 836 beta = cg->zeta*tmp/(dnorm*dnorm); 837 gamma = 0.0; 838 } 839 } else { /* original KD 2015 implementation */ 840 if (beta < cg->zeta*tmp/(dnorm*dnorm)) { 841 beta = cg->zeta*tmp/(dnorm*dnorm); 842 gamma = 0.0; 843 } else { 844 gamma = cg->xi*gd/dk_yk; 845 } 846 } 847 /* Do the update in two steps */ 848 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr); 849 ierr = VecAXPY(tao->stepdirection, gamma, cg->y_work);CHKERRQ(ierr); 850 } 851 } 852 break; 853 854 case CG_SSML_BFGS: 855 /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory." 856 Discussion Papers 269 (1977). */ 857 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 858 ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 859 snorm = step*dnorm; 860 cg->yts = dk_yk*step; 861 cg->yty = ynorm2; 862 cg->sts = snorm*snorm; 863 if (!cg->diag_scaling){ 864 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 865 ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 866 tmp = gd/dk_yk; 867 beta = tau_k*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*tmp; 868 /* d <- -t*g + beta*t*d + t*tmp*yk */ 869 ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 870 } else { 871 /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */ 872 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 873 ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 874 /* compute scalar gamma */ 875 ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 876 ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 877 gamma = gd/dk_yk; 878 /* Compute scalar beta */ 879 beta = (gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk; 880 /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 881 ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 882 } 883 break; 884 885 case CG_SSML_DFP: 886 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 887 ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 888 snorm = step*dnorm; 889 cg->yts = dk_yk*step; 890 cg->yty = ynorm2; 891 cg->sts = snorm*snorm; 892 if (!cg->diag_scaling){ 893 /* Instead of a regular convex combination, we will solve a quadratic formula. */ 894 ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr); 895 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 896 tau_k = cg->dfp_scale*tau_k; 897 tmp = tau_k*gkp1_yk/cg->yty; 898 beta = -step*gd/dk_yk; 899 /* d <- -t*g + beta*d + tmp*yk */ 900 ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 901 } else { 902 /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */ 903 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 904 ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 905 /* compute scalar gamma */ 906 ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 907 ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 908 gamma = (gkp1_yk/tmp); 909 /* Compute scalar beta */ 910 beta = -step*gd/dk_yk; 911 /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 912 ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 913 } 914 break; 915 916 case CG_SSML_BROYDEN: 917 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 918 ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr); 919 snorm = step*dnorm; 920 cg->yts = step*dk_yk; 921 cg->yty = ynorm2; 922 cg->sts = snorm*snorm; 923 if (!cg->diag_scaling){ 924 /* Instead of a regular convex combination, we will solve a quadratic formula. */ 925 ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_bfgs, cg->bfgs_scale);CHKERRQ(ierr); 926 ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_dfp, cg->dfp_scale);CHKERRQ(ierr); 927 ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr); 928 tau_k = cg->theta*tau_bfgs + (1.0-cg->theta)*tau_dfp; 929 /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0, 930 it should reproduce the tau_dfp scaling. Same with dfp_scale. */ 931 tmp = cg->theta*tau_bfgs*gd/dk_yk + (1-cg->theta)*tau_dfp*gkp1_yk/cg->yty; 932 beta = cg->theta*tau_bfgs*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*gd/dk_yk; 933 /* d <- -t*g + beta*d + tmp*yk */ 934 ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr); 935 } else { 936 /* We have diagonal scaling enabled */ 937 ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr); 938 ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr); 939 /* compute scalar gamma */ 940 ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr); 941 ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr); 942 gamma = cg->theta*gd/dk_yk + (1-cg->theta)*(gkp1_yk/tmp); 943 /* Compute scalar beta */ 944 beta = cg->theta*(gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk; 945 /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */ 946 ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr); 947 } 948 break; 949 950 default: 951 break; 952 953 } 954 } 955 PetscFunctionReturn(0); 956 } 957 958 PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm) 959 { 960 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 961 PetscErrorCode ierr; 962 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 963 PetscReal step=1.0,gnorm2,gd,dnorm=0.0; 964 PetscReal gnorm2_old,f_old,resnorm, gnorm_old; 965 PetscBool pcgd_fallback = PETSC_FALSE; 966 967 PetscFunctionBegin; 968 /* We are now going to perform a line search along the direction. */ 969 /* Store solution and gradient info before it changes */ 970 ierr = VecCopy(tao->solution, cg->X_old);CHKERRQ(ierr); 971 ierr = VecCopy(tao->gradient, cg->G_old);CHKERRQ(ierr); 972 ierr = VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old);CHKERRQ(ierr); 973 974 gnorm_old = gnorm; 975 gnorm2_old = gnorm_old*gnorm_old; 976 f_old = cg->f; 977 /* Perform bounded line search. If we are recycling a solution from a previous */ 978 /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */ 979 if (!(cg->recycle && 0 == tao->niter)){ 980 /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */ 981 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 982 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 983 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 984 985 /* Check linesearch failure */ 986 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 987 ++cg->ls_fails; 988 if (cg->cg_type == CG_GradientDescent){ 989 /* Nothing left to do but fail out of the optimization */ 990 step = 0.0; 991 tao->reason = TAO_DIVERGED_LS_FAILURE; 992 } else { 993 /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */ 994 ierr = VecCopy(cg->X_old, tao->solution);CHKERRQ(ierr); 995 ierr = VecCopy(cg->G_old, tao->gradient);CHKERRQ(ierr); 996 ierr = VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient);CHKERRQ(ierr); 997 gnorm = gnorm_old; 998 gnorm2 = gnorm2_old; 999 cg->f = f_old; 1000 1001 /* Fall back on preconditioned CG (so long as you're not already using it) */ 1002 if (cg->cg_type != CG_PCGradientDescent && cg->diag_scaling){ 1003 pcgd_fallback = PETSC_TRUE; 1004 ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback);CHKERRQ(ierr); 1005 1006 ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1007 ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1008 1009 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 1010 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1011 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1012 1013 pcgd_fallback = PETSC_FALSE; 1014 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 1015 /* Going to perform a regular gradient descent step. */ 1016 ++cg->ls_fails; 1017 step = 0.0; 1018 } 1019 } 1020 /* Fall back on the scaled gradient step */ 1021 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 1022 ++cg->ls_fails; 1023 ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1024 ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1025 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 1026 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 1027 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 1028 } 1029 1030 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 1031 /* Nothing left to do but fail out of the optimization */ 1032 ++cg->ls_fails; 1033 step = 0.0; 1034 tao->reason = TAO_DIVERGED_LS_FAILURE; 1035 } else { 1036 /* One of the fallbacks worked. Set them both back equal to false. */ 1037 pcgd_fallback = PETSC_FALSE; 1038 } 1039 } 1040 } 1041 /* Convergence test for line search failure */ 1042 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1043 1044 /* Standard convergence test */ 1045 ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 1046 ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 1047 if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 1048 ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 1049 ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 1050 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 1051 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 1052 } 1053 /* Assert we have an updated step and we need at least one more iteration. */ 1054 /* Calculate the next direction */ 1055 /* Estimate the active set at the new solution */ 1056 ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 1057 /* Compute the projected gradient and its norm */ 1058 ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 1059 ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 1060 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 1061 gnorm2 = gnorm*gnorm; 1062 1063 /* Calculate some quantities used in the StepDirectionUpdate. */ 1064 ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 1065 /* Update the step direction. */ 1066 ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback);CHKERRQ(ierr); 1067 ++tao->niter; 1068 ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1069 1070 if (cg->cg_type != CG_GradientDescent) { 1071 /* Figure out which previously active variables became inactive this iteration */ 1072 ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 1073 if (cg->inactive_idx && cg->inactive_old) { 1074 ierr = ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives);CHKERRQ(ierr); 1075 } 1076 /* Selectively reset the CG step those freshly inactive variables */ 1077 if (cg->new_inactives) { 1078 ierr = VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1079 ierr = VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1080 ierr = VecCopy(cg->inactive_grad, cg->inactive_step);CHKERRQ(ierr); 1081 ierr = VecScale(cg->inactive_step, -1.0);CHKERRQ(ierr); 1082 ierr = VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 1083 ierr = VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 1084 } 1085 /* Verify that this is a descent direction */ 1086 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 1087 ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 1088 if (PetscIsInfOrNanReal(gd) || (gd/(dnorm*dnorm) <= -1e10 || gd/(dnorm*dnorm) >= -1e-10)) { 1089 /* Not a descent direction, so we reset back to projected gradient descent */ 1090 ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr); 1091 ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr); 1092 ++cg->descent_error; 1093 } else { 1094 } 1095 } 1096 PetscFunctionReturn(0); 1097 } 1098 1099 PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0) 1100 { 1101 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 1102 PetscErrorCode ierr; 1103 1104 PetscFunctionBegin; 1105 ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr); 1106 cg->pc = H0; 1107 PetscFunctionReturn(0); 1108 } 1109