1 #include <petsctaolinesearch.h> 2 #include <../src/tao/bound/impls/bncg/bncg.h> 3 4 #define CG_FletcherReeves 0 5 #define CG_PolakRibiere 1 6 #define CG_PolakRibierePlus 2 7 #define CG_HestenesStiefel 3 8 #define CG_DaiYuan 4 9 #define CG_Types 5 10 11 static const char *CG_Table[64] = {"fr", "pr", "prp", "hs", "dy"}; 12 13 #define CG_AS_NONE 0 14 #define CG_AS_BERTSEKAS 1 15 #define CG_AS_SIZE 2 16 17 static const char *CG_AS_TYPE[64] = {"none", "bertsekas"}; 18 19 PetscErrorCode TaoBNCGSetRecycleFlag(Tao tao, PetscBool recycle) 20 { 21 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 22 23 PetscFunctionBegin; 24 cg->recycle = recycle; 25 PetscFunctionReturn(0); 26 } 27 28 PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType) 29 { 30 PetscErrorCode ierr; 31 TAO_BNCG *cg = (TAO_BNCG *)tao->data; 32 33 PetscFunctionBegin; 34 ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr); 35 if (cg->inactive_idx) { 36 ierr = ISDuplicate(cg->inactive_idx, &cg->inactive_old);CHKERRQ(ierr); 37 ierr = ISCopy(cg->inactive_idx, cg->inactive_old);CHKERRQ(ierr); 38 } 39 switch (asType) { 40 case CG_AS_NONE: 41 ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr); 42 ierr = VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx);CHKERRQ(ierr); 43 ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr); 44 ierr = ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx);CHKERRQ(ierr); 45 break; 46 47 case CG_AS_BERTSEKAS: 48 /* Use gradient descent to estimate the active set */ 49 ierr = VecCopy(cg->unprojected_gradient, cg->W);CHKERRQ(ierr); 50 ierr = VecScale(cg->W, -1.0);CHKERRQ(ierr); 51 ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->as_step, &cg->as_tol, &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx);CHKERRQ(ierr); 52 53 default: 54 break; 55 } 56 PetscFunctionReturn(0); 57 } 58 59 PetscErrorCode TaoBNCGBoundStep(Tao tao, Vec step) 60 { 61 PetscErrorCode ierr; 62 TAO_BNCG *cg = (TAO_BNCG *)tao->data; 63 64 PetscFunctionBegin; 65 switch (cg->as_type) { 66 case CG_AS_NONE: 67 if (cg->active_idx) {ierr = VecISSet(step, cg->active_idx, 0.0);CHKERRQ(ierr);} 68 break; 69 70 case CG_AS_BERTSEKAS: 71 ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, step);CHKERRQ(ierr); 72 break; 73 74 default: 75 break; 76 } 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode TaoSolve_BNCG(Tao tao) 81 { 82 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 83 PetscErrorCode ierr; 84 TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING; 85 PetscReal step=1.0,gnorm,gnorm2,delta,gd,ginner,beta,dnorm,resnorm; 86 PetscReal gd_old,gnorm2_old,f_old; 87 PetscBool cg_restart; 88 89 PetscFunctionBegin; 90 /* Project the current point onto the feasible set */ 91 ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 92 ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr); 93 94 /* Project the initial point onto the feasible region */ 95 ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr); 96 97 if (!cg->recycle) { 98 /* Solver is not being recycled so just compute the objective function and criteria */ 99 ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient);CHKERRQ(ierr); 100 } else { 101 /* We are recycling, so we have to compute ||g_old||^2 for use in the CG step calculation */ 102 ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr); 103 } 104 ierr = VecNorm(cg->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr); 105 if (PetscIsInfOrNanReal(cg->f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN"); 106 107 /* Estimate the active set and compute the projected gradient */ 108 ierr = VecCopy(cg->unprojected_gradient, cg->W);CHKERRQ(ierr); 109 ierr = VecScale(cg->W, -1.0);CHKERRQ(ierr); 110 ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 111 112 /* Project the gradient and calculate the norm */ 113 ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 114 ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 115 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 116 gnorm2 = gnorm*gnorm; 117 118 /* Convergence check */ 119 tao->niter = 0; 120 tao->reason = TAO_CONTINUE_ITERATING; 121 ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 122 ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 123 ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 124 ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 125 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 126 if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0); 127 128 /* Start optimization iterations */ 129 cg->ls_fails = cg->broken_ortho = cg->descent_error = 0; 130 cg->resets = -1; 131 while (tao->reason == TAO_CONTINUE_ITERATING) { 132 /* Check restart conditions for using steepest descent */ 133 cg_restart = PETSC_FALSE; 134 ierr = VecDot(tao->gradient, cg->G_old, &ginner);CHKERRQ(ierr); 135 ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr); 136 if (tao->niter == 0 && !cg->recycle && dnorm != 0.0) { 137 /* 1) First iteration, with recycle disabled, and a non-zero previous step */ 138 cg_restart = PETSC_TRUE; 139 } else if (PetscAbsScalar(ginner) >= cg->eta * gnorm2) { 140 /* 2) Gradients are far from orthogonal */ 141 cg_restart = PETSC_TRUE; 142 cg->broken_ortho++; 143 } 144 145 /* Compute CG step */ 146 if (cg_restart) { 147 beta = 0.0; 148 cg->resets++; 149 } else { 150 switch (cg->cg_type) { 151 case CG_FletcherReeves: 152 beta = gnorm2 / gnorm2_old; 153 break; 154 155 case CG_PolakRibiere: 156 beta = (gnorm2 - ginner) / gnorm2_old; 157 break; 158 159 case CG_PolakRibierePlus: 160 beta = PetscMax((gnorm2-ginner)/gnorm2_old, 0.0); 161 break; 162 163 case CG_HestenesStiefel: 164 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 165 ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 166 beta = (gnorm2 - ginner) / (gd - gd_old); 167 break; 168 169 case CG_DaiYuan: 170 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 171 ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr); 172 beta = gnorm2 / (gd - gd_old); 173 break; 174 175 default: 176 beta = 0.0; 177 break; 178 } 179 } 180 181 /* Compute the direction d=-g + beta*d */ 182 ierr = VecAXPBY(tao->stepdirection, -1.0, beta, tao->gradient);CHKERRQ(ierr); 183 ierr = TaoBNCGBoundStep(tao, tao->stepdirection);CHKERRQ(ierr); 184 if (cg->inactive_old) { 185 /* Compute which new indexes that were active before became inactive this iteration */ 186 ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr); 187 ierr = ISDifference(cg->inactive_old, cg->inactive_idx, &cg->new_inactives); 188 /* Selectively reset the CG step those freshly inactive variables to be the gradient descent direction */ 189 ierr = VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 190 ierr = VecGetSubVector(tao->gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 191 ierr = VecCopy(cg->inactive_grad, cg->inactive_step);CHKERRQ(ierr); 192 ierr = VecScale(cg->inactive_step, -1.0);CHKERRQ(ierr); 193 ierr = VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr); 194 ierr = VecRestoreSubVector(tao->gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr); 195 } 196 197 /* Verify that this is a descent direction */ 198 ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr); 199 ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm); 200 if (gd > -cg->rho*PetscPowReal(dnorm, cg->pow)) { 201 /* Not a descent direction, so we reset back to projected gradient descent */ 202 ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, tao->gradient);CHKERRQ(ierr); 203 cg->resets++; 204 cg->descent_error++; 205 } 206 207 /* update initial steplength choice */ 208 delta = 1.0; 209 delta = PetscMax(delta, cg->delta_min); 210 delta = PetscMin(delta, cg->delta_max); 211 212 /* Store solution and gradient info before it changes */ 213 ierr = VecCopy(tao->solution, cg->X_old);CHKERRQ(ierr); 214 ierr = VecCopy(tao->gradient, cg->G_old);CHKERRQ(ierr); 215 ierr = VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old);CHKERRQ(ierr); 216 gnorm2_old = gnorm2; 217 f_old = cg->f; 218 219 /* Perform bounded line search */ 220 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta);CHKERRQ(ierr); 221 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 222 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 223 224 /* Check linesearch failure */ 225 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) { 226 cg->ls_fails++; 227 /* Restore previous point */ 228 gnorm2 = gnorm2_old; 229 cg->f = f_old; 230 ierr = VecCopy(cg->X_old, tao->solution);CHKERRQ(ierr); 231 ierr = VecCopy(cg->G_old, tao->gradient);CHKERRQ(ierr); 232 ierr = VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient);CHKERRQ(ierr); 233 234 /* Fall back on the unscaled gradient step */ 235 delta = 1.0; 236 ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr); 237 ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 238 ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, tao->stepdirection);CHKERRQ(ierr); 239 240 ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,delta);CHKERRQ(ierr); 241 ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr); 242 ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr); 243 244 if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER){ 245 cg->ls_fails++; 246 /* Restore previous point */ 247 gnorm2 = gnorm2_old; 248 cg->f = f_old; 249 ierr = VecCopy(cg->X_old, tao->solution);CHKERRQ(ierr); 250 ierr = VecCopy(cg->G_old, tao->gradient);CHKERRQ(ierr); 251 ierr = VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient);CHKERRQ(ierr); 252 253 /* Nothing left to do but fail out of the optimization */ 254 step = 0.0; 255 tao->reason = TAO_DIVERGED_LS_FAILURE; 256 } 257 } 258 259 /* Estimate the active set at the new solution */ 260 ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr); 261 262 /* Compute the projected gradient and its norm */ 263 ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr); 264 ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr); 265 ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr); 266 gnorm2 = gnorm*gnorm; 267 268 /* Convergence test */ 269 tao->niter++; 270 ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr); 271 ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr); 272 ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr); 273 ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr); 274 ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 275 } 276 PetscFunctionReturn(0); 277 } 278 279 static PetscErrorCode TaoSetUp_BNCG(Tao tao) 280 { 281 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 282 PetscErrorCode ierr; 283 284 PetscFunctionBegin; 285 if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);} 286 if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);} 287 if (!cg->W) {ierr = VecDuplicate(tao->solution,&cg->W);CHKERRQ(ierr);} 288 if (!cg->X_old) {ierr = VecDuplicate(tao->solution,&cg->X_old);CHKERRQ(ierr);} 289 if (!cg->G_old) {ierr = VecDuplicate(tao->gradient,&cg->G_old);CHKERRQ(ierr);} 290 if (!cg->unprojected_gradient) {ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient);CHKERRQ(ierr);} 291 if (!cg->unprojected_gradient_old) {ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient_old);CHKERRQ(ierr);} 292 PetscFunctionReturn(0); 293 } 294 295 static PetscErrorCode TaoDestroy_BNCG(Tao tao) 296 { 297 TAO_BNCG *cg = (TAO_BNCG*) tao->data; 298 PetscErrorCode ierr; 299 300 PetscFunctionBegin; 301 if (tao->setupcalled) { 302 ierr = VecDestroy(&cg->W);CHKERRQ(ierr); 303 ierr = VecDestroy(&cg->X_old);CHKERRQ(ierr); 304 ierr = VecDestroy(&cg->G_old);CHKERRQ(ierr); 305 ierr = VecDestroy(&cg->unprojected_gradient);CHKERRQ(ierr); 306 ierr = VecDestroy(&cg->unprojected_gradient_old);CHKERRQ(ierr); 307 } 308 ierr = PetscFree(tao->data);CHKERRQ(ierr); 309 PetscFunctionReturn(0); 310 } 311 312 static PetscErrorCode TaoSetFromOptions_BNCG(PetscOptionItems *PetscOptionsObject,Tao tao) 313 { 314 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 315 PetscErrorCode ierr; 316 317 PetscFunctionBegin; 318 ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 319 ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");CHKERRQ(ierr); 320 ierr = PetscOptionsReal("-tao_bncg_eta","restart tolerance", "", cg->eta,&cg->eta,NULL);CHKERRQ(ierr); 321 ierr = PetscOptionsReal("-tao_bncg_rho","descent direction tolerance", "", cg->rho,&cg->rho,NULL);CHKERRQ(ierr); 322 ierr = PetscOptionsReal("-tao_bncg_pow","descent direction exponent", "", cg->pow,&cg->pow,NULL);CHKERRQ(ierr); 323 ierr = PetscOptionsEList("-tao_bncg_type","cg formula", "", CG_Table, CG_Types, CG_Table[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr); 324 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); 325 ierr = PetscOptionsReal("-tao_bncg_delta_min","minimum delta value", "", cg->delta_min,&cg->delta_min,NULL);CHKERRQ(ierr); 326 ierr = PetscOptionsReal("-tao_bncg_delta_max","maximum delta value", "", cg->delta_max,&cg->delta_max,NULL);CHKERRQ(ierr); 327 ierr = PetscOptionsBool("-tao_bncg_recycle","enable recycling the existing solution and gradient at the start of a new solve","",cg->recycle,&cg->recycle,NULL);CHKERRQ(ierr); 328 ierr = PetscOptionsReal("-tao_bncg_as_tol", "initial tolerance used when estimating actively bounded variables", "", cg->as_tol, &cg->as_tol,NULL);CHKERRQ(ierr); 329 ierr = PetscOptionsReal("-tao_bncg_as_step", "step length used when estimating actively bounded variables", "", cg->as_step, &cg->as_step,NULL);CHKERRQ(ierr); 330 ierr = PetscOptionsTail();CHKERRQ(ierr); 331 PetscFunctionReturn(0); 332 } 333 334 static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer) 335 { 336 PetscBool isascii; 337 TAO_BNCG *cg = (TAO_BNCG*)tao->data; 338 PetscErrorCode ierr; 339 340 PetscFunctionBegin; 341 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 342 if (isascii) { 343 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 344 ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]);CHKERRQ(ierr); 345 ierr = PetscViewerASCIIPrintf(viewer, "Resets: %i\n", cg->resets);CHKERRQ(ierr); 346 ierr = PetscViewerASCIIPrintf(viewer, " Broken ortho: %i\n", cg->broken_ortho);CHKERRQ(ierr); 347 ierr = PetscViewerASCIIPrintf(viewer, " Not a descent dir.: %i\n", cg->descent_error);CHKERRQ(ierr); 348 ierr = PetscViewerASCIIPrintf(viewer, "Line search fails: %i\n", cg->ls_fails);CHKERRQ(ierr); 349 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 350 } 351 PetscFunctionReturn(0); 352 } 353 354 /*MC 355 TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method. 356 357 Options Database Keys: 358 + -tao_bncg_eta <r> - restart tolerance 359 . -tao_bncg_type <taocg_type> - cg formula 360 . -tao_bncg_delta_min <r> - minimum delta value 361 - -tao_bncg_delta_max <r> - maximum delta value 362 363 Notes: 364 CG formulas are: 365 "fr" - Fletcher-Reeves 366 "pr" - Polak-Ribiere 367 "prp" - Polak-Ribiere-Plus 368 "hs" - Hestenes-Steifel 369 "dy" - Dai-Yuan 370 Level: beginner 371 M*/ 372 373 374 PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao) 375 { 376 TAO_BNCG *cg; 377 const char *morethuente_type = TAOLINESEARCHMT; 378 PetscErrorCode ierr; 379 380 PetscFunctionBegin; 381 tao->ops->setup = TaoSetUp_BNCG; 382 tao->ops->solve = TaoSolve_BNCG; 383 tao->ops->view = TaoView_BNCG; 384 tao->ops->setfromoptions = TaoSetFromOptions_BNCG; 385 tao->ops->destroy = TaoDestroy_BNCG; 386 387 /* Override default settings (unless already changed) */ 388 if (!tao->max_it_changed) tao->max_it = 2000; 389 if (!tao->max_funcs_changed) tao->max_funcs = 4000; 390 391 /* Note: nondefault values should be used for nonlinear conjugate gradient */ 392 /* method. In particular, gtol should be less that 0.5; the value used in */ 393 /* Nocedal and Wright is 0.10. We use the default values for the */ 394 /* linesearch because it seems to work better. */ 395 ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 396 ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 397 ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr); 398 ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr); 399 ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 400 401 ierr = PetscNewLog(tao,&cg);CHKERRQ(ierr); 402 tao->data = (void*)cg; 403 cg->rho = 1e-4; 404 cg->pow = 2.1; 405 cg->eta = 0.5; 406 cg->delta_min = 1e-7; 407 cg->delta_max = 100; 408 cg->as_step = 0.001; 409 cg->as_tol = 0.001; 410 cg->as_type = CG_AS_BERTSEKAS; 411 cg->cg_type = CG_DaiYuan; 412 cg->recycle = PETSC_FALSE; 413 PetscFunctionReturn(0); 414 } 415