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