1 #include <../src/snes/impls/tr/trimpl.h> /*I "petscsnes.h" I*/ 2 3 typedef struct { 4 SNES snes; 5 PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 6 PetscErrorCode (*convdestroy)(void *); 7 void *convctx; 8 } SNES_TR_KSPConverged_Ctx; 9 10 const char *const SNESNewtonTRFallbackTypes[] = {"NEWTON", "CAUCHY", "DOGLEG", "SNESNewtonTRFallbackType", "SNES_TR_FALLBACK_", NULL}; 11 12 static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 13 { 14 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 15 SNES snes = ctx->snes; 16 SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 17 Vec x; 18 PetscReal nrm; 19 20 PetscFunctionBegin; 21 PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 22 if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 23 /* Determine norm of solution */ 24 PetscCall(KSPBuildSolution(ksp, NULL, &x)); 25 PetscCall(VecNorm(x, NORM_2, &nrm)); 26 if (nrm >= neP->delta) { 27 PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 28 *reason = KSP_CONVERGED_STEP_LENGTH; 29 } 30 PetscFunctionReturn(PETSC_SUCCESS); 31 } 32 33 static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 34 { 35 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 36 37 PetscFunctionBegin; 38 PetscCall((*ctx->convdestroy)(ctx->convctx)); 39 PetscCall(PetscFree(ctx)); 40 PetscFunctionReturn(PETSC_SUCCESS); 41 } 42 43 static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 44 { 45 SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 46 47 PetscFunctionBegin; 48 *reason = SNES_CONVERGED_ITERATING; 49 if (neP->delta < snes->deltatol) { 50 PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)snes->deltatol)); 51 *reason = SNES_DIVERGED_TR_DELTA; 52 } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 53 PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 54 *reason = SNES_DIVERGED_FUNCTION_COUNT; 55 } 56 PetscFunctionReturn(PETSC_SUCCESS); 57 } 58 59 /*@ 60 SNESNewtonTRSetFallbackType - Set the type of fallback if the solution of the trust region subproblem is outside the radius 61 62 Input Parameters: 63 + snes - the nonlinear solver object 64 - ftype - the fallback type, see `SNESNewtonTRFallbackType` 65 66 Level: intermediate 67 68 .seealso: `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`, 69 `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 70 @*/ 71 PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype) 72 { 73 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 74 PetscBool flg; 75 76 PetscFunctionBegin; 77 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 78 PetscValidLogicalCollectiveEnum(snes, ftype, 2); 79 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 80 if (flg) tr->fallback = ftype; 81 PetscFunctionReturn(PETSC_SUCCESS); 82 } 83 84 /*@C 85 SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 86 Allows the user a chance to change or override the trust region decision. 87 88 Logically Collective 89 90 Input Parameters: 91 + snes - the nonlinear solver object 92 . func - [optional] function evaluation routine, see `SNESNewtonTRPreCheck()` for the calling sequence 93 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 94 95 Level: intermediate 96 97 Note: 98 This function is called BEFORE the function evaluation within the solver. 99 100 .seealso: `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 101 @*/ 102 PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 103 { 104 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 105 PetscBool flg; 106 107 PetscFunctionBegin; 108 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 109 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 110 if (flg) { 111 if (func) tr->precheck = func; 112 if (ctx) tr->precheckctx = ctx; 113 } 114 PetscFunctionReturn(PETSC_SUCCESS); 115 } 116 117 /*@C 118 SNESNewtonTRGetPreCheck - Gets the pre-check function 119 120 Not collective 121 122 Input Parameter: 123 . snes - the nonlinear solver context 124 125 Output Parameters: 126 + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPreCheck()` 127 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 128 129 Level: intermediate 130 131 .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 132 @*/ 133 PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 134 { 135 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 136 PetscBool flg; 137 138 PetscFunctionBegin; 139 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 140 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 141 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 142 if (func) *func = tr->precheck; 143 if (ctx) *ctx = tr->precheckctx; 144 PetscFunctionReturn(PETSC_SUCCESS); 145 } 146 147 /*@C 148 SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 149 function evaluation. Allows the user a chance to change or override the internal decision of the solver 150 151 Logically Collective 152 153 Input Parameters: 154 + snes - the nonlinear solver object 155 . func - [optional] function evaluation routine, see `SNESNewtonTRPostCheck()` for the calling sequence 156 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 157 158 Level: intermediate 159 160 Note: 161 This function is called BEFORE the function evaluation within the solver while the function set in 162 `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 163 164 .seealso: `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 165 @*/ 166 PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 167 { 168 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 169 PetscBool flg; 170 171 PetscFunctionBegin; 172 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 173 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 174 if (flg) { 175 if (func) tr->postcheck = func; 176 if (ctx) tr->postcheckctx = ctx; 177 } 178 PetscFunctionReturn(PETSC_SUCCESS); 179 } 180 181 /*@C 182 SNESNewtonTRGetPostCheck - Gets the post-check function 183 184 Not collective 185 186 Input Parameter: 187 . snes - the nonlinear solver context 188 189 Output Parameters: 190 + func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPostCheck()` 191 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 192 193 Level: intermediate 194 195 .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 196 @*/ 197 PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 198 { 199 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 200 PetscBool flg; 201 202 PetscFunctionBegin; 203 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 204 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 205 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 206 if (func) *func = tr->postcheck; 207 if (ctx) *ctx = tr->postcheckctx; 208 PetscFunctionReturn(PETSC_SUCCESS); 209 } 210 211 /*@C 212 SNESNewtonTRPreCheck - Runs the precheck routine 213 214 Logically Collective 215 216 Input Parameters: 217 + snes - the solver 218 . X - The last solution 219 - Y - The step direction 220 221 Output Parameters: 222 . changed_Y - Indicator that the step direction Y has been changed. 223 224 Level: intermediate 225 226 .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()` 227 @*/ 228 PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 229 { 230 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 231 PetscBool flg; 232 233 PetscFunctionBegin; 234 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 235 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 236 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 237 *changed_Y = PETSC_FALSE; 238 if (tr->precheck) { 239 PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 240 PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 241 } 242 PetscFunctionReturn(PETSC_SUCCESS); 243 } 244 245 /*@C 246 SNESNewtonTRPostCheck - Runs the postcheck routine 247 248 Logically Collective 249 250 Input Parameters: 251 + snes - the solver 252 . X - The last solution 253 . Y - The full step direction 254 - W - The updated solution, W = X - Y 255 256 Output Parameters: 257 + changed_Y - indicator if step has been changed 258 - changed_W - Indicator if the new candidate solution W has been changed. 259 260 Note: 261 If Y is changed then W is recomputed as X - Y 262 263 Level: intermediate 264 265 .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck() 266 @*/ 267 PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 268 { 269 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 270 PetscBool flg; 271 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 274 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 275 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 276 *changed_Y = PETSC_FALSE; 277 *changed_W = PETSC_FALSE; 278 if (tr->postcheck) { 279 PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 280 PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 281 PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 282 } 283 PetscFunctionReturn(PETSC_SUCCESS); 284 } 285 286 static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 287 { 288 PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 289 PetscReal x1 = temp / a; 290 PetscReal x2 = c / temp; 291 *xm = PetscMin(x1, x2); 292 *xp = PetscMax(x1, x2); 293 } 294 295 /* 296 SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 297 (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 298 nonlinear equations 299 300 */ 301 static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 302 { 303 SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 304 Vec X, F, Y, G, W, GradF, YU; 305 PetscInt maxits, lits; 306 PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm; 307 PetscReal deltaM, fk, fkp1, deltaqm, gTy, yTHy; 308 PetscReal auk, gfnorm, ycnorm, gTBg; 309 KSP ksp; 310 PetscBool already_done = PETSC_FALSE; 311 PetscBool clear_converged_test, rho_satisfied; 312 PetscVoidFunction ksp_has_radius; 313 SNES_TR_KSPConverged_Ctx *ctx; 314 void *convctx; 315 PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 316 PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 317 318 PetscFunctionBegin; 319 PetscCall(SNESGetObjective(snes, &objective, NULL)); 320 321 maxits = snes->max_its; /* maximum number of iterations */ 322 X = snes->vec_sol; /* solution vector */ 323 F = snes->vec_func; /* residual vector */ 324 Y = snes->vec_sol_update; /* update vector */ 325 G = snes->work[0]; /* updated residual */ 326 W = snes->work[1]; /* temporary vector */ 327 GradF = !objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */ 328 YU = snes->work[3]; /* work vector for dogleg method */ 329 330 PetscCheck(!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); 331 332 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 333 snes->iter = 0; 334 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 335 336 /* Set the linear stopping criteria to use the More' trick if needed */ 337 clear_converged_test = PETSC_FALSE; 338 PetscCall(SNESGetKSP(snes, &ksp)); 339 PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 340 PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &ksp_has_radius)); 341 if (convtest != SNESTR_KSPConverged_Private && !ksp_has_radius) { 342 clear_converged_test = PETSC_TRUE; 343 PetscCall(PetscNew(&ctx)); 344 ctx->snes = snes; 345 PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 346 PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 347 PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 348 } 349 350 if (!snes->vec_func_init_set) { 351 PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 352 } else snes->vec_func_init_set = PETSC_FALSE; 353 354 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 355 SNESCheckFunctionNorm(snes, fnorm); 356 PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 357 358 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 359 snes->norm = fnorm; 360 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 361 delta = neP->delta0; 362 deltaM = neP->deltaM; 363 neP->delta = delta; 364 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 365 PetscCall(SNESMonitor(snes, 0, fnorm)); 366 367 /* test convergence */ 368 rho_satisfied = PETSC_FALSE; 369 PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 370 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 371 372 if (objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 373 else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 374 375 while (snes->iter < maxits) { 376 PetscBool changed_y; 377 PetscBool changed_w; 378 379 /* solve trust-region subproblem */ 380 if (!already_done) { 381 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 382 SNESCheckJacobianDomainerror(snes); 383 } 384 PetscCall(KSPCGSetRadius(snes->ksp, delta)); 385 PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre)); 386 PetscCall(KSPSolve(snes->ksp, F, Y)); 387 SNESCheckKSPSolve(snes); 388 PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 389 390 /* calculating GradF of minimization function only once */ 391 if (!already_done) { 392 if (objective) gfnorm = fnorm; 393 else { 394 PetscCall(MatMultTranspose(snes->jacobian, F, GradF)); /* grad f = J^T F */ 395 PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); 396 } 397 already_done = PETSC_TRUE; 398 } 399 400 /* decide what to do when the update is outside of trust region */ 401 PetscCall(VecNorm(Y, NORM_2, &ynorm)); 402 if (ynorm > delta) { 403 switch (neP->fallback) { 404 case SNES_TR_FALLBACK_NEWTON: 405 auk = delta / ynorm; 406 PetscCall(VecScale(Y, auk)); 407 break; 408 case SNES_TR_FALLBACK_CAUCHY: 409 case SNES_TR_FALLBACK_DOGLEG: 410 PetscCall(MatMult(snes->jacobian, GradF, W)); 411 if (objective) PetscCall(VecDotRealPart(GradF, W, &gTBg)); 412 else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */ 413 /* Eqs 4.7 and 4.8 in Nocedal and Wright */ 414 auk = delta / gfnorm; 415 if (gTBg > 0.0) auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1); 416 ycnorm = auk * gfnorm; 417 if (neP->fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) { 418 /* Cauchy solution */ 419 PetscCall(VecAXPBY(Y, auk, 0.0, GradF)); 420 PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg)); 421 } else { /* take linear combination of Cauchy and Newton direction and step */ 422 PetscReal c0, c1, c2, tau = 0.0, tpos, tneg; 423 PetscBool noroots; 424 425 auk = gfnorm * gfnorm / gTBg; 426 PetscCall(VecAXPBY(YU, auk, 0.0, GradF)); 427 PetscCall(VecAXPY(Y, -1.0, YU)); 428 PetscCall(VecNorm(Y, NORM_2, &c0)); 429 PetscCall(VecDotRealPart(YU, Y, &c1)); 430 c0 = PetscSqr(c0); 431 c2 = PetscSqr(ycnorm) - PetscSqr(delta); 432 PetscQuadraticRoots(c0, c1, c2, &tneg, &tpos); 433 434 noroots = PetscIsInfOrNanReal(tneg); 435 if (noroots) { /* No roots, select Cauchy point */ 436 auk = delta / gfnorm; 437 auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1); 438 PetscCall(VecAXPBY(Y, auk, 0.0, GradF)); 439 } else { /* Here roots corresponds to tau-1 in Nocedal and Wright */ 440 tpos += 1.0; 441 tneg += 1.0; 442 tau = PetscClipInterval(tpos, 0.0, 2.0); /* clip to tau [0,2] */ 443 if (tau < 1.0) { 444 PetscCall(VecAXPBY(Y, tau, 0.0, YU)); 445 } else { 446 PetscCall(VecAXPBY(Y, 1.0, tau - 1, YU)); 447 } 448 } 449 PetscCall(VecNorm(Y, NORM_2, &c0)); /* this norm will be cached and reused later */ 450 PetscCall(PetscInfo(snes, "%s evaluated. roots: (%g, %g), tau %g, ynorm: %g, ycnorm: %g, ydlnorm %g, gTBg: %g\n", noroots ? "CP" : "DL", (double)tneg, (double)tpos, (double)tau, (double)ynorm, (double)ycnorm, (double)c0, (double)gTBg)); 451 } 452 break; 453 default: 454 SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode"); 455 break; 456 } 457 } 458 459 /* Evaluate the solution to meet the improvement ratio criteria */ 460 461 /* compute the final ynorm */ 462 PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 463 PetscCall(VecNorm(Y, NORM_2, &ynorm)); 464 465 /* the quadratic model difference */ 466 PetscCall(MatMult(snes->jacobian, Y, W)); 467 if (objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 468 else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */ 469 PetscCall(VecDotRealPart(GradF, Y, &gTy)); 470 deltaqm = -(-gTy + 0.5 * yTHy); /* difference in quadratic model, -gTy because SNES solves it this way */ 471 472 /* update */ 473 PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */ 474 PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 475 if (changed_y) { 476 /* Need to recompute the quadratic model difference */ 477 PetscCall(MatMult(snes->jacobian, Y, W)); 478 if (objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 479 else PetscCall(VecDotRealPart(W, W, &yTHy)); 480 PetscCall(VecDotRealPart(GradF, Y, &gTy)); 481 deltaqm = -(-gTy + 0.5 * yTHy); 482 /* User changed Y but not W */ 483 if (!changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X)); 484 } 485 486 /* Compute new objective function */ 487 PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */ 488 PetscCall(VecNorm(G, NORM_2, &gnorm)); 489 if (objective) PetscCall(SNESComputeObjective(snes, W, &fkp1)); 490 else fkp1 = 0.5 * PetscSqr(gnorm); 491 SNESCheckFunctionNorm(snes, fkp1); 492 493 if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */ 494 else rho = -1.0; /* no reduction in quadratic model, step must be rejected */ 495 PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy)); 496 497 if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */ 498 else if (rho > neP->eta3) delta *= neP->t2; /* expand the region */ 499 delta = PetscMin(delta, deltaM); /* but not greater than deltaM */ 500 501 neP->delta = delta; 502 if (rho >= neP->eta1) { 503 rho_satisfied = PETSC_TRUE; 504 } else { 505 rho_satisfied = PETSC_FALSE; 506 PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 507 /* check to see if progress is hopeless */ 508 PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP)); 509 if (!snes->reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 510 if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_INNER; 511 snes->numFailures++; 512 /* We're not progressing, so return with the current iterate */ 513 if (snes->reason) break; 514 } 515 if (rho_satisfied) { 516 /* Update function values */ 517 already_done = PETSC_FALSE; 518 fnorm = gnorm; 519 fk = fkp1; 520 521 /* New residual and linearization point */ 522 PetscCall(VecCopy(G, F)); 523 PetscCall(VecCopy(W, X)); 524 525 /* Monitor convergence */ 526 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 527 snes->iter++; 528 snes->norm = fnorm; 529 snes->xnorm = xnorm; 530 snes->ynorm = ynorm; 531 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 532 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 533 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 534 535 /* Test for convergence, xnorm = || X || */ 536 PetscCall(VecNorm(X, NORM_2, &xnorm)); 537 PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP); 538 if (snes->reason) break; 539 } 540 } 541 542 if (snes->iter == maxits) { 543 PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 544 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 545 } 546 if (clear_converged_test) { 547 PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 548 PetscCall(PetscFree(ctx)); 549 PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 550 } 551 PetscFunctionReturn(PETSC_SUCCESS); 552 } 553 554 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 555 { 556 PetscFunctionBegin; 557 PetscCall(SNESSetWorkVecs(snes, 4)); 558 PetscCall(SNESSetUpMatrices(snes)); 559 PetscFunctionReturn(PETSC_SUCCESS); 560 } 561 562 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 563 { 564 PetscFunctionBegin; 565 PetscCall(PetscFree(snes->data)); 566 PetscFunctionReturn(PETSC_SUCCESS); 567 } 568 569 static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) 570 { 571 SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 572 573 PetscFunctionBegin; 574 PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 575 PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL)); 576 PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 577 PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 578 PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 579 PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 580 PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 581 PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 582 PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 583 PetscCall(PetscOptionsEnum("-snes_tr_fallback_type", "Type of fallback if subproblem solution is outside of the trust region", "SNESNewtonTRSetFallbackType", SNESNewtonTRFallbackTypes, (PetscEnum)ctx->fallback, (PetscEnum *)&ctx->fallback, NULL)); 584 PetscOptionsHeadEnd(); 585 PetscFunctionReturn(PETSC_SUCCESS); 586 } 587 588 static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 589 { 590 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 591 PetscBool iascii; 592 593 PetscFunctionBegin; 594 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 595 if (iascii) { 596 PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)snes->deltatol)); 597 PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 598 PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM)); 599 PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback])); 600 } 601 PetscFunctionReturn(PETSC_SUCCESS); 602 } 603 604 /*MC 605 SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 606 607 Options Database Keys: 608 + -snes_tr_tol <tol> - trust region tolerance 609 . -snes_tr_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 610 . -snes_tr_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 611 . -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 612 . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 613 . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 614 . -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL) 615 . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2) 616 - -snes_tr_fallback_type <newton,cauchy,dogleg> - Solution strategy to test reduction when step is outside of trust region. Can use scaled Newton direction, Cauchy point (Steepest Descent direction) or dogleg method. 617 618 Reference: 619 . * - "Numerical Optimization" by Nocedal and Wright, chapter 4. 620 621 Level: intermediate 622 623 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 624 `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 625 `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()` 626 M*/ 627 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 628 { 629 SNES_NEWTONTR *neP; 630 631 PetscFunctionBegin; 632 snes->ops->setup = SNESSetUp_NEWTONTR; 633 snes->ops->solve = SNESSolve_NEWTONTR; 634 snes->ops->destroy = SNESDestroy_NEWTONTR; 635 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 636 snes->ops->view = SNESView_NEWTONTR; 637 638 snes->usesksp = PETSC_TRUE; 639 snes->usesnpc = PETSC_FALSE; 640 641 snes->alwayscomputesfinalresidual = PETSC_TRUE; 642 643 PetscCall(PetscNew(&neP)); 644 snes->data = (void *)neP; 645 neP->delta = 0.0; 646 neP->delta0 = 0.2; 647 neP->eta1 = 0.001; 648 neP->eta2 = 0.25; 649 neP->eta3 = 0.75; 650 neP->t1 = 0.25; 651 neP->t2 = 2.0; 652 neP->deltaM = 1.e10; 653 neP->fallback = SNES_TR_FALLBACK_NEWTON; 654 PetscFunctionReturn(PETSC_SUCCESS); 655 } 656