1 #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 2 3 typedef struct { 4 SNES snes; 5 /* Information on the regular SNES convergence test; which may have been user provided */ 6 PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 7 PetscErrorCode (*convdestroy)(void*); 8 void *convctx; 9 } SNES_TR_KSPConverged_Ctx; 10 11 static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx) 12 { 13 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 14 SNES snes = ctx->snes; 15 SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 16 Vec x; 17 PetscReal nrm; 18 PetscErrorCode ierr; 19 20 PetscFunctionBegin; 21 ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr); 22 if (*reason) { 23 ierr = PetscInfo(snes,"Default or user provided convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr); 24 } 25 /* Determine norm of solution */ 26 ierr = KSPBuildSolution(ksp,NULL,&x);CHKERRQ(ierr); 27 ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr); 28 if (nrm >= neP->delta) { 29 ierr = PetscInfo(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr); 30 *reason = KSP_CONVERGED_STEP_LENGTH; 31 } 32 PetscFunctionReturn(0); 33 } 34 35 static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 36 { 37 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx; 38 PetscErrorCode ierr; 39 40 PetscFunctionBegin; 41 ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr); 42 ierr = PetscFree(ctx);CHKERRQ(ierr); 43 PetscFunctionReturn(0); 44 } 45 46 /* ---------------------------------------------------------------- */ 47 /* 48 SNESTR_Converged_Private -test convergence JUST for 49 the trust region tolerance. 50 51 */ 52 static PetscErrorCode SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy) 53 { 54 SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 55 PetscErrorCode ierr; 56 57 PetscFunctionBegin; 58 *reason = SNES_CONVERGED_ITERATING; 59 if (neP->delta < xnorm * snes->deltatol) { 60 ierr = PetscInfo(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr); 61 *reason = SNES_DIVERGED_TR_DELTA; 62 } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 63 ierr = PetscInfo(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr); 64 *reason = SNES_DIVERGED_FUNCTION_COUNT; 65 } 66 PetscFunctionReturn(0); 67 } 68 69 /*@C 70 SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 71 Allows the user a chance to change or override the decision of the line search routine. 72 73 Logically Collective on snes 74 75 Input Parameters: 76 + snes - the nonlinear solver object 77 . func - [optional] function evaluation routine, see SNESNewtonTRPreCheck() for the calling sequence 78 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 79 80 Level: intermediate 81 82 Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver. 83 84 .seealso: SNESNewtonTRPreCheck(), SNESNewtonTRGetPreCheck(), SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck() 85 @*/ 86 PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,PetscBool*,void*),void *ctx) 87 { 88 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 89 90 PetscFunctionBegin; 91 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 92 if (func) tr->precheck = func; 93 if (ctx) tr->precheckctx = ctx; 94 PetscFunctionReturn(0); 95 } 96 97 /*@C 98 SNESNewtonTRGetPreCheck - Gets the pre-check function 99 100 Not collective 101 102 Input Parameter: 103 . snes - the nonlinear solver context 104 105 Output Parameters: 106 + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPreCheck() 107 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 108 109 Level: intermediate 110 111 .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRPreCheck() 112 @*/ 113 PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,PetscBool*,void*),void **ctx) 114 { 115 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 116 117 PetscFunctionBegin; 118 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 119 if (func) *func = tr->precheck; 120 if (ctx) *ctx = tr->precheckctx; 121 PetscFunctionReturn(0); 122 } 123 124 /*@C 125 SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 126 function evaluation. Allows the user a chance to change or override the decision of the line search routine 127 128 Logically Collective on snes 129 130 Input Parameters: 131 + snes - the nonlinear solver object 132 . func - [optional] function evaluation routine, see SNESNewtonTRPostCheck() for the calling sequence 133 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 134 135 Level: intermediate 136 137 Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver while the function set in 138 SNESLineSearchSetPostCheck() is called AFTER the function evaluation. 139 140 .seealso: SNESNewtonTRPostCheck(), SNESNewtonTRGetPostCheck() 141 @*/ 142 PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void *ctx) 143 { 144 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 145 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 148 if (func) tr->postcheck = func; 149 if (ctx) tr->postcheckctx = ctx; 150 PetscFunctionReturn(0); 151 } 152 153 /*@C 154 SNESNewtonTRGetPostCheck - Gets the post-check function 155 156 Not collective 157 158 Input Parameter: 159 . snes - the nonlinear solver context 160 161 Output Parameters: 162 + func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck() 163 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 164 165 Level: intermediate 166 167 .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRPostCheck() 168 @*/ 169 PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes,PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,PetscBool*,void*),void **ctx) 170 { 171 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 172 173 PetscFunctionBegin; 174 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 175 if (func) *func = tr->postcheck; 176 if (ctx) *ctx = tr->postcheckctx; 177 PetscFunctionReturn(0); 178 } 179 180 /*@C 181 SNESNewtonTRPreCheck - Called before the step has been determined in SNESNEWTONTR 182 183 Logically Collective on snes 184 185 Input Parameters: 186 + snes - the solver 187 . X - The last solution 188 - Y - The step direction 189 190 Output Parameters: 191 . changed_Y - Indicator that the step direction Y has been changed. 192 193 Level: developer 194 195 .seealso: SNESNewtonTRSetPreCheck(), SNESNewtonTRGetPreCheck() 196 @*/ 197 static PetscErrorCode SNESNewtonTRPreCheck(SNES snes,Vec X,Vec Y,PetscBool *changed_Y) 198 { 199 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 *changed_Y = PETSC_FALSE; 204 if (tr->precheck) { 205 ierr = (*tr->precheck)(snes,X,Y,changed_Y,tr->precheckctx);CHKERRQ(ierr); 206 PetscValidLogicalCollectiveBool(snes,*changed_Y,4); 207 } 208 PetscFunctionReturn(0); 209 } 210 211 /*@C 212 SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation 213 214 Logically Collective on snes 215 216 Input Parameters: 217 + snes - the solver 218 . X - The last solution 219 . Y - The full step direction 220 - W - The updated solution, W = X - Y 221 222 Output Parameters: 223 + changed_Y - indicator if step has been changed 224 - changed_W - Indicator if the new candidate solution W has been changed. 225 226 Notes: 227 If Y is changed then W is recomputed as X - Y 228 229 Level: developer 230 231 .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck() 232 @*/ 233 static PetscErrorCode SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_Y,PetscBool *changed_W) 234 { 235 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 236 PetscErrorCode ierr; 237 238 PetscFunctionBegin; 239 *changed_Y = PETSC_FALSE; 240 *changed_W = PETSC_FALSE; 241 if (tr->postcheck) { 242 ierr = (*tr->postcheck)(snes,X,Y,W,changed_Y,changed_W,tr->postcheckctx);CHKERRQ(ierr); 243 PetscValidLogicalCollectiveBool(snes,*changed_Y,5); 244 PetscValidLogicalCollectiveBool(snes,*changed_W,6); 245 } 246 PetscFunctionReturn(0); 247 } 248 249 /* 250 SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust 251 region approach for solving systems of nonlinear equations. 252 253 */ 254 static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 255 { 256 SNES_NEWTONTR *neP = (SNES_NEWTONTR*)snes->data; 257 Vec X,F,Y,G,Ytmp,W; 258 PetscErrorCode ierr; 259 PetscInt maxits,i,lits; 260 PetscReal rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1; 261 PetscScalar cnorm; 262 KSP ksp; 263 SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 264 PetscBool breakout = PETSC_FALSE; 265 SNES_TR_KSPConverged_Ctx *ctx; 266 PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),(*convdestroy)(void*); 267 void *convctx; 268 269 PetscFunctionBegin; 270 PetscCheckFalse(snes->xl || snes->xu || snes->ops->computevariablebounds,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 271 272 maxits = snes->max_its; /* maximum number of iterations */ 273 X = snes->vec_sol; /* solution vector */ 274 F = snes->vec_func; /* residual vector */ 275 Y = snes->work[0]; /* work vectors */ 276 G = snes->work[1]; 277 Ytmp = snes->work[2]; 278 W = snes->work[3]; 279 280 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 281 snes->iter = 0; 282 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 283 284 /* Set the linear stopping criteria to use the More' trick. */ 285 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 286 ierr = KSPGetConvergenceTest(ksp,&convtest,&convctx,&convdestroy);CHKERRQ(ierr); 287 if (convtest != SNESTR_KSPConverged_Private) { 288 ierr = PetscNew(&ctx);CHKERRQ(ierr); 289 ctx->snes = snes; 290 ierr = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr); 291 ierr = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr); 292 ierr = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr); 293 } 294 295 if (!snes->vec_func_init_set) { 296 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); /* F(X) */ 297 } else snes->vec_func_init_set = PETSC_FALSE; 298 299 ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- || F || */ 300 SNESCheckFunctionNorm(snes,fnorm); 301 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); /* xnorm <- || X || */ 302 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 303 snes->norm = fnorm; 304 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 305 delta = xnorm ? neP->delta0*xnorm : neP->delta0; 306 neP->delta = delta; 307 ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 308 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 309 310 /* test convergence */ 311 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 312 if (snes->reason) PetscFunctionReturn(0); 313 314 for (i=0; i<maxits; i++) { 315 316 /* Call general purpose update function */ 317 if (snes->ops->update) { 318 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 319 } 320 321 /* Solve J Y = F, where J is Jacobian matrix */ 322 ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 323 SNESCheckJacobianDomainerror(snes); 324 ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 325 ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr); 326 ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr); 327 snes->linear_its += lits; 328 329 ierr = PetscInfo(snes,"iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n",snes->iter,lits);CHKERRQ(ierr); 330 ierr = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr); 331 norm1 = nrm; 332 333 while (1) { 334 PetscBool changed_y; 335 PetscBool changed_w; 336 ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr); 337 nrm = norm1; 338 339 /* Scale Y if need be and predict new value of F norm */ 340 if (nrm >= delta) { 341 nrm = delta/nrm; 342 gpnorm = (1.0 - nrm)*fnorm; 343 cnorm = nrm; 344 ierr = PetscInfo(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr); 345 ierr = VecScale(Y,cnorm);CHKERRQ(ierr); 346 nrm = gpnorm; 347 ynorm = delta; 348 } else { 349 gpnorm = 0.0; 350 ierr = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr); 351 ynorm = nrm; 352 } 353 /* PreCheck() allows for updates to Y prior to W <- X - Y */ 354 ierr = SNESNewtonTRPreCheck(snes,X,Y,&changed_y);CHKERRQ(ierr); 355 ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr); /* W <- X - Y */ 356 ierr = SNESNewtonTRPostCheck(snes,X,Y,W,&changed_y,&changed_w);CHKERRQ(ierr); 357 if (changed_y) {ierr = VecWAXPY(W,-1.0,Y,X);CHKERRQ(ierr);} 358 ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr); 359 ierr = SNESComputeFunction(snes,W,G);CHKERRQ(ierr); /* F(X-Y) = G */ 360 ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr); /* gnorm <- || g || */ 361 SNESCheckFunctionNorm(snes,gnorm); 362 if (fnorm == gpnorm) rho = 0.0; 363 else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm); 364 365 /* Update size of trust region */ 366 if (rho < neP->mu) delta *= neP->delta1; 367 else if (rho < neP->eta) delta *= neP->delta2; 368 else delta *= neP->delta3; 369 ierr = PetscInfo(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr); 370 ierr = PetscInfo(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr); 371 372 neP->delta = delta; 373 if (rho > neP->sigma) break; 374 ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr); 375 376 /* check to see if progress is hopeless */ 377 neP->itflag = PETSC_FALSE; 378 ierr = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 379 if (!reason) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);} 380 if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 381 if (reason) { 382 /* We're not progressing, so return with the current iterate */ 383 ierr = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr); 384 breakout = PETSC_TRUE; 385 break; 386 } 387 snes->numFailures++; 388 } 389 if (!breakout) { 390 /* Update function and solution vectors */ 391 fnorm = gnorm; 392 ierr = VecCopy(G,F);CHKERRQ(ierr); 393 ierr = VecCopy(W,X);CHKERRQ(ierr); 394 /* Monitor convergence */ 395 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 396 snes->iter = i+1; 397 snes->norm = fnorm; 398 snes->xnorm = xnorm; 399 snes->ynorm = ynorm; 400 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 401 ierr = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr); 402 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 403 /* Test for convergence, xnorm = || X || */ 404 neP->itflag = PETSC_TRUE; 405 if (snes->ops->converged != SNESConvergedSkip) {ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);} 406 ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr); 407 if (reason) break; 408 } else break; 409 } 410 411 if (i == maxits) { 412 ierr = PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",maxits);CHKERRQ(ierr); 413 if (!reason) reason = SNES_DIVERGED_MAX_IT; 414 } 415 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 416 snes->reason = reason; 417 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 418 if (convtest != SNESTR_KSPConverged_Private) { 419 ierr = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr); 420 ierr = PetscFree(ctx);CHKERRQ(ierr); 421 ierr = KSPSetConvergenceTest(ksp,convtest,convctx,convdestroy);CHKERRQ(ierr); 422 } 423 PetscFunctionReturn(0); 424 } 425 426 /*------------------------------------------------------------*/ 427 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 428 { 429 PetscErrorCode ierr; 430 431 PetscFunctionBegin; 432 ierr = SNESSetWorkVecs(snes,4);CHKERRQ(ierr); 433 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 434 PetscFunctionReturn(0); 435 } 436 437 PetscErrorCode SNESReset_NEWTONTR(SNES snes) 438 { 439 440 PetscFunctionBegin; 441 PetscFunctionReturn(0); 442 } 443 444 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 445 { 446 PetscErrorCode ierr; 447 448 PetscFunctionBegin; 449 ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr); 450 ierr = PetscFree(snes->data);CHKERRQ(ierr); 451 PetscFunctionReturn(0); 452 } 453 /*------------------------------------------------------------*/ 454 455 static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes) 456 { 457 SNES_NEWTONTR *ctx = (SNES_NEWTONTR*)snes->data; 458 PetscErrorCode ierr; 459 460 PetscFunctionBegin; 461 ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr); 462 ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr); 463 ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr); 464 ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr); 465 ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr); 466 ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr); 467 ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr); 468 ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr); 469 ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr); 470 ierr = PetscOptionsTail();CHKERRQ(ierr); 471 PetscFunctionReturn(0); 472 } 473 474 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer) 475 { 476 SNES_NEWTONTR *tr = (SNES_NEWTONTR*)snes->data; 477 PetscErrorCode ierr; 478 PetscBool iascii; 479 480 PetscFunctionBegin; 481 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 482 if (iascii) { 483 ierr = PetscViewerASCIIPrintf(viewer," Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr); 484 ierr = PetscViewerASCIIPrintf(viewer," mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr); 485 ierr = PetscViewerASCIIPrintf(viewer," delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);CHKERRQ(ierr); 486 } 487 PetscFunctionReturn(0); 488 } 489 /* ------------------------------------------------------------ */ 490 /*MC 491 SNESNEWTONTR - Newton based nonlinear solver that uses a trust region 492 493 Options Database: 494 + -snes_trtol <tol> - trust region tolerance 495 . -snes_tr_mu <mu> - trust region parameter 496 . -snes_tr_eta <eta> - trust region parameter 497 . -snes_tr_sigma <sigma> - trust region parameter 498 . -snes_tr_delta0 <delta0> - initial size of the trust region is delta0*norm2(x) 499 . -snes_tr_delta1 <delta1> - trust region parameter 500 . -snes_tr_delta2 <delta2> - trust region parameter 501 - -snes_tr_delta3 <delta3> - trust region parameter 502 503 The basic algorithm is taken from "The Minpack Project", by More', 504 Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development 505 of Mathematical Software", Wayne Cowell, editor. 506 507 Level: intermediate 508 509 .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance() 510 511 M*/ 512 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 513 { 514 SNES_NEWTONTR *neP; 515 PetscErrorCode ierr; 516 517 PetscFunctionBegin; 518 snes->ops->setup = SNESSetUp_NEWTONTR; 519 snes->ops->solve = SNESSolve_NEWTONTR; 520 snes->ops->destroy = SNESDestroy_NEWTONTR; 521 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 522 snes->ops->view = SNESView_NEWTONTR; 523 snes->ops->reset = SNESReset_NEWTONTR; 524 525 snes->usesksp = PETSC_TRUE; 526 snes->usesnpc = PETSC_FALSE; 527 528 snes->alwayscomputesfinalresidual = PETSC_TRUE; 529 530 ierr = PetscNewLog(snes,&neP);CHKERRQ(ierr); 531 snes->data = (void*)neP; 532 neP->mu = 0.25; 533 neP->eta = 0.75; 534 neP->delta = 0.0; 535 neP->delta0 = 0.2; 536 neP->delta1 = 0.3; 537 neP->delta2 = 0.75; 538 neP->delta3 = 2.0; 539 neP->sigma = 0.0001; 540 neP->itflag = PETSC_FALSE; 541 neP->rnorm0 = 0.0; 542 neP->ttol = 0.0; 543 PetscFunctionReturn(0); 544 } 545