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 const char *const SNESNewtonTRQNTypes[] = {"NONE", "SAME", "DIFFERENT", "SNESNewtonTRQNType", "SNES_TR_QN_", NULL}; 12 13 static PetscErrorCode SNESComputeJacobian_MATLMVM(SNES snes, Vec X, Mat J, Mat B, void *dummy) 14 { 15 PetscFunctionBegin; 16 // PetscCall(MatLMVMSymBroydenSetDelta(B, _some_delta)); 17 PetscCall(MatLMVMUpdate(B, X, snes->vec_func)); 18 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 19 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 20 if (J != B) { 21 // PetscCall(MatLMVMSymBroydenSetDelta(J, _some_delta)); 22 PetscCall(MatLMVMUpdate(J, X, snes->vec_func)); 23 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 24 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 25 } 26 PetscFunctionReturn(PETSC_SUCCESS); 27 } 28 29 static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 30 { 31 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 32 SNES snes = ctx->snes; 33 SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 34 Vec x; 35 PetscReal nrm; 36 37 PetscFunctionBegin; 38 PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx)); 39 if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); 40 /* Determine norm of solution */ 41 PetscCall(KSPBuildSolution(ksp, NULL, &x)); 42 PetscCall(VecNorm(x, neP->norm, &nrm)); 43 if (nrm >= neP->delta) { 44 PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm)); 45 *reason = KSP_CONVERGED_STEP_LENGTH; 46 } 47 PetscFunctionReturn(PETSC_SUCCESS); 48 } 49 50 static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) 51 { 52 SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx; 53 54 PetscFunctionBegin; 55 PetscCall((*ctx->convdestroy)(ctx->convctx)); 56 PetscCall(PetscFree(ctx)); 57 PetscFunctionReturn(PETSC_SUCCESS); 58 } 59 60 static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 61 { 62 SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 63 64 PetscFunctionBegin; 65 *reason = SNES_CONVERGED_ITERATING; 66 if (neP->delta < snes->deltatol) { 67 PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)snes->deltatol)); 68 *reason = SNES_DIVERGED_TR_DELTA; 69 } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 70 PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 71 *reason = SNES_DIVERGED_FUNCTION_COUNT; 72 } 73 PetscFunctionReturn(PETSC_SUCCESS); 74 } 75 76 /*@ 77 SNESNewtonTRSetNormType - Specify the type of norm to use for the computation of the trust region. 78 79 Input Parameters: 80 + snes - the nonlinear solver object 81 - norm - the norm type 82 83 Level: intermediate 84 85 .seealso: `SNESNEWTONTR`, `NormType` 86 @*/ 87 PetscErrorCode SNESNewtonTRSetNormType(SNES snes, NormType norm) 88 { 89 PetscBool flg; 90 91 PetscFunctionBegin; 92 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 93 PetscValidLogicalCollectiveEnum(snes, norm, 2); 94 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 95 if (flg) { 96 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 97 98 tr->norm = norm; 99 } 100 PetscFunctionReturn(PETSC_SUCCESS); 101 } 102 103 /*@ 104 SNESNewtonTRSetQNType - Specify to use a quasi-Newton model. 105 106 Input Parameters: 107 + snes - the nonlinear solver object 108 - use - the type of approximations to be used 109 110 Level: intermediate 111 112 Notes: 113 Options for the approximations can be set with the snes_tr_qn_ and snes_tr_qn_pre_ prefixes. 114 115 .seealso: `SNESNEWTONTR`, `SNESNewtonTRQNType`, `MATLMVM` 116 @*/ 117 PetscErrorCode SNESNewtonTRSetQNType(SNES snes, SNESNewtonTRQNType use) 118 { 119 PetscBool flg; 120 121 PetscFunctionBegin; 122 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 123 PetscValidLogicalCollectiveEnum(snes, use, 2); 124 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 125 if (flg) { 126 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 127 128 tr->qn = use; 129 } 130 PetscFunctionReturn(PETSC_SUCCESS); 131 } 132 133 /*@ 134 SNESNewtonTRSetFallbackType - Set the type of fallback to use if the solution of the trust region subproblem is outside the radius 135 136 Input Parameters: 137 + snes - the nonlinear solver object 138 - ftype - the fallback type, see `SNESNewtonTRFallbackType` 139 140 Level: intermediate 141 142 .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`, 143 `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()` 144 @*/ 145 PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype) 146 { 147 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 148 PetscBool flg; 149 150 PetscFunctionBegin; 151 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 152 PetscValidLogicalCollectiveEnum(snes, ftype, 2); 153 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 154 if (flg) tr->fallback = ftype; 155 PetscFunctionReturn(PETSC_SUCCESS); 156 } 157 158 /*@C 159 SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined. 160 Allows the user a chance to change or override the trust region decision. 161 162 Logically Collective 163 164 Input Parameters: 165 + snes - the nonlinear solver object 166 . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 167 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 168 169 Level: deprecated (since 3.19) 170 171 Note: 172 This function is called BEFORE the function evaluation within the solver. 173 174 .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 175 @*/ 176 PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 177 { 178 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 179 PetscBool flg; 180 181 PetscFunctionBegin; 182 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 183 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 184 if (flg) { 185 if (func) tr->precheck = func; 186 if (ctx) tr->precheckctx = ctx; 187 } 188 PetscFunctionReturn(PETSC_SUCCESS); 189 } 190 191 /*@C 192 SNESNewtonTRGetPreCheck - Gets the pre-check function 193 194 Deprecated use `SNESNEWTONDCTRDC` 195 196 Not Collective 197 198 Input Parameter: 199 . snes - the nonlinear solver context 200 201 Output Parameters: 202 + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()` 203 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 204 205 Level: deprecated (since 3.19) 206 207 .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()` 208 @*/ 209 PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 210 { 211 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 212 PetscBool flg; 213 214 PetscFunctionBegin; 215 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 216 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 217 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 218 if (func) *func = tr->precheck; 219 if (ctx) *ctx = tr->precheckctx; 220 PetscFunctionReturn(PETSC_SUCCESS); 221 } 222 223 /*@C 224 SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 225 function evaluation. Allows the user a chance to change or override the internal decision of the solver 226 227 Logically Collective 228 229 Input Parameters: 230 + snes - the nonlinear solver object 231 . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 232 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 233 234 Level: deprecated (since 3.19) 235 236 Note: 237 This function is called BEFORE the function evaluation within the solver while the function set in 238 `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 239 240 .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()` 241 @*/ 242 PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 243 { 244 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 245 PetscBool flg; 246 247 PetscFunctionBegin; 248 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 249 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 250 if (flg) { 251 if (func) tr->postcheck = func; 252 if (ctx) tr->postcheckctx = ctx; 253 } 254 PetscFunctionReturn(PETSC_SUCCESS); 255 } 256 257 /*@C 258 SNESNewtonTRGetPostCheck - Gets the post-check function 259 260 Not Collective 261 262 Input Parameter: 263 . snes - the nonlinear solver context 264 265 Output Parameters: 266 + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()` 267 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 268 269 Level: intermediate 270 271 .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()` 272 @*/ 273 PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 274 { 275 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 276 PetscBool flg; 277 278 PetscFunctionBegin; 279 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 280 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 281 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 282 if (func) *func = tr->postcheck; 283 if (ctx) *ctx = tr->postcheckctx; 284 PetscFunctionReturn(PETSC_SUCCESS); 285 } 286 287 /*@C 288 SNESNewtonTRPreCheck - Runs the precheck routine 289 290 Logically Collective 291 292 Input Parameters: 293 + snes - the solver 294 . X - The last solution 295 - Y - The step direction 296 297 Output Parameter: 298 . changed_Y - Indicator that the step direction `Y` has been changed. 299 300 Level: intermediate 301 302 .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()` 303 @*/ 304 PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 305 { 306 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 307 PetscBool flg; 308 309 PetscFunctionBegin; 310 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 311 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 312 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 313 *changed_Y = PETSC_FALSE; 314 if (tr->precheck) { 315 PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 316 PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 317 } 318 PetscFunctionReturn(PETSC_SUCCESS); 319 } 320 321 /*@C 322 SNESNewtonTRPostCheck - Runs the postcheck routine 323 324 Logically Collective 325 326 Input Parameters: 327 + snes - the solver 328 . X - The last solution 329 . Y - The full step direction 330 - W - The updated solution, W = X - Y 331 332 Output Parameters: 333 + changed_Y - indicator if step has been changed 334 - changed_W - Indicator if the new candidate solution W has been changed. 335 336 Note: 337 If Y is changed then W is recomputed as X - Y 338 339 Level: intermediate 340 341 .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck()` 342 @*/ 343 PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 344 { 345 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 346 PetscBool flg; 347 348 PetscFunctionBegin; 349 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 350 PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg)); 351 PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name); 352 *changed_Y = PETSC_FALSE; 353 *changed_W = PETSC_FALSE; 354 if (tr->postcheck) { 355 PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 356 PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 357 PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 358 } 359 PetscFunctionReturn(PETSC_SUCCESS); 360 } 361 362 /* stable implementation of roots of a*x^2 + b*x + c = 0 */ 363 static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp) 364 { 365 PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c)); 366 PetscReal x1 = temp / a; 367 PetscReal x2 = c / temp; 368 *xm = PetscMin(x1, x2); 369 *xp = PetscMax(x1, x2); 370 } 371 372 /* Computes the quadratic model difference */ 373 static PetscErrorCode SNESNewtonTRQuadraticDelta(SNES snes, Mat J, PetscBool has_objective, Vec Y, Vec GradF, Vec W, PetscReal *yTHy_, PetscReal *gTy_, PetscReal *deltaqm) 374 { 375 PetscReal yTHy, gTy; 376 377 PetscFunctionBegin; 378 PetscCall(MatMult(J, Y, W)); 379 if (has_objective) PetscCall(VecDotRealPart(Y, W, &yTHy)); 380 else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */ 381 PetscCall(VecDotRealPart(GradF, Y, &gTy)); 382 *deltaqm = -(-(gTy) + 0.5 * (yTHy)); /* difference in quadratic model, -gTy because SNES solves it this way */ 383 if (yTHy_) *yTHy_ = yTHy; 384 if (gTy_) *gTy_ = gTy; 385 PetscFunctionReturn(PETSC_SUCCESS); 386 } 387 388 /* Computes the new objective given X = Xk, Y = direction 389 W work vector, on output W = X - Y 390 G work vector, on output G = SNESFunction(W) */ 391 static PetscErrorCode SNESNewtonTRObjective(SNES snes, PetscBool has_objective, Vec X, Vec Y, Vec W, Vec G, PetscReal *gnorm, PetscReal *fkp1) 392 { 393 PetscBool changed_y, changed_w; 394 395 PetscFunctionBegin; 396 /* TODO: we can add a linesearch here */ 397 PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y)); 398 PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */ 399 PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 400 if (changed_y && !changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X)); 401 402 PetscCall(SNESComputeFunction(snes, W, G)); /* F(Xkp1) = G */ 403 PetscCall(VecNorm(G, NORM_2, gnorm)); 404 if (has_objective) PetscCall(SNESComputeObjective(snes, W, fkp1)); 405 else *fkp1 = 0.5 * PetscSqr(*gnorm); 406 PetscFunctionReturn(PETSC_SUCCESS); 407 } 408 409 static PetscErrorCode SNESSetUpQN_NEWTONTR(SNES snes) 410 { 411 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 412 413 PetscFunctionBegin; 414 PetscCall(MatDestroy(&tr->qnB)); 415 PetscCall(MatDestroy(&tr->qnB_pre)); 416 if (tr->qn) { 417 PetscInt n, N; 418 const char *optionsprefix; 419 Mat B; 420 421 PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 422 PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 423 PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_")); 424 PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 425 PetscCall(MatSetType(B, MATLMVMBFGS)); 426 PetscCall(VecGetLocalSize(snes->vec_sol, &n)); 427 PetscCall(VecGetSize(snes->vec_sol, &N)); 428 PetscCall(MatSetSizes(B, n, n, N, N)); 429 PetscCall(MatSetUp(B)); 430 PetscCall(MatSetFromOptions(B)); 431 PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 432 tr->qnB = B; 433 if (tr->qn == SNES_TR_QN_DIFFERENT) { 434 PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B)); 435 PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 436 PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_pre_")); 437 PetscCall(MatAppendOptionsPrefix(B, optionsprefix)); 438 PetscCall(MatSetType(B, MATLMVMBFGS)); 439 PetscCall(MatSetSizes(B, n, n, N, N)); 440 PetscCall(MatSetUp(B)); 441 PetscCall(MatSetFromOptions(B)); 442 PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func)); 443 tr->qnB_pre = B; 444 } else { 445 PetscCall(PetscObjectReference((PetscObject)tr->qnB)); 446 tr->qnB_pre = tr->qnB; 447 } 448 } 449 PetscFunctionReturn(PETSC_SUCCESS); 450 } 451 452 /* 453 SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 454 (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 455 nonlinear equations 456 457 */ 458 static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) 459 { 460 SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data; 461 Vec X, F, Y, G, W, GradF, YU, Yc; 462 PetscInt maxits, lits; 463 PetscReal rho, fnorm, gnorm = 0.0, xnorm = 0.0, delta, ynorm; 464 PetscReal deltaM, fk, fkp1, deltaqm = 0.0, gTy = 0.0, yTHy = 0.0; 465 PetscReal auk, tauk, gfnorm, gfnorm_k, ycnorm, gTBg, objmin = 0.0, beta_k = 1.0; 466 PC pc; 467 Mat J, Jp; 468 PetscBool already_done = PETSC_FALSE, on_boundary; 469 PetscBool clear_converged_test, rho_satisfied, has_objective; 470 SNES_TR_KSPConverged_Ctx *ctx; 471 void *convctx; 472 473 PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 474 PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *); 475 476 PetscFunctionBegin; 477 PetscCall(SNESGetObjective(snes, &objective, NULL)); 478 has_objective = objective ? PETSC_TRUE : PETSC_FALSE; 479 480 maxits = snes->max_its; /* maximum number of iterations */ 481 X = snes->vec_sol; /* solution vector */ 482 F = snes->vec_func; /* residual vector */ 483 Y = snes->vec_sol_update; /* update vector */ 484 G = snes->work[0]; /* updated residual */ 485 W = snes->work[1]; /* temporary vector */ 486 GradF = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */ 487 YU = snes->work[3]; /* work vector for dogleg method */ 488 Yc = snes->work[4]; /* Cauchy point */ 489 490 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); 491 492 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 493 snes->iter = 0; 494 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 495 496 /* setup QN matrices if needed */ 497 PetscCall(SNESSetUpQN_NEWTONTR(snes)); 498 499 /* Set the linear stopping criteria to use the More' trick if needed */ 500 clear_converged_test = PETSC_FALSE; 501 PetscCall(SNESGetKSP(snes, &snes->ksp)); 502 PetscCall(KSPGetConvergenceTest(snes->ksp, &convtest, &convctx, &convdestroy)); 503 if (convtest != SNESTR_KSPConverged_Private) { 504 clear_converged_test = PETSC_TRUE; 505 PetscCall(PetscNew(&ctx)); 506 ctx->snes = snes; 507 PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 508 PetscCall(KSPSetConvergenceTest(snes->ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy)); 509 PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n")); 510 } 511 512 if (!snes->vec_func_init_set) { 513 PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 514 } else snes->vec_func_init_set = PETSC_FALSE; 515 516 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 517 SNESCheckFunctionNorm(snes, fnorm); 518 PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 519 520 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 521 snes->norm = fnorm; 522 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 523 delta = neP->delta0; 524 deltaM = neP->deltaM; 525 neP->delta = delta; 526 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 527 528 /* test convergence */ 529 rho_satisfied = PETSC_FALSE; 530 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 531 PetscCall(SNESMonitor(snes, 0, fnorm)); 532 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 533 534 if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 535 else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 536 537 /* hook state vector to BFGS preconditioner */ 538 PetscCall(KSPGetPC(snes->ksp, &pc)); 539 PetscCall(PCLMVMSetUpdateVec(pc, X)); 540 541 if (neP->kmdc) PetscCall(KSPSetComputeEigenvalues(snes->ksp, PETSC_TRUE)); 542 543 while (snes->iter < maxits) { 544 /* calculating Jacobian and GradF of minimization function only once */ 545 if (!already_done) { 546 /* Call general purpose update function */ 547 PetscTryTypeMethod(snes, update, snes->iter); 548 549 /* apply the nonlinear preconditioner */ 550 if (snes->npc && snes->npcside == PC_RIGHT) { 551 SNESConvergedReason reason; 552 553 PetscCall(SNESSetInitialFunction(snes->npc, F)); 554 PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 555 PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X)); 556 PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0)); 557 PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 558 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) { 559 snes->reason = SNES_DIVERGED_INNER; 560 PetscFunctionReturn(PETSC_SUCCESS); 561 } 562 // XXX 563 PetscCall(SNESGetNPCFunction(snes, F, &fnorm)); 564 } else if (snes->ops->update) { /* if update is present, recompute objective function and function norm */ 565 PetscCall(SNESComputeFunction(snes, X, F)); 566 } 567 568 /* Jacobian */ 569 J = NULL; 570 Jp = NULL; 571 if (!neP->qnB) { 572 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 573 J = snes->jacobian; 574 Jp = snes->jacobian_pre; 575 } else { /* QN model */ 576 PetscCall(SNESComputeJacobian_MATLMVM(snes, X, neP->qnB, neP->qnB_pre, NULL)); 577 J = neP->qnB; 578 Jp = neP->qnB_pre; 579 } 580 SNESCheckJacobianDomainerror(snes); 581 582 /* objective function */ 583 PetscCall(VecNorm(F, NORM_2, &fnorm)); 584 if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk)); 585 else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */ 586 587 /* GradF */ 588 if (has_objective) gfnorm = fnorm; 589 else { 590 PetscCall(MatMultTranspose(J, F, GradF)); /* grad f = J^T F */ 591 PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); 592 } 593 PetscCall(VecNorm(GradF, neP->norm, &gfnorm_k)); 594 } 595 already_done = PETSC_TRUE; 596 597 /* solve trust-region subproblem */ 598 599 /* first compute Cauchy Point */ 600 PetscCall(MatMult(J, GradF, W)); 601 if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg)); 602 else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */ 603 /* Eqs 4.11 and 4.12 in Nocedal and Wright 2nd Edition (4.7 and 4.8 in 1st Edition) */ 604 auk = delta / gfnorm_k; 605 if (gTBg < 0.0) tauk = 1.0; 606 else tauk = PetscMin(gfnorm * gfnorm * gfnorm_k / (delta * gTBg), 1); 607 auk *= tauk; 608 ycnorm = auk * gfnorm; 609 PetscCall(VecAXPBY(Yc, auk, 0.0, GradF)); 610 611 on_boundary = PETSC_FALSE; 612 if (tauk != 1.0) { 613 KSPConvergedReason reason; 614 615 /* sufficient decrease (see 6.3.27 in Conn, Gould, Toint "Trust Region Methods") 616 beta_k the largest eigenvalue of the Hessian. Here we use the previous estimated value */ 617 objmin = -neP->kmdc * gnorm * PetscMin(gnorm / beta_k, delta); 618 PetscCall(KSPCGSetObjectiveTarget(snes->ksp, objmin)); 619 620 /* specify radius if looking for Newton step and trust region norm is the l2 norm */ 621 PetscCall(KSPCGSetRadius(snes->ksp, neP->fallback == SNES_TR_FALLBACK_NEWTON && neP->norm == NORM_2 ? delta : 0.0)); 622 PetscCall(KSPSetOperators(snes->ksp, J, Jp)); 623 PetscCall(KSPSolve(snes->ksp, F, Y)); 624 SNESCheckKSPSolve(snes); 625 PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 626 PetscCall(KSPGetConvergedReason(snes->ksp, &reason)); 627 on_boundary = (PetscBool)(reason == KSP_CONVERGED_STEP_LENGTH); 628 PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits)); 629 if (neP->kmdc) { /* update estimated Hessian largest eigenvalue */ 630 PetscReal emax, emin; 631 PetscCall(KSPComputeExtremeSingularValues(snes->ksp, &emax, &emin)); 632 if (emax > 0.0) beta_k = emax + 1; 633 } 634 } else { /* Cauchy point is on the boundary, accept it */ 635 on_boundary = PETSC_TRUE; 636 PetscCall(VecCopy(Yc, Y)); 637 PetscCall(PetscInfo(snes, "CP evaluated on boundary. delta: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ycnorm, (double)gTBg)); 638 } 639 PetscCall(VecNorm(Y, neP->norm, &ynorm)); 640 641 /* decide what to do when the update is outside of trust region */ 642 if (ynorm > delta || ynorm == 0.0) { 643 SNESNewtonTRFallbackType fallback = ynorm > 0.0 ? neP->fallback : SNES_TR_FALLBACK_CAUCHY; 644 645 PetscCheck(neP->norm == NORM_2 || fallback != SNES_TR_FALLBACK_DOGLEG, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "DOGLEG without l2 norm not implemented"); 646 switch (fallback) { 647 case SNES_TR_FALLBACK_NEWTON: 648 auk = delta / ynorm; 649 PetscCall(VecScale(Y, auk)); 650 PetscCall(PetscInfo(snes, "SN evaluated. delta: %g, ynorm: %g\n", (double)delta, (double)ynorm)); 651 break; 652 case SNES_TR_FALLBACK_CAUCHY: 653 case SNES_TR_FALLBACK_DOGLEG: 654 if (fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) { 655 PetscCall(VecCopy(Yc, Y)); 656 PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg)); 657 } else { /* take linear combination of Cauchy and Newton direction and step */ 658 auk = gfnorm * gfnorm / gTBg; 659 if (gfnorm_k * auk >= delta) { /* first leg: Cauchy point outside of trust region */ 660 PetscCall(VecAXPBY(Y, delta / gfnorm_k, 0.0, GradF)); 661 PetscCall(PetscInfo(snes, "CP evaluated (outside region). delta: %g, ynorm: %g, ycnorm: %g\n", (double)delta, (double)ynorm, (double)ycnorm)); 662 } else { /* second leg */ 663 PetscReal c0, c1, c2, tau = 0.0, tpos, tneg; 664 PetscBool noroots; 665 666 /* Find solutions of (Eq. 4.16 in Nocedal and Wright) 667 ||p_U + lambda * (p_B - p_U)||^2 - delta^2 = 0, 668 where p_U the Cauchy direction, p_B the Newton direction */ 669 PetscCall(VecAXPBY(YU, auk, 0.0, GradF)); 670 PetscCall(VecAXPY(Y, -1.0, YU)); 671 PetscCall(VecNorm(Y, NORM_2, &c0)); 672 PetscCall(VecDotRealPart(YU, Y, &c1)); 673 c0 = PetscSqr(c0); 674 c2 = PetscSqr(ycnorm) - PetscSqr(delta); 675 PetscQuadraticRoots(c0, 2 * c1, c2, &tneg, &tpos); 676 677 /* In principle the DL strategy as a unique solution in [0,1] 678 here we check that for some reason we numerically failed 679 to compute it. In that case, we use the Cauchy point */ 680 noroots = PetscIsInfOrNanReal(tneg); 681 if (!noroots) { 682 if (tpos > 1) { 683 if (tneg >= 0 && tneg <= 1) { 684 tau = tneg; 685 } else noroots = PETSC_TRUE; 686 } else if (tpos >= 0) { 687 tau = tpos; 688 } else noroots = PETSC_TRUE; 689 } 690 if (noroots) { /* No roots, select Cauchy point */ 691 PetscCall(VecCopy(Yc, Y)); 692 } else { 693 PetscCall(VecAXPBY(Y, 1.0, tau, YU)); 694 } 695 PetscCall(PetscInfo(snes, "%s evaluated. roots: (%g, %g), tau %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", noroots ? "CP" : "DL", (double)tneg, (double)tpos, (double)tau, (double)ynorm, (double)ycnorm, (double)gTBg)); 696 } 697 } 698 break; 699 default: 700 SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode"); 701 break; 702 } 703 } 704 705 /* compute the quadratic model difference */ 706 PetscCall(SNESNewtonTRQuadraticDelta(snes, J, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm)); 707 708 /* Compute new objective function */ 709 PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1)); 710 if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1; 711 else { 712 if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */ 713 else rho = neP->eta1; /* no reduction in quadratic model, step must be rejected */ 714 } 715 716 PetscCall(VecNorm(Y, neP->norm, &ynorm)); 717 PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g, ynormk=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy, (double)ynorm)); 718 719 /* update the size of the trust region */ 720 if (rho < neP->eta2) delta *= neP->t1; /* shrink the region */ 721 else if (rho > neP->eta3 && on_boundary) delta *= neP->t2; /* expand the region */ 722 delta = PetscMin(delta, deltaM); /* but not greater than deltaM */ 723 724 /* log 2-norm of update for moniroting routines */ 725 PetscCall(VecNorm(Y, NORM_2, &ynorm)); 726 727 /* decide on new step */ 728 neP->delta = delta; 729 if (rho > neP->eta1) { 730 rho_satisfied = PETSC_TRUE; 731 } else { 732 rho_satisfied = PETSC_FALSE; 733 PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 734 /* check to see if progress is hopeless */ 735 PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP)); 736 if (!snes->reason) PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 737 if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_TR_DELTA; 738 snes->numFailures++; 739 /* We're not progressing, so return with the current iterate */ 740 if (snes->reason) break; 741 } 742 if (rho_satisfied) { 743 /* Update function values */ 744 already_done = PETSC_FALSE; 745 fnorm = gnorm; 746 fk = fkp1; 747 748 /* New residual and linearization point */ 749 PetscCall(VecCopy(G, F)); 750 PetscCall(VecCopy(W, X)); 751 752 /* Monitor convergence */ 753 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 754 snes->iter++; 755 snes->norm = fnorm; 756 snes->xnorm = xnorm; 757 snes->ynorm = ynorm; 758 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 759 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 760 761 /* Test for convergence, xnorm = || X || */ 762 PetscCall(VecNorm(X, NORM_2, &xnorm)); 763 PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm)); 764 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 765 if (snes->reason) break; 766 } 767 } 768 769 if (clear_converged_test) { 770 PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 771 PetscCall(PetscFree(ctx)); 772 PetscCall(KSPSetConvergenceTest(snes->ksp, convtest, convctx, convdestroy)); 773 } 774 PetscFunctionReturn(PETSC_SUCCESS); 775 } 776 777 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) 778 { 779 PetscFunctionBegin; 780 PetscCall(SNESSetWorkVecs(snes, 5)); 781 PetscCall(SNESSetUpMatrices(snes)); 782 PetscFunctionReturn(PETSC_SUCCESS); 783 } 784 785 static PetscErrorCode SNESReset_NEWTONTR(SNES snes) 786 { 787 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 788 789 PetscFunctionBegin; 790 PetscCall(MatDestroy(&tr->qnB)); 791 PetscCall(MatDestroy(&tr->qnB_pre)); 792 PetscFunctionReturn(PETSC_SUCCESS); 793 } 794 795 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) 796 { 797 PetscFunctionBegin; 798 PetscCall(SNESReset_NEWTONTR(snes)); 799 PetscCall(PetscFree(snes->data)); 800 PetscFunctionReturn(PETSC_SUCCESS); 801 } 802 803 static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) 804 { 805 SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data; 806 SNESNewtonTRQNType qn; 807 SNESNewtonTRFallbackType fallback; 808 NormType norm; 809 PetscReal deltatol; 810 PetscBool flg; 811 812 PetscFunctionBegin; 813 PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 814 PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 815 PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 816 PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 817 PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 818 PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 819 PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 820 PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 821 PetscCall(PetscOptionsReal("-snes_tr_kmdc", "sufficient decrease parameter", "None", ctx->kmdc, &ctx->kmdc, NULL)); 822 823 deltatol = snes->deltatol; 824 PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", deltatol, &deltatol, &flg)); 825 if (flg) PetscCall(SNESSetTrustRegionTolerance(snes, deltatol)); 826 827 fallback = ctx->fallback; 828 PetscCall(PetscOptionsEnum("-snes_tr_fallback_type", "Type of fallback if subproblem solution is outside of the trust region", "SNESNewtonTRSetFallbackType", SNESNewtonTRFallbackTypes, (PetscEnum)fallback, (PetscEnum *)&fallback, &flg)); 829 if (flg) PetscCall(SNESNewtonTRSetFallbackType(snes, fallback)); 830 831 qn = ctx->qn; 832 PetscCall(PetscOptionsEnum("-snes_tr_qn", "Use Quasi-Newton approximations for the model", "SNESNewtonTRSetQNType", SNESNewtonTRQNTypes, (PetscEnum)qn, (PetscEnum *)&qn, &flg)); 833 if (flg) PetscCall(SNESNewtonTRSetQNType(snes, qn)); 834 835 norm = ctx->norm; 836 PetscCall(PetscOptionsEnum("-snes_tr_norm_type", "Type of norm for trust region bounds", "SNESNewtonTRSetNormType", NormTypes, (PetscEnum)norm, (PetscEnum *)&norm, &flg)); 837 if (flg) PetscCall(SNESNewtonTRSetNormType(snes, norm)); 838 839 PetscOptionsHeadEnd(); 840 PetscFunctionReturn(PETSC_SUCCESS); 841 } 842 843 static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) 844 { 845 SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data; 846 PetscBool iascii; 847 848 PetscFunctionBegin; 849 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 850 if (iascii) { 851 PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)snes->deltatol)); 852 PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 853 PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM)); 854 PetscCall(PetscViewerASCIIPrintf(viewer, " kmdc=%g\n", (double)tr->kmdc)); 855 PetscCall(PetscViewerASCIIPrintf(viewer, " fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback])); 856 if (tr->qn) PetscCall(PetscViewerASCIIPrintf(viewer, " qn=%s\n", SNESNewtonTRQNTypes[tr->qn])); 857 if (tr->norm != NORM_2) PetscCall(PetscViewerASCIIPrintf(viewer, " norm=%s\n", NormTypes[tr->norm])); 858 } 859 PetscFunctionReturn(PETSC_SUCCESS); 860 } 861 862 /*MC 863 SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction {cite}`nocedal2006numerical` 864 865 Options Database Keys: 866 + -snes_tr_tol <tol> - trust region tolerance 867 . -snes_tr_eta1 <eta1> - trust region parameter eta1 <= eta2, rho > eta1 breaks out of the inner iteration (default: eta1=0.001) 868 . -snes_tr_eta2 <eta2> - trust region parameter, rho <= eta2 shrinks the trust region (default: eta2=0.25) 869 . -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 870 . -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 871 . -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 872 . -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL) 873 . -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2) 874 - -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. 875 Level: deprecated 876 877 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, 878 `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, 879 `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()`, `SNESNewtonTRSetQNType()` 880 M*/ 881 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) 882 { 883 SNES_NEWTONTR *neP; 884 885 PetscFunctionBegin; 886 snes->ops->setup = SNESSetUp_NEWTONTR; 887 snes->ops->solve = SNESSolve_NEWTONTR; 888 snes->ops->reset = SNESReset_NEWTONTR; 889 snes->ops->destroy = SNESDestroy_NEWTONTR; 890 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR; 891 snes->ops->view = SNESView_NEWTONTR; 892 893 snes->stol = 0.0; 894 snes->usesksp = PETSC_TRUE; 895 snes->npcside = PC_RIGHT; 896 snes->usesnpc = PETSC_TRUE; 897 898 snes->alwayscomputesfinalresidual = PETSC_TRUE; 899 900 PetscCall(PetscNew(&neP)); 901 snes->data = (void *)neP; 902 neP->delta = 0.0; 903 neP->delta0 = 0.2; 904 neP->eta1 = 0.001; 905 neP->eta2 = 0.25; 906 neP->eta3 = 0.75; 907 neP->t1 = 0.25; 908 neP->t2 = 2.0; 909 neP->deltaM = 1.e10; 910 neP->norm = NORM_2; 911 neP->fallback = SNES_TR_FALLBACK_NEWTON; 912 neP->kmdc = 0.0; /* by default do not use sufficient decrease */ 913 PetscFunctionReturn(PETSC_SUCCESS); 914 } 915