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