1 #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I "petscsnes.h" I*/ 2 3 typedef struct { 4 SNES snes; 5 /* Information on the regular SNES convergence test; which may have been user provided 6 Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho 7 Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private 8 */ 9 10 PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 11 PetscErrorCode (*convdestroy)(void *); 12 void *convctx; 13 } SNES_TRDC_KSPConverged_Ctx; 14 15 static PetscErrorCode SNESNewtonTRSetTolerances_TRDC(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0) 16 { 17 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 18 19 PetscFunctionBegin; 20 if (delta_min == PETSC_DETERMINE) delta_min = 1.e-12; 21 if (delta_max == PETSC_DETERMINE) delta_max = 0.5; 22 if (delta_0 == PETSC_DETERMINE) delta_0 = 0.1; 23 if (delta_min != PETSC_CURRENT) tr->deltatol = delta_min; 24 if (delta_max != PETSC_CURRENT) tr->deltaM = delta_max; 25 if (delta_0 != PETSC_CURRENT) tr->delta0 = delta_0; 26 PetscFunctionReturn(PETSC_SUCCESS); 27 } 28 29 static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) 30 { 31 SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx; 32 SNES snes = ctx->snes; 33 SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)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, NORM_2, &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 SNESTRDC_KSPConverged_Destroy(void *cctx) 51 { 52 SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx; 53 54 PetscFunctionBegin; 55 PetscCall((*ctx->convdestroy)(ctx->convctx)); 56 PetscCall(PetscFree(ctx)); 57 PetscFunctionReturn(PETSC_SUCCESS); 58 } 59 60 /* 61 SNESTRDC_Converged_Private -test convergence JUST for 62 the trust region tolerance. 63 64 */ 65 static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) 66 { 67 SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 68 69 PetscFunctionBegin; 70 *reason = SNES_CONVERGED_ITERATING; 71 if (neP->delta < xnorm * neP->deltatol) { 72 PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)neP->deltatol)); 73 *reason = SNES_DIVERGED_TR_DELTA; 74 } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) { 75 PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs)); 76 *reason = SNES_DIVERGED_FUNCTION_COUNT; 77 } 78 PetscFunctionReturn(PETSC_SUCCESS); 79 } 80 81 /*@ 82 SNESNewtonTRDCGetRhoFlag - Get whether the current solution update is within the trust-region. 83 84 Logically Collective 85 86 Input Parameter: 87 . snes - the nonlinear solver object 88 89 Output Parameter: 90 . rho_flag - `PETSC_FALSE` or `PETSC_TRUE` 91 92 Level: developer 93 94 .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPreCheck()`, 95 `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()` 96 @*/ 97 PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag) 98 { 99 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 100 101 PetscFunctionBegin; 102 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 103 PetscAssertPointer(rho_flag, 2); 104 *rho_flag = tr->rho_satisfied; 105 PetscFunctionReturn(PETSC_SUCCESS); 106 } 107 108 /*@C 109 SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined. 110 Allows the user a chance to change or override the trust region decision. 111 112 Logically Collective 113 114 Input Parameters: 115 + snes - the nonlinear solver object 116 . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()` 117 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 118 119 Level: intermediate 120 121 Note: 122 This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver. 123 124 .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 125 `SNESNewtonTRDCGetRhoFlag()` 126 @*/ 127 PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) 128 { 129 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 130 131 PetscFunctionBegin; 132 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 133 if (func) tr->precheck = func; 134 if (ctx) tr->precheckctx = ctx; 135 PetscFunctionReturn(PETSC_SUCCESS); 136 } 137 138 /*@C 139 SNESNewtonTRDCGetPreCheck - Gets the pre-check function optionally set with `SNESNewtonTRDCSetPreCheck()` 140 141 Not Collective 142 143 Input Parameter: 144 . snes - the nonlinear solver context 145 146 Output Parameters: 147 + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPreCheck()` 148 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 149 150 Level: intermediate 151 152 .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()` 153 @*/ 154 PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) 155 { 156 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 157 158 PetscFunctionBegin; 159 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 160 if (func) *func = tr->precheck; 161 if (ctx) *ctx = tr->precheckctx; 162 PetscFunctionReturn(PETSC_SUCCESS); 163 } 164 165 /*@C 166 SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next 167 function evaluation. Allows the user a chance to change or override the decision of the line search routine 168 169 Logically Collective 170 171 Input Parameters: 172 + snes - the nonlinear solver object 173 . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()` 174 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 175 176 Level: intermediate 177 178 Note: 179 This function is called BEFORE the function evaluation within the `SNESNEWTONTRDC` solver while the function set in 180 `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation. 181 182 .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 183 @*/ 184 PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) 185 { 186 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 187 188 PetscFunctionBegin; 189 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 190 if (func) tr->postcheck = func; 191 if (ctx) tr->postcheckctx = ctx; 192 PetscFunctionReturn(PETSC_SUCCESS); 193 } 194 195 /*@C 196 SNESNewtonTRDCGetPostCheck - Gets the post-check function optionally set with `SNESNewtonTRDCSetPostCheck()` 197 198 Not Collective 199 200 Input Parameter: 201 . snes - the nonlinear solver context 202 203 Output Parameters: 204 + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRDCPostCheck()` 205 - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 206 207 Level: intermediate 208 209 .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()` 210 @*/ 211 PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) 212 { 213 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 214 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 217 if (func) *func = tr->postcheck; 218 if (ctx) *ctx = tr->postcheckctx; 219 PetscFunctionReturn(PETSC_SUCCESS); 220 } 221 222 // PetscClangLinter pragma disable: -fdoc-internal-linkage 223 /*@C 224 SNESNewtonTRDCPreCheck - Called before the step has been determined in `SNESNEWTONTRDC` 225 226 Logically Collective 227 228 Input Parameters: 229 + snes - the solver 230 . X - The last solution 231 - Y - The step direction 232 233 Output Parameter: 234 . changed_Y - Indicator that the step direction `Y` has been changed. 235 236 Level: developer 237 238 .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCPostCheck()` 239 @*/ 240 static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) 241 { 242 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 243 244 PetscFunctionBegin; 245 *changed_Y = PETSC_FALSE; 246 if (tr->precheck) { 247 PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx)); 248 PetscValidLogicalCollectiveBool(snes, *changed_Y, 4); 249 } 250 PetscFunctionReturn(PETSC_SUCCESS); 251 } 252 253 // PetscClangLinter pragma disable: -fdoc-internal-linkage 254 /*@C 255 SNESNewtonTRDCPostCheck - Called after the step has been determined in `SNESNEWTONTRDC` but before the function evaluation at that step 256 257 Logically Collective 258 259 Input Parameters: 260 + snes - the solver 261 . X - The last solution 262 . Y - The full step direction 263 - W - The updated solution, W = X - Y 264 265 Output Parameters: 266 + changed_Y - indicator if step has been changed 267 - changed_W - Indicator if the new candidate solution `W` has been changed. 268 269 Level: developer 270 271 Note: 272 If `Y` is changed then `W` is recomputed as `X` - `Y` 273 274 .seealso: [](ch_snes), `SNES`, `SNESNEWTONTRDC`, `SNESNEWTONTRDC`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, `SNESNewtonTRDCPreCheck() 275 @*/ 276 static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) 277 { 278 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 279 280 PetscFunctionBegin; 281 *changed_Y = PETSC_FALSE; 282 *changed_W = PETSC_FALSE; 283 if (tr->postcheck) { 284 PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx)); 285 PetscValidLogicalCollectiveBool(snes, *changed_Y, 5); 286 PetscValidLogicalCollectiveBool(snes, *changed_W, 6); 287 } 288 PetscFunctionReturn(PETSC_SUCCESS); 289 } 290 291 /* 292 SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy 293 (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of 294 nonlinear equations 295 296 */ 297 static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) 298 { 299 SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data; 300 Vec X, F, Y, G, W, GradF, YNtmp; 301 Vec YCtmp; 302 Mat jac; 303 PetscInt maxits, i, j, lits, inner_count, bs; 304 PetscReal rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */ 305 PetscReal inorms[99]; /* need to make it dynamic eventually, fixed max block size of 99 for now */ 306 PetscReal deltaM, ynnorm, f0, mp, gTy, g, yTHy; /* rho calculation */ 307 PetscReal auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg; /* Cauchy Point */ 308 KSP ksp; 309 SNESConvergedReason reason = SNES_CONVERGED_ITERATING; 310 PetscBool breakout = PETSC_FALSE; 311 SNES_TRDC_KSPConverged_Ctx *ctx; 312 PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *); 313 void *convctx; 314 315 PetscFunctionBegin; 316 maxits = snes->max_its; /* maximum number of iterations */ 317 X = snes->vec_sol; /* solution vector */ 318 F = snes->vec_func; /* residual vector */ 319 Y = snes->work[0]; /* update vector */ 320 G = snes->work[1]; /* updated residual */ 321 W = snes->work[2]; /* temporary vector */ 322 GradF = snes->work[3]; /* grad f = J^T F */ 323 YNtmp = snes->work[4]; /* Newton solution */ 324 YCtmp = snes->work[5]; /* Cauchy solution */ 325 326 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); 327 328 PetscCall(VecGetBlockSize(YNtmp, &bs)); 329 330 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 331 snes->iter = 0; 332 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 333 334 /* Set the linear stopping criteria to use the More' trick. From tr.c */ 335 PetscCall(SNESGetKSP(snes, &ksp)); 336 PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy)); 337 if (convtest != SNESTRDC_KSPConverged_Private) { 338 PetscCall(PetscNew(&ctx)); 339 ctx->snes = snes; 340 PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 341 PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy)); 342 PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n")); 343 } 344 345 if (!snes->vec_func_init_set) { 346 PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */ 347 } else snes->vec_func_init_set = PETSC_FALSE; 348 349 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */ 350 SNESCheckFunctionNorm(snes, fnorm); 351 PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */ 352 353 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 354 snes->norm = fnorm; 355 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 356 delta = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */ 357 deltaM = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */ 358 neP->delta = delta; 359 PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 360 PetscCall(SNESMonitor(snes, 0, fnorm)); 361 362 neP->rho_satisfied = PETSC_FALSE; 363 364 /* test convergence */ 365 PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 366 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 367 368 for (i = 0; i < maxits; i++) { 369 PetscBool changed_y; 370 PetscBool changed_w; 371 372 /* dogleg method */ 373 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 374 SNESCheckJacobianDomainerror(snes); 375 PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian)); 376 PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */ 377 SNESCheckKSPSolve(snes); /* this is necessary but old tr.c did not have it*/ 378 PetscCall(KSPGetIterationNumber(snes->ksp, &lits)); 379 PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL)); 380 381 /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable) 382 for inner iteration and Cauchy direction calculation 383 */ 384 if (bs > 1 && neP->auto_scale_multiphase) { 385 PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms)); 386 for (j = 0; j < bs; j++) { 387 if (neP->auto_scale_max > 1.0) { 388 if (inorms[j] < 1.0 / neP->auto_scale_max) inorms[j] = 1.0 / neP->auto_scale_max; 389 } 390 PetscCall(VecStrideSet(W, j, inorms[j])); 391 PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j])); 392 PetscCall(VecStrideScale(X, j, 1.0 / inorms[j])); 393 } 394 PetscCall(VecNorm(X, NORM_2, &xnorm)); 395 if (i == 0) { 396 delta = neP->delta0 * xnorm; 397 } else { 398 delta = neP->delta * xnorm; 399 } 400 deltaM = neP->deltaM * xnorm; 401 PetscCall(MatDiagonalScale(jac, NULL, W)); 402 } 403 404 /* calculating GradF of minimization function */ 405 PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */ 406 PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */ 407 408 inner_count = 0; 409 neP->rho_satisfied = PETSC_FALSE; 410 while (1) { 411 if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */ 412 PetscCall(VecCopy(YNtmp, Y)); 413 } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */ 414 PetscCall(MatMult(jac, GradF, W)); 415 PetscCall(VecDotRealPart(W, W, &gTBg)); /* completes GradF^T J^T J GradF */ 416 PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */ 417 if (gTBg <= 0.0) { 418 auk = PETSC_MAX_REAL; 419 } else { 420 auk = PetscSqr(gfnorm) / gTBg; 421 } 422 auk = PetscMin(delta / gfnorm, auk); 423 PetscCall(VecCopy(GradF, YCtmp)); /* this could be improved */ 424 PetscCall(VecScale(YCtmp, auk)); /* YCtmp, Cauchy solution*/ 425 PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */ 426 if (ycnorm >= delta) { /* see if the Cauchy solution meets the criteria */ 427 PetscCall(VecCopy(YCtmp, Y)); 428 PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm)); 429 } else { /* take ratio, tau, of Cauchy and Newton direction and step */ 430 PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */ 431 PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */ 432 c0 = PetscSqr(c0); 433 PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1)); 434 c1 = 2.0 * c1; 435 PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */ 436 c2 = PetscSqr(c2) - PetscSqr(delta); 437 tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */ 438 tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); 439 tau = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */ 440 PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm)); 441 PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp)); 442 PetscCall(VecAXPY(W, -tau, YCtmp)); 443 PetscCall(VecCopy(W, Y)); /* this could be improved */ 444 } 445 } else { 446 /* if Cauchy is disabled, only use Newton direction */ 447 auk = delta / ynnorm; 448 PetscCall(VecScale(YNtmp, auk)); 449 PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/ 450 } 451 452 PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm */ 453 f0 = 0.5 * PetscSqr(fnorm); /* minimizing function f(X) */ 454 PetscCall(MatMult(jac, Y, W)); 455 PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */ 456 PetscCall(VecDotRealPart(GradF, Y, &gTy)); 457 mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/ 458 459 /* scale back solution update */ 460 if (bs > 1 && neP->auto_scale_multiphase) { 461 for (j = 0; j < bs; j++) { 462 PetscCall(VecStrideScale(Y, j, inorms[j])); 463 if (inner_count == 0) { 464 /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */ 465 /* need to scale back X to match Y and provide proper update to the external code */ 466 PetscCall(VecStrideScale(X, j, inorms[j])); 467 } 468 } 469 if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */ 470 PetscCall(VecNorm(Y, NORM_2, &temp_ynorm)); 471 } else { 472 temp_xnorm = xnorm; 473 temp_ynorm = ynorm; 474 } 475 inner_count++; 476 477 /* Evaluate the solution to meet the improvement ratio criteria */ 478 PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y)); 479 PetscCall(VecWAXPY(W, -1.0, Y, X)); 480 PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w)); 481 if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X)); 482 PetscCall(VecCopy(Y, snes->vec_sol_update)); 483 PetscCall(SNESComputeFunction(snes, W, G)); /* F(X-Y) = G */ 484 PetscCall(VecNorm(G, NORM_2, &gnorm)); /* gnorm <- || g || */ 485 SNESCheckFunctionNorm(snes, gnorm); 486 g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */ 487 if (f0 == mp) rho = 0.0; 488 else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */ 489 490 if (rho < neP->eta2) { 491 delta *= neP->t1; /* shrink the region */ 492 } else if (rho > neP->eta3) { 493 delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */ 494 } 495 496 neP->delta = delta; 497 if (rho >= neP->eta1) { 498 /* unscale delta and xnorm before going to the next outer iteration */ 499 if (bs > 1 && neP->auto_scale_multiphase) { 500 neP->delta = delta / xnorm; 501 xnorm = temp_xnorm; 502 ynorm = temp_ynorm; 503 } 504 neP->rho_satisfied = PETSC_TRUE; 505 break; /* the improvement ratio is satisfactory */ 506 } 507 PetscCall(PetscInfo(snes, "Trying again in smaller region\n")); 508 509 /* check to see if progress is hopeless */ 510 neP->itflag = PETSC_FALSE; 511 /* both delta, ynorm, and xnorm are either scaled or unscaled */ 512 PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP)); 513 /* if multiphase state changes, break out inner iteration */ 514 if (reason == SNES_BREAKOUT_INNER_ITER) { 515 if (bs > 1 && neP->auto_scale_multiphase) { 516 /* unscale delta and xnorm before going to the next outer iteration */ 517 neP->delta = delta / xnorm; 518 xnorm = temp_xnorm; 519 ynorm = temp_ynorm; 520 } 521 reason = SNES_CONVERGED_ITERATING; 522 break; 523 } 524 if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER; 525 if (reason) { 526 if (reason < 0) { 527 /* We're not progressing, so return with the current iterate */ 528 PetscCall(SNESMonitor(snes, i + 1, fnorm)); 529 breakout = PETSC_TRUE; 530 break; 531 } else if (reason > 0) { 532 /* We're converged, so return with the current iterate and update solution */ 533 PetscCall(SNESMonitor(snes, i + 1, fnorm)); 534 breakout = PETSC_FALSE; 535 break; 536 } 537 } 538 snes->numFailures++; 539 } 540 if (!breakout) { 541 /* Update function and solution vectors */ 542 fnorm = gnorm; 543 PetscCall(VecCopy(G, F)); 544 PetscCall(VecCopy(W, X)); 545 /* Monitor convergence */ 546 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 547 snes->iter = i + 1; 548 snes->norm = fnorm; 549 snes->xnorm = xnorm; 550 snes->ynorm = ynorm; 551 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 552 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits)); 553 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 554 /* Test for convergence, xnorm = || X || */ 555 neP->itflag = PETSC_TRUE; 556 if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm)); 557 PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP); 558 if (reason) break; 559 } else break; 560 } 561 562 /* PetscCall(PetscFree(inorms)); */ 563 if (i == maxits) { 564 PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits)); 565 if (!reason) reason = SNES_DIVERGED_MAX_IT; 566 } 567 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 568 snes->reason = reason; 569 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 570 if (convtest != SNESTRDC_KSPConverged_Private) { 571 PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy)); 572 PetscCall(PetscFree(ctx)); 573 PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy)); 574 } 575 PetscFunctionReturn(PETSC_SUCCESS); 576 } 577 578 static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) 579 { 580 PetscFunctionBegin; 581 PetscCall(SNESSetWorkVecs(snes, 6)); 582 PetscCall(SNESSetUpMatrices(snes)); 583 PetscFunctionReturn(PETSC_SUCCESS); 584 } 585 586 static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) 587 { 588 PetscFunctionBegin; 589 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", NULL)); 590 PetscCall(PetscFree(snes->data)); 591 PetscFunctionReturn(PETSC_SUCCESS); 592 } 593 594 static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems PetscOptionsObject) 595 { 596 SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data; 597 598 PetscFunctionBegin; 599 PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations"); 600 PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESNewtonTRSetTolerances", ctx->deltatol, &ctx->deltatol, NULL)); 601 PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL)); 602 PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL)); 603 PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL)); 604 PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL)); 605 PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL)); 606 PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL)); 607 PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL)); 608 PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL)); 609 PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL)); 610 PetscCall(PetscOptionsBool("-snes_trdc_auto_scale_multiphase", "auto_scale_multiphase", "Auto scaling for proper cauchy direction", ctx->auto_scale_multiphase, &ctx->auto_scale_multiphase, NULL)); 611 PetscOptionsHeadEnd(); 612 PetscFunctionReturn(PETSC_SUCCESS); 613 } 614 615 static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer) 616 { 617 SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data; 618 PetscBool iascii; 619 620 PetscFunctionBegin; 621 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 622 if (iascii) { 623 PetscCall(PetscViewerASCIIPrintf(viewer, " Trust region tolerance %g\n", (double)tr->deltatol)); 624 PetscCall(PetscViewerASCIIPrintf(viewer, " eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3)); 625 PetscCall(PetscViewerASCIIPrintf(viewer, " delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM)); 626 } 627 PetscFunctionReturn(PETSC_SUCCESS); 628 } 629 630 /*MC 631 SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction 632 633 Options Database Keys: 634 + -snes_trdc_tol <tol> - trust region tolerance 635 . -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001) 636 . -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25) 637 . -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75) 638 . -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25) 639 . -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0) 640 . -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, $deltaM*norm2(x)$ (default: 0.5) 641 . -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, $delta0*norm2(x)$ (default: 0.1) 642 . -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor 643 . -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm 644 - -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region 645 646 Level: advanced 647 648 Notes: 649 `SNESNEWTONTRDC` only works for root-finding problems and does not support objective functions. 650 The main difference with respect to `SNESNEWTONTR` is that `SNESNEWTONTRDC` scales the trust region by the norm of the current linearization point. 651 Future version may extend the `SNESNEWTONTR` code and deprecate `SNESNEWTONTRDC`. 652 653 For details, see {cite}`park2021linear` 654 655 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNewtonTRSetTolerances()`, 656 `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`, 657 `SNESNewtonTRDCGetRhoFlag()`, `SNESNewtonTRDCSetPreCheck()` 658 M*/ 659 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) 660 { 661 SNES_NEWTONTRDC *neP; 662 663 PetscFunctionBegin; 664 snes->ops->setup = SNESSetUp_NEWTONTRDC; 665 snes->ops->solve = SNESSolve_NEWTONTRDC; 666 snes->ops->destroy = SNESDestroy_NEWTONTRDC; 667 snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC; 668 snes->ops->view = SNESView_NEWTONTRDC; 669 670 snes->usesksp = PETSC_TRUE; 671 snes->usesnpc = PETSC_FALSE; 672 673 snes->alwayscomputesfinalresidual = PETSC_TRUE; 674 675 PetscCall(SNESParametersInitialize(snes)); 676 677 PetscCall(PetscNew(&neP)); 678 snes->data = (void *)neP; 679 neP->eta1 = 0.001; 680 neP->eta2 = 0.25; 681 neP->eta3 = 0.75; 682 neP->t1 = 0.25; 683 neP->t2 = 2.0; 684 neP->sigma = 0.0001; 685 neP->itflag = PETSC_FALSE; 686 neP->rnorm0 = 0.0; 687 neP->ttol = 0.0; 688 neP->use_cauchy = PETSC_TRUE; 689 neP->auto_scale_multiphase = PETSC_FALSE; 690 neP->auto_scale_max = -1.0; 691 neP->rho_satisfied = PETSC_FALSE; 692 neP->delta = 0.0; 693 neP->deltaM = 0.5; 694 neP->delta0 = 0.1; 695 neP->deltatol = 1.e-12; 696 697 /* for multiphase (multivariable) scaling */ 698 /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13 699 on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now. 700 PetscCall(VecGetBlockSize(snes->work[0],&neP->bs)); 701 PetscCall(PetscCalloc1(neP->bs,&neP->inorms)); 702 */ 703 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", SNESNewtonTRSetTolerances_TRDC)); 704 PetscFunctionReturn(PETSC_SUCCESS); 705 } 706