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