Actual source code: tr.c

  1: #include <../src/snes/impls/tr/trimpl.h>

  3: typedef struct {
  4:   SNES                  snes;
  5:   KSPConvergenceTestFn *convtest;
  6:   PetscCtxDestroyFn    *convdestroy;
  7:   void                 *convctx;
  8: } SNES_TR_KSPConverged_Ctx;

 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};

 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;

 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: }

 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;

 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: }

 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: }

 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;

 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: }

 76: static PetscErrorCode SNESTR_KSPConverged_Destroy(PetscCtxRt cctx)
 77: {
 78:   SNES_TR_KSPConverged_Ctx *ctx = *(SNES_TR_KSPConverged_Ctx **)cctx;

 80:   PetscFunctionBegin;
 81:   PetscCall((*ctx->convdestroy)(&ctx->convctx));
 82:   PetscCall(PetscFree(ctx));
 83:   PetscFunctionReturn(PETSC_SUCCESS);
 84: }

 86: static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscReal xnorm, PetscReal ynorm, PetscReal fnorm)
 87: {
 88:   SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data;

 90:   PetscFunctionBegin;
 91:   /* only check convergence for iter > 0 */
 92:   if (snes->iter) PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
 93:   if (snes->reason > 0) PetscFunctionReturn(PETSC_SUCCESS);
 94:   if (neP->delta < neP->deltam) {
 95:     PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)neP->deltam));
 96:     snes->reason = SNES_DIVERGED_TR_DELTA;
 97:   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
 98:     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
 99:     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
100:   }
101:   PetscFunctionReturn(PETSC_SUCCESS);
102: }

104: /*@
105:   SNESNewtonTRSetNormType - Specify the type of norm to use for the computation of the trust region.

107:   Input Parameters:
108: + snes - the nonlinear solver object
109: - norm - the norm type

111:   Level: intermediate

113: .seealso: `SNESNEWTONTR`, `NormType`
114: @*/
115: PetscErrorCode SNESNewtonTRSetNormType(SNES snes, NormType norm)
116: {
117:   PetscBool flg;

119:   PetscFunctionBegin;
122:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
123:   if (flg) {
124:     SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;

126:     tr->norm = norm;
127:   }
128:   PetscFunctionReturn(PETSC_SUCCESS);
129: }

131: /*@
132:   SNESNewtonTRSetQNType - Specify to use a quasi-Newton model.

134:   Input Parameters:
135: + snes - the nonlinear solver object
136: - use  - the type of approximations to be used

138:   Level: intermediate

140:   Notes:
141:   Options for the approximations can be set with the `snes_tr_qn_` and `snes_tr_qn_pre_` prefixes.

143: .seealso: `SNESNEWTONTR`, `SNESNewtonTRQNType`, `MATLMVM`
144: @*/
145: PetscErrorCode SNESNewtonTRSetQNType(SNES snes, SNESNewtonTRQNType use)
146: {
147:   PetscBool flg;

149:   PetscFunctionBegin;
152:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
153:   if (flg) {
154:     SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;

156:     tr->qn = use;
157:   }
158:   PetscFunctionReturn(PETSC_SUCCESS);
159: }

161: /*@
162:   SNESNewtonTRSetFallbackType - Set the type of fallback to use if the solution of the trust region subproblem is outside the radius

164:   Input Parameters:
165: + snes  - the nonlinear solver object
166: - ftype - the fallback type, see `SNESNewtonTRFallbackType`

168:   Level: intermediate

170: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`,
171:           `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`
172: @*/
173: PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype)
174: {
175:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
176:   PetscBool      flg;

178:   PetscFunctionBegin;
181:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
182:   if (flg) tr->fallback = ftype;
183:   PetscFunctionReturn(PETSC_SUCCESS);
184: }

186: /*@C
187:   SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined.
188:   Allows the user a chance to change or override the trust region decision.

190:   Logically Collective

192:   Input Parameters:
193: + snes - the nonlinear solver object
194: . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()`
195: - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

197:   Calling sequence of `func`:
198: + snes    - the nonlinear solver object
199: . X       - the current solution value
200: . Y       - the tentative update step
201: . changed - output, flag indicating `Y` has been changed by the pre-check
202: - ctx     - the optional application context

204:   Level: intermediate

206:   Note:
207:   This function is called BEFORE the function evaluation within the solver.

209: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`
210: @*/
211: PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES snes, Vec X, Vec Y, PetscBool *changed, PetscCtx ctx), PetscCtx ctx)
212: {
213:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
214:   PetscBool      flg;

216:   PetscFunctionBegin;
218:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
219:   if (flg) {
220:     if (func) tr->precheck = func;
221:     if (ctx) tr->precheckctx = ctx;
222:   }
223:   PetscFunctionReturn(PETSC_SUCCESS);
224: }

226: /*@C
227:   SNESNewtonTRGetPreCheck - Gets the pre-check function

229:   Not Collective

231:   Input Parameter:
232: . snes - the nonlinear solver context

234:   Output Parameters:
235: + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPreCheck()`
236: - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

238:   Calling sequence of `func`:
239: + snes    - the nonlinear solver object
240: . X       - the current solution value
241: . Y       - the tentative update step
242: . changed - output, flag indicating `Y` has been changed by the pre-check
243: - ctx     - the optional application context

245:   Level: intermediate

247: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()`
248: @*/
249: PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES snes, Vec X, Vec Y, PetscBool *changed, PetscCtx ctx), PetscCtxRt ctx)
250: {
251:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
252:   PetscBool      flg;

254:   PetscFunctionBegin;
256:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
257:   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
258:   if (func) *func = tr->precheck;
259:   if (ctx) *(void **)ctx = tr->precheckctx;
260:   PetscFunctionReturn(PETSC_SUCCESS);
261: }

263: /*@C
264:   SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
265:   function evaluation. Allows the user a chance to change or override the internal decision of the solver

267:   Logically Collective

269:   Input Parameters:
270: + snes - the nonlinear solver object
271: . func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()`
272: - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

274:   Calling sequence of `func`:
275: + snes      - the nonlinear solver object
276: . X         - the current solution value
277: . Y         - the tentative update step
278: . W         - the tentative new solution value
279: . changed_Y - output, flag indicated `Y` has been changed by the post-check
280: . changed_W - output, flag indicated `W` has been changed by the post-check
281: - ctx       - the optional application context

283:   Level: intermediate

285:   Note:
286:   This function is called BEFORE the function evaluation within the solver while the function set in
287:   `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation.

289: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`
290: @*/
291: PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W, PetscCtx ctx), PetscCtx ctx)
292: {
293:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
294:   PetscBool      flg;

296:   PetscFunctionBegin;
298:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
299:   if (flg) {
300:     if (func) tr->postcheck = func;
301:     if (ctx) tr->postcheckctx = ctx;
302:   }
303:   PetscFunctionReturn(PETSC_SUCCESS);
304: }

306: /*@C
307:   SNESNewtonTRGetPostCheck - Gets the post-check function

309:   Not Collective

311:   Input Parameter:
312: . snes - the nonlinear solver context

314:   Output Parameters:
315: + func - [optional] function evaluation routine, for the calling sequence see `SNESNewtonTRPostCheck()`
316: - ctx  - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`)

318:   Calling sequence of `func`:
319: + snes      - the nonlinear solver object
320: . X         - the current solution value
321: . Y         - the tentative update step
322: . W         - the tentative new solution value
323: . changed_Y - output, flag indicated `Y` has been changed by the post-check
324: . changed_W - output, flag indicated `W` has been changed by the post-check
325: - ctx       - the optional application context

327:   Level: intermediate

329: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()`
330: @*/
331: PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W, PetscCtx ctx), PetscCtxRt ctx)
332: {
333:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
334:   PetscBool      flg;

336:   PetscFunctionBegin;
338:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
339:   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
340:   if (func) *func = tr->postcheck;
341:   if (ctx) *(void **)ctx = tr->postcheckctx;
342:   PetscFunctionReturn(PETSC_SUCCESS);
343: }

345: /*@C
346:   SNESNewtonTRPreCheck - Runs the precheck routine

348:   Logically Collective

350:   Input Parameters:
351: + snes - the solver
352: . X    - The last solution
353: - Y    - The step direction

355:   Output Parameter:
356: . changed_Y - Indicator that the step direction `Y` has been changed.

358:   Level: intermediate

360: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()`
361: @*/
362: PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y)
363: {
364:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
365:   PetscBool      flg;

367:   PetscFunctionBegin;
369:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
370:   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
371:   *changed_Y = PETSC_FALSE;
372:   if (tr->precheck) {
373:     PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
375:   }
376:   PetscFunctionReturn(PETSC_SUCCESS);
377: }

379: /*@C
380:   SNESNewtonTRPostCheck - Runs the postcheck routine

382:   Logically Collective

384:   Input Parameters:
385: + snes - the solver
386: . X    - The last solution
387: . Y    - The full step direction
388: - W    - The updated solution, W = X - Y

390:   Output Parameters:
391: + changed_Y - indicator if step has been changed
392: - changed_W - Indicator if the new candidate solution W has been changed.

394:   Note:
395:   If Y is changed then W is recomputed as X - Y

397:   Level: intermediate

399: .seealso: [](ch_snes), `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck()`
400: @*/
401: PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
402: {
403:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
404:   PetscBool      flg;

406:   PetscFunctionBegin;
408:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
409:   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
410:   *changed_Y = PETSC_FALSE;
411:   *changed_W = PETSC_FALSE;
412:   if (tr->postcheck) {
413:     PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
416:   }
417:   PetscFunctionReturn(PETSC_SUCCESS);
418: }

420: /* stable implementation of roots of a*x^2 + b*x + c = 0 */
421: static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp)
422: {
423:   PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c));
424:   PetscReal x1   = temp / a;
425:   PetscReal x2   = c / temp;
426:   *xm            = PetscMin(x1, x2);
427:   *xp            = PetscMax(x1, x2);
428: }

430: /* Computes the quadratic model difference */
431: static PetscErrorCode SNESNewtonTRQuadraticDelta(SNES snes, Mat J, PetscBool has_objective, Vec Y, Vec GradF, Vec W, PetscReal *yTHy_, PetscReal *gTy_, PetscReal *deltaqm)
432: {
433:   PetscReal yTHy, gTy;

435:   PetscFunctionBegin;
436:   PetscCall(MatMult(J, Y, W));
437:   if (has_objective) PetscCall(VecDotRealPart(Y, W, &yTHy));
438:   else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */
439:   PetscCall(VecDotRealPart(GradF, Y, &gTy));
440:   *deltaqm = -(-(gTy) + 0.5 * (yTHy)); /* difference in quadratic model, -gTy because SNES solves it this way */
441:   if (yTHy_) *yTHy_ = yTHy;
442:   if (gTy_) *gTy_ = gTy;
443:   PetscFunctionReturn(PETSC_SUCCESS);
444: }

446: /* Computes the new objective given X = Xk, Y = direction
447:    W work vector, on output W = X - Y
448:    G work vector, on output G = SNESFunction(W) */
449: static PetscErrorCode SNESNewtonTRObjective(SNES snes, PetscBool has_objective, Vec X, Vec Y, Vec W, Vec G, PetscReal *gnorm, PetscReal *fkp1)
450: {
451:   PetscBool changed_Y, changed_W;

453:   PetscFunctionBegin;
454:   /* TODO: we can add a linesearch here */
455:   PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_Y));
456:   PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */
457:   PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_Y, &changed_W));
458:   if (changed_Y && !changed_W) PetscCall(VecWAXPY(W, -1.0, Y, X));

460:   PetscCall(SNESComputeFunction(snes, W, G)); /*  F(Xkp1) = G */
461:   PetscCall(VecNorm(G, NORM_2, gnorm));
462:   SNESCheckFunctionDomainError(snes, *gnorm);
463:   if (has_objective) {
464:     PetscCall(SNESComputeObjective(snes, W, fkp1));
465:     SNESCheckObjectiveDomainError(snes, *fkp1);
466:   } else *fkp1 = 0.5 * PetscSqr(*gnorm);
467:   PetscFunctionReturn(PETSC_SUCCESS);
468: }

470: static PetscErrorCode SNESSetUpQN_NEWTONTR(SNES snes)
471: {
472:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;

474:   PetscFunctionBegin;
475:   PetscCall(MatDestroy(&tr->qnB));
476:   PetscCall(MatDestroy(&tr->qnB_pre));
477:   if (tr->qn) {
478:     PetscInt    n, N;
479:     const char *optionsprefix;
480:     Mat         B;

482:     PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B));
483:     PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
484:     PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_"));
485:     PetscCall(MatAppendOptionsPrefix(B, optionsprefix));
486:     PetscCall(MatSetType(B, MATLMVMBFGS));
487:     PetscCall(VecGetLocalSize(snes->vec_sol, &n));
488:     PetscCall(VecGetSize(snes->vec_sol, &N));
489:     PetscCall(MatSetSizes(B, n, n, N, N));
490:     PetscCall(MatSetUp(B));
491:     PetscCall(MatSetFromOptions(B));
492:     PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func));
493:     tr->qnB = B;
494:     if (tr->qn == SNES_TR_QN_DIFFERENT) {
495:       PetscCall(MatCreate(PetscObjectComm((PetscObject)snes), &B));
496:       PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
497:       PetscCall(MatSetOptionsPrefix(B, "snes_tr_qn_pre_"));
498:       PetscCall(MatAppendOptionsPrefix(B, optionsprefix));
499:       PetscCall(MatSetType(B, MATLMVMBFGS));
500:       PetscCall(MatSetSizes(B, n, n, N, N));
501:       PetscCall(MatSetUp(B));
502:       PetscCall(MatSetFromOptions(B));
503:       PetscCall(MatLMVMAllocate(B, snes->vec_sol, snes->vec_func));
504:       tr->qnB_pre = B;
505:     } else {
506:       PetscCall(PetscObjectReference((PetscObject)tr->qnB));
507:       tr->qnB_pre = tr->qnB;
508:     }
509:   }
510:   PetscFunctionReturn(PETSC_SUCCESS);
511: }

513: /*
514:    SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
515:    (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
516:    nonlinear equations

518: */
519: static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
520: {
521:   SNES_NEWTONTR            *neP = (SNES_NEWTONTR *)snes->data;
522:   Vec                       X, F, Y, G, W, GradF, YU, Yc;
523:   PetscInt                  maxits, lits;
524:   PetscReal                 rho, fnorm, gnorm = 0.0, xnorm = 0.0, delta, ynorm;
525:   PetscReal                 fk, fkp1, deltaqm = 0.0, gTy = 0.0, yTHy = 0.0;
526:   PetscReal                 auk, tauk, gfnorm, gfnorm_k, ycnorm, gTBg, objmin = 0.0, beta_k = 1.0;
527:   PC                        pc;
528:   Mat                       J, Jp;
529:   PetscBool                 already_done = PETSC_FALSE, on_boundary, use_cauchy;
530:   PetscBool                 clear_converged_test, step_ok, has_objective;
531:   SNES_TR_KSPConverged_Ctx *ctx;
532:   void                     *convctx;
533:   SNESObjectiveFn          *objective;
534:   KSPConvergenceTestFn     *convtest;
535:   PetscCtxDestroyFn        *convdestroy;

537:   PetscFunctionBegin;
538:   PetscCall(SNESGetObjective(snes, &objective, NULL));
539:   has_objective = objective ? PETSC_TRUE : PETSC_FALSE;

541:   maxits = snes->max_its;                                   /* maximum number of iterations */
542:   X      = snes->vec_sol;                                   /* solution vector */
543:   F      = snes->vec_func;                                  /* residual vector */
544:   Y      = snes->vec_sol_update;                            /* update vector */
545:   G      = snes->work[0];                                   /* updated residual */
546:   W      = snes->work[1];                                   /* temporary vector */
547:   GradF  = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */
548:   YU     = snes->work[3];                                   /* work vector for dogleg method */
549:   Yc     = snes->work[4];                                   /* Cauchy point */

551:   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);

553:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
554:   snes->iter = 0;
555:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));

557:   /* setup QN matrices if needed */
558:   PetscCall(SNESSetUpQN_NEWTONTR(snes));

560:   /* Set the linear stopping criteria to use the More' trick if needed */
561:   clear_converged_test = PETSC_FALSE;
562:   PetscCall(SNESGetKSP(snes, &snes->ksp));
563:   PetscCall(KSPGetConvergenceTest(snes->ksp, &convtest, &convctx, &convdestroy));
564:   if (convtest != SNESTR_KSPConverged_Private) {
565:     clear_converged_test = PETSC_TRUE;
566:     PetscCall(PetscNew(&ctx));
567:     ctx->snes = snes;
568:     PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
569:     PetscCall(KSPSetConvergenceTest(snes->ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy));
570:     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n"));
571:   }

573:   if (!snes->vec_func_init_set) {
574:     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
575:   } else snes->vec_func_init_set = PETSC_FALSE;

577:   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
578:   SNESCheckFunctionDomainError(snes, fnorm);
579:   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */

581:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
582:   snes->norm = fnorm;
583:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
584:   delta      = neP->delta0;
585:   neP->delta = delta;
586:   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));

588:   /* test convergence */
589:   step_ok = PETSC_FALSE;
590:   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
591:   PetscCall(SNESMonitor(snes, 0, fnorm));
592:   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);

594:   if (has_objective) {
595:     PetscCall(SNESComputeObjective(snes, X, &fk));
596:     SNESCheckObjectiveDomainError(snes, fk);
597:   } else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */

599:   /* hook state vector to BFGS preconditioner */
600:   PetscCall(KSPGetPC(snes->ksp, &pc));
601:   PetscCall(PCLMVMSetUpdateVec(pc, X));

603:   if (neP->kmdc) PetscCall(KSPSetComputeEigenvalues(snes->ksp, PETSC_TRUE));

605:   while (snes->iter < maxits) {
606:     /* calculating Jacobian and GradF of minimization function only once */
607:     if (!already_done) {
608:       /* Call general purpose update function */
609:       PetscTryTypeMethod(snes, update, snes->iter);

611:       /* apply the nonlinear preconditioner */
612:       if (snes->npc && snes->npcside == PC_RIGHT) {
613:         SNESConvergedReason reason;

615:         PetscCall(SNESSetInitialFunction(snes->npc, F));
616:         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
617:         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
618:         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
619:         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
620:         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
621:           snes->reason = SNES_DIVERGED_INNER;
622:           PetscFunctionReturn(PETSC_SUCCESS);
623:         }
624:         // XXX
625:         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
626:       }

628:       /* Jacobian */
629:       J  = NULL;
630:       Jp = NULL;
631:       if (!neP->qnB) {
632:         PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
633:         SNESCheckJacobianDomainError(snes);
634:         J  = snes->jacobian;
635:         Jp = snes->jacobian_pre;
636:       } else { /* QN model */
637:         PetscCall(SNESComputeJacobian_MATLMVM(snes, X, neP->qnB, neP->qnB_pre, NULL));
638:         J  = neP->qnB;
639:         Jp = neP->qnB_pre;
640:       }
641:       SNESCheckJacobianDomainError(snes);

643:       /* objective function */
644:       PetscCall(VecNorm(F, NORM_2, &fnorm));
645:       SNESCheckFunctionDomainError(snes, fnorm);
646:       if (has_objective) {
647:         PetscCall(SNESComputeObjective(snes, X, &fk));
648:         SNESCheckObjectiveDomainError(snes, fk);
649:       } else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */

651:       /* GradF */
652:       if (has_objective) gfnorm = fnorm;
653:       else {
654:         PetscCall(MatMultTranspose(J, F, GradF)); /* grad f = J^T F */
655:         PetscCall(VecNorm(GradF, NORM_2, &gfnorm));
656:       }
657:       PetscCall(VecNorm(GradF, neP->norm, &gfnorm_k));
658:     }
659:     already_done = PETSC_TRUE;

661:     /* solve trust-region subproblem */

663:     /* first compute Cauchy Point */
664:     PetscCall(MatMult(J, GradF, W));
665:     if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg));
666:     else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */
667:     /* Eqs 4.11 and 4.12 in Nocedal and Wright 2nd Edition (4.7 and 4.8 in 1st Edition) */
668:     auk = delta / gfnorm_k;
669:     if (gTBg < 0.0) tauk = 1.0;
670:     else tauk = PetscMin(gfnorm * gfnorm * gfnorm_k / (delta * gTBg), 1);
671:     auk *= tauk;
672:     ycnorm = auk * gfnorm;
673:     PetscCall(VecAXPBY(Yc, auk, 0.0, GradF));

675:     on_boundary = PETSC_FALSE;
676:     use_cauchy  = (PetscBool)(tauk == 1.0 && has_objective);
677:     if (!use_cauchy) {
678:       KSPConvergedReason reason;

680:       /* sufficient decrease (see 6.3.27 in Conn, Gould, Toint "Trust Region Methods")
681:          beta_k the largest eigenvalue of the Hessian. Here we use the previous estimated value */
682:       objmin = -neP->kmdc * gnorm * PetscMin(gnorm / beta_k, delta);
683:       PetscCall(KSPCGSetObjectiveTarget(snes->ksp, objmin));

685:       /* specify radius if looking for Newton step and trust region norm is the l2 norm */
686:       PetscCall(KSPCGSetRadius(snes->ksp, neP->fallback == SNES_TR_FALLBACK_NEWTON && neP->norm == NORM_2 ? delta : 0.0));
687:       PetscCall(KSPSetOperators(snes->ksp, J, Jp));
688:       PetscCall(KSPSolve(snes->ksp, F, Y));
689:       SNESCheckKSPSolve(snes);
690:       PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
691:       PetscCall(KSPGetConvergedReason(snes->ksp, &reason));
692:       on_boundary = (PetscBool)(reason == KSP_CONVERGED_STEP_LENGTH);
693:       PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
694:       if (neP->kmdc) { /* update estimated Hessian largest eigenvalue */
695:         PetscReal emax, emin;
696:         PetscCall(KSPComputeExtremeSingularValues(snes->ksp, &emax, &emin));
697:         if (emax > 0.0) beta_k = emax + 1;
698:       }
699:     } else { /* Cauchy point is on the boundary, accept it */
700:       on_boundary = PETSC_TRUE;
701:       PetscCall(VecCopy(Yc, Y));
702:       PetscCall(PetscInfo(snes, "CP evaluated on boundary. delta: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ycnorm, (double)gTBg));
703:     }
704:     PetscCall(VecNorm(Y, neP->norm, &ynorm));

706:     /* decide what to do when the update is outside of trust region */
707:     if (!use_cauchy && (ynorm > delta || ynorm == 0.0)) {
708:       SNESNewtonTRFallbackType fallback = ynorm > 0.0 ? neP->fallback : SNES_TR_FALLBACK_CAUCHY;

710:       PetscCheck(neP->norm == NORM_2 || fallback != SNES_TR_FALLBACK_DOGLEG, PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "DOGLEG without l2 norm not implemented");
711:       switch (fallback) {
712:       case SNES_TR_FALLBACK_NEWTON:
713:         auk = delta / ynorm;
714:         PetscCall(VecScale(Y, auk));
715:         PetscCall(PetscInfo(snes, "SN evaluated. delta: %g, ynorm: %g\n", (double)delta, (double)ynorm));
716:         break;
717:       case SNES_TR_FALLBACK_CAUCHY:
718:       case SNES_TR_FALLBACK_DOGLEG:
719:         if (fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) {
720:           PetscCall(VecCopy(Yc, Y));
721:           PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg));
722:         } else { /* take linear combination of Cauchy and Newton direction and step */
723:           auk = gfnorm * gfnorm / gTBg;
724:           if (gfnorm_k * auk >= delta) { /* first leg: Cauchy point outside of trust region */
725:             PetscCall(VecAXPBY(Y, delta / gfnorm_k, 0.0, GradF));
726:             PetscCall(PetscInfo(snes, "CP evaluated (outside region). delta: %g, ynorm: %g, ycnorm: %g\n", (double)delta, (double)ynorm, (double)ycnorm));
727:           } else { /* second leg */
728:             PetscReal c0, c1, c2, tau = 0.0, tpos, tneg;
729:             PetscBool noroots;

731:             /* Find solutions of (Eq. 4.16 in Nocedal and Wright)
732:                  ||p_U + lambda * (p_B - p_U)||^2 - delta^2 = 0,
733:                where p_U  the Cauchy direction, p_B the Newton direction */
734:             PetscCall(VecAXPBY(YU, auk, 0.0, GradF));
735:             PetscCall(VecAXPY(Y, -1.0, YU));
736:             PetscCall(VecNorm(Y, NORM_2, &c0));
737:             PetscCall(VecDotRealPart(YU, Y, &c1));
738:             c0 = PetscSqr(c0);
739:             c2 = PetscSqr(ycnorm) - PetscSqr(delta);
740:             PetscQuadraticRoots(c0, 2 * c1, c2, &tneg, &tpos);

742:             /* In principle the DL strategy as a unique solution in [0,1]
743:                here we check that for some reason we numerically failed
744:                to compute it. In that case, we use the Cauchy point */
745:             noroots = PetscIsInfOrNanReal(tneg);
746:             if (!noroots) {
747:               if (tpos > 1) {
748:                 if (tneg >= 0 && tneg <= 1) {
749:                   tau = tneg;
750:                 } else noroots = PETSC_TRUE;
751:               } else if (tpos >= 0) {
752:                 tau = tpos;
753:               } else noroots = PETSC_TRUE;
754:             }
755:             if (noroots) { /* No roots, select Cauchy point */
756:               PetscCall(VecCopy(Yc, Y));
757:             } else {
758:               PetscCall(VecAXPBY(Y, 1.0, tau, YU));
759:             }
760:             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));
761:           }
762:         }
763:         break;
764:       default:
765:         SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode");
766:         break;
767:       }
768:     }

770:     /* compute the quadratic model difference */
771:     PetscCall(SNESNewtonTRQuadraticDelta(snes, J, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm));

773:     /* Compute new objective function */
774:     PetscCall(SNESNewtonTRObjective(snes, has_objective, X, Y, W, G, &gnorm, &fkp1));
775:     if (PetscIsInfOrNanReal(fkp1)) rho = neP->eta1;
776:     else {
777:       if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */
778:       else rho = neP->eta1;                           /*  no reduction in quadratic model, step must be rejected */
779:     }

781:     PetscCall(VecNorm(Y, neP->norm, &ynorm));
782:     PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g, ynormk=%g, gnorm=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy, (double)ynorm, (double)gnorm));

784:     /* update the size of the trust region */
785:     if (rho < neP->eta2) delta *= neP->t1;                     /* shrink the region */
786:     else if (rho > neP->eta3 && on_boundary) delta *= neP->t2; /* expand the region */
787:     delta = PetscMin(delta, neP->deltaM);                      /* but not greater than deltaM */

789:     /* log 2-norm of update for monitoring routines */
790:     PetscCall(VecNorm(Y, NORM_2, &ynorm));

792:     /* decide on new step */
793:     neP->delta = delta;
794:     if (rho > neP->eta1) {
795:       step_ok = PETSC_TRUE;
796:     } else {
797:       step_ok = PETSC_FALSE;
798:       /* check convergence */
799:       PetscCall(SNESTR_Converged_Private(snes, xnorm, ynorm, gnorm));
800:       if (snes->reason < 0) {
801:         snes->numFailures++;
802:         break;
803:       } else if (!snes->reason) {
804:         snes->numFailures++;
805:         PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
806:       } else step_ok = PETSC_TRUE;
807:     }
808:     if (step_ok) {
809:       /* Update function values */
810:       already_done = PETSC_FALSE;
811:       fnorm        = gnorm;
812:       fk           = fkp1;

814:       /* New residual and linearization point */
815:       PetscCall(VecCopy(G, F));
816:       PetscCall(VecCopy(W, X));

818:       /* Monitor convergence */
819:       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
820:       snes->iter++;
821:       snes->norm  = fnorm;
822:       snes->xnorm = xnorm;
823:       snes->ynorm = ynorm;
824:       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
825:       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));

827:       /* Test for convergence, xnorm = || X || */
828:       PetscCall(VecNorm(X, NORM_2, &xnorm));
829:       PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
830:       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
831:       if (snes->reason) break;
832:     }
833:   }

835:   if (clear_converged_test) {
836:     PetscCall(KSPGetAndClearConvergenceTest(snes->ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
837:     PetscCall(PetscFree(ctx));
838:     PetscCall(KSPSetConvergenceTest(snes->ksp, convtest, convctx, convdestroy));
839:   }
840:   PetscFunctionReturn(PETSC_SUCCESS);
841: }

843: static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
844: {
845:   PetscFunctionBegin;
846:   PetscCall(SNESSetWorkVecs(snes, 5));
847:   PetscCall(SNESSetUpMatrices(snes));
848:   PetscFunctionReturn(PETSC_SUCCESS);
849: }

851: static PetscErrorCode SNESReset_NEWTONTR(SNES snes)
852: {
853:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;

855:   PetscFunctionBegin;
856:   PetscCall(MatDestroy(&tr->qnB));
857:   PetscCall(MatDestroy(&tr->qnB_pre));
858:   PetscFunctionReturn(PETSC_SUCCESS);
859: }

861: static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
862: {
863:   PetscFunctionBegin;
864:   PetscCall(SNESReset_NEWTONTR(snes));
865:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", NULL));
866:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", NULL));
867:   PetscCall(PetscFree(snes->data));
868:   PetscFunctionReturn(PETSC_SUCCESS);
869: }

871: static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems PetscOptionsObject)
872: {
873:   SNES_NEWTONTR           *ctx = (SNES_NEWTONTR *)snes->data;
874:   SNESNewtonTRQNType       qn;
875:   SNESNewtonTRFallbackType fallback;
876:   NormType                 norm;
877:   PetscBool                flg;

879:   PetscFunctionBegin;
880:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
881:   PetscCall(PetscOptionsDeprecated("-snes_tr_deltaM", "-snes_tr_deltamax", "3.22", NULL));
882:   PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "SNESNewtonTRSetUpdateParameters", ctx->eta1, &ctx->eta1, NULL));
883:   PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "SNESNewtonTRSetUpdateParameters", ctx->eta2, &ctx->eta2, NULL));
884:   PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "SNESNewtonTRSetUpdateParameters", ctx->eta3, &ctx->eta3, NULL));
885:   PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "SNESNewtonTRSetUpdateParameters", ctx->t1, &ctx->t1, NULL));
886:   PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "SNESNewtonTRSetUpdateParameters", ctx->t2, &ctx->t2, NULL));
887:   PetscCall(PetscOptionsReal("-snes_tr_delta0", "Initial trust region size", "SNESNewtonTRSetTolerances", ctx->delta0, &ctx->delta0, NULL));
888:   PetscCall(PetscOptionsReal("-snes_tr_deltamin", "Minimum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltam, &ctx->deltam, NULL));
889:   PetscCall(PetscOptionsReal("-snes_tr_deltamax", "Maximum allowed trust region size", "SNESNewtonTRSetTolerances", ctx->deltaM, &ctx->deltaM, NULL));
890:   PetscCall(PetscOptionsReal("-snes_tr_kmdc", "sufficient decrease parameter", "None", ctx->kmdc, &ctx->kmdc, NULL));

892:   fallback = ctx->fallback;
893:   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));
894:   if (flg) PetscCall(SNESNewtonTRSetFallbackType(snes, fallback));

896:   qn = ctx->qn;
897:   PetscCall(PetscOptionsEnum("-snes_tr_qn", "Use Quasi-Newton approximations for the model", "SNESNewtonTRSetQNType", SNESNewtonTRQNTypes, (PetscEnum)qn, (PetscEnum *)&qn, &flg));
898:   if (flg) PetscCall(SNESNewtonTRSetQNType(snes, qn));

900:   norm = ctx->norm;
901:   PetscCall(PetscOptionsEnum("-snes_tr_norm_type", "Type of norm for trust region bounds", "SNESNewtonTRSetNormType", NormTypes, (PetscEnum)norm, (PetscEnum *)&norm, &flg));
902:   if (flg) PetscCall(SNESNewtonTRSetNormType(snes, norm));

904:   PetscOptionsHeadEnd();
905:   PetscFunctionReturn(PETSC_SUCCESS);
906: }

908: static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer)
909: {
910:   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
911:   PetscBool      isascii;

913:   PetscFunctionBegin;
914:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
915:   if (isascii) {
916:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Trust region parameters:\n"));
917:     PetscCall(PetscViewerASCIIPrintf(viewer, "    eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
918:     PetscCall(PetscViewerASCIIPrintf(viewer, "    t1=%g, t2=%g\n", (double)tr->t1, (double)tr->t2));
919:     PetscCall(PetscViewerASCIIPrintf(viewer, "    delta_min=%g, delta_0=%g, delta_max=%g\n", (double)tr->deltam, (double)tr->delta0, (double)tr->deltaM));
920:     PetscCall(PetscViewerASCIIPrintf(viewer, "    kmdc=%g\n", (double)tr->kmdc));
921:     PetscCall(PetscViewerASCIIPrintf(viewer, "    fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback]));
922:     if (tr->qn) PetscCall(PetscViewerASCIIPrintf(viewer, "    qn=%s\n", SNESNewtonTRQNTypes[tr->qn]));
923:     if (tr->norm != NORM_2) PetscCall(PetscViewerASCIIPrintf(viewer, "    norm=%s\n", NormTypes[tr->norm]));
924:   }
925:   PetscFunctionReturn(PETSC_SUCCESS);
926: }

928: /*@
929:   SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.

931:   Logically Collective

933:   Input Parameters:
934: + snes - the `SNES` context
935: - tol  - tolerance

937:   Level: deprecated

939: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetTolerances()`
940: @*/
941: PetscErrorCode SNESSetTrustRegionTolerance(SNES snes, PetscReal tol)
942: {
943:   return SNESNewtonTRSetTolerances(snes, tol, PETSC_CURRENT, PETSC_CURRENT);
944: }

946: /*@
947:   SNESNewtonTRSetTolerances - Sets the trust region parameter tolerances.

949:   Logically Collective

951:   Input Parameters:
952: + snes      - the `SNES` context
953: . delta_min - minimum allowed trust region size
954: . delta_max - maximum allowed trust region size
955: - delta_0   - initial trust region size

957:   Options Database Key:
958: + -snes_tr_deltamin tol - Set minimum size
959: . -snes_tr_deltamax tol - Set maximum size
960: - -snes_tr_delta0   tol - Set initial size

962:   Note:
963:   Use `PETSC_DETERMINE` to use the default value for the given `SNES`.
964:   Use `PETSC_CURRENT` to retain a value.

966:   Fortran Note:
967:   Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL`

969:   Level: intermediate

971: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRGetTolerances()`
972: @*/
973: PetscErrorCode SNESNewtonTRSetTolerances(SNES snes, PetscReal delta_min, PetscReal delta_max, PetscReal delta_0)
974: {
975:   PetscFunctionBegin;
980:   PetscTryMethod(snes, "SNESNewtonTRSetTolerances_C", (SNES, PetscReal, PetscReal, PetscReal), (snes, delta_min, delta_max, delta_0));
981:   PetscFunctionReturn(PETSC_SUCCESS);
982: }

984: /*@
985:   SNESNewtonTRGetTolerances - Gets the trust region parameter tolerances.

987:   Not Collective

989:   Input Parameter:
990: . snes - the `SNES` context

992:   Output Parameters:
993: + delta_min - minimum allowed trust region size or `NULL`
994: . delta_max - maximum allowed trust region size or `NULL`
995: - delta_0   - initial trust region size or `NULL`

997:   Level: intermediate

999: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetTolerances()`
1000: @*/
1001: PetscErrorCode SNESNewtonTRGetTolerances(SNES snes, PetscReal *delta_min, PetscReal *delta_max, PetscReal *delta_0)
1002: {
1003:   PetscFunctionBegin;
1005:   if (delta_min) PetscAssertPointer(delta_min, 2);
1006:   if (delta_max) PetscAssertPointer(delta_max, 3);
1007:   if (delta_0) PetscAssertPointer(delta_0, 4);
1008:   PetscUseMethod(snes, "SNESNewtonTRGetTolerances_C", (SNES, PetscReal *, PetscReal *, PetscReal *), (snes, delta_min, delta_max, delta_0));
1009:   PetscFunctionReturn(PETSC_SUCCESS);
1010: }

1012: /*@
1013:   SNESNewtonTRSetUpdateParameters - Sets the trust region update parameters.

1015:   Logically Collective

1017:   Input Parameters:
1018: + snes - the `SNES` context
1019: . eta1 - acceptance tolerance
1020: . eta2 - shrinking tolerance
1021: . eta3 - enlarging tolerance
1022: . t1   - shrink factor
1023: - t2   - enlarge factor

1025:   Options Database Key:
1026: + -snes_tr_eta1 tol - Set `eta1`
1027: . -snes_tr_eta2 tol - Set `eta2`
1028: . -snes_tr_eta3 tol - Set `eta3`
1029: . -snes_tr_t1   tol - Set `t1`
1030: - -snes_tr_t2   tol - Set `t2`

1032:   Notes:
1033:   Given the ratio $\rho = \frac{f(x_k) - f(x_k+s_k)}{m(0) - m(s_k)}$, with $x_k$ the current iterate,
1034:   $s_k$ the computed step, $f$ the objective function, and $m$ the quadratic model, the trust region
1035:   radius is modified as follows

1037:   $
1038:   \delta =
1039:   \begin{cases}
1040:   \delta * t_1 ,& \rho < \eta_2 \\
1041:   \delta * t_2 ,& \rho > \eta_3 \\
1042:   \end{cases}
1043:   $

1045:   The step is accepted if $\rho > \eta_1$.
1046:   Use `PETSC_DETERMINE` to use the default value for the given `SNES`.
1047:   Use `PETSC_CURRENT` to retain a value.

1049:   Fortran Note:
1050:   Use `PETSC_DETERMINE_REAL`, `PETSC_CURRENT_REAL`

1052:   Level: intermediate

1054: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESSetObjective()`, `SNESNewtonTRGetUpdateParameters()`
1055: @*/
1056: PetscErrorCode SNESNewtonTRSetUpdateParameters(SNES snes, PetscReal eta1, PetscReal eta2, PetscReal eta3, PetscReal t1, PetscReal t2)
1057: {
1058:   PetscBool flg;

1060:   PetscFunctionBegin;
1067:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
1068:   if (flg) {
1069:     SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;

1071:     if (eta1 == PETSC_DETERMINE) eta1 = tr->default_eta1;
1072:     if (eta2 == PETSC_DETERMINE) eta2 = tr->default_eta2;
1073:     if (eta3 == PETSC_DETERMINE) eta3 = tr->default_eta3;
1074:     if (t1 == PETSC_DETERMINE) t1 = tr->default_t1;
1075:     if (t2 == PETSC_DETERMINE) t2 = tr->default_t2;
1076:     if (eta1 != PETSC_CURRENT) tr->eta1 = eta1;
1077:     if (eta2 != PETSC_CURRENT) tr->eta2 = eta2;
1078:     if (eta3 != PETSC_CURRENT) tr->eta3 = eta3;
1079:     if (t1 != PETSC_CURRENT) tr->t1 = t1;
1080:     if (t2 != PETSC_CURRENT) tr->t2 = t2;
1081:   }
1082:   PetscFunctionReturn(PETSC_SUCCESS);
1083: }

1085: /*@
1086:   SNESNewtonTRGetUpdateParameters - Gets the trust region update parameters.

1088:   Not Collective

1090:   Input Parameter:
1091: . snes - the `SNES` context

1093:   Output Parameters:
1094: + eta1 - acceptance tolerance
1095: . eta2 - shrinking tolerance
1096: . eta3 - enlarging tolerance
1097: . t1   - shrink factor
1098: - t2   - enlarge factor

1100:   Level: intermediate

1102: .seealso: [](ch_snes), `SNES`, `SNESNEWTONTR`, `SNESNewtonTRSetUpdateParameters()`
1103: @*/
1104: PetscErrorCode SNESNewtonTRGetUpdateParameters(SNES snes, PetscReal *eta1, PetscReal *eta2, PetscReal *eta3, PetscReal *t1, PetscReal *t2)
1105: {
1106:   SNES_NEWTONTR *tr;
1107:   PetscBool      flg;

1109:   PetscFunctionBegin;
1111:   if (eta1) PetscAssertPointer(eta1, 2);
1112:   if (eta2) PetscAssertPointer(eta2, 3);
1113:   if (eta3) PetscAssertPointer(eta3, 4);
1114:   if (t1) PetscAssertPointer(t1, 5);
1115:   if (t2) PetscAssertPointer(t2, 6);
1116:   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
1117:   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
1118:   tr = (SNES_NEWTONTR *)snes->data;
1119:   if (eta1) *eta1 = tr->eta1;
1120:   if (eta2) *eta2 = tr->eta2;
1121:   if (eta3) *eta3 = tr->eta3;
1122:   if (t1) *t1 = tr->t1;
1123:   if (t2) *t2 = tr->t2;
1124:   PetscFunctionReturn(PETSC_SUCCESS);
1125: }

1127: /*MC
1128:    SNESNEWTONTR - Newton based nonlinear solver that uses a trust-region strategy

1130:    Options Database Keys:
1131: +  -snes_tr_deltamin deltamin                    - trust region parameter, minimum size of trust region
1132: .  -snes_tr_deltamax deltamax                    - trust region parameter, max size of trust region (default: 1e10)
1133: .  -snes_tr_delta0 delta0                        - trust region parameter, initial size of trust region (default: 0.2)
1134: .  -snes_tr_eta1 eta1                            - trust region parameter $eta1 \le eta2$, $rho > eta1$ breaks out of the inner iteration (default: 0.001)
1135: .  -snes_tr_eta2 eta2                            - trust region parameter, $rho \le eta2$ shrinks the trust region (default: 0.25)
1136: .  -snes_tr_eta3 eta3                            - trust region parameter $eta3 > eta2$, $rho \ge eta3$ expands the trust region (default: 0.75)
1137: .  -snes_tr_t1 t1                                - trust region parameter, shrinking factor of trust region (default: 0.25)
1138: .  -snes_tr_t2 t2                                - trust region parameter, expanding factor of trust region (default: 2.0)
1139: .  -snes_tr_norm_type (1|2|infinity)             - Type of norm for trust region bounds (default: 2)
1140: -  -snes_tr_fallback_type (newton,cauchy,dogleg) - Solution strategy to test reduction when step is outside of trust region.
1141:                                                    Can use scaled Newton direction, Cauchy point (Steepest Descent direction) or dogleg method.

1143:    Level: beginner

1145:    Notes:
1146:    The code is largely based on the book {cite}`nocedal2006numerical` and supports minimizing objective functions using a quadratic model.
1147:    Quasi-Newton models are also supported.

1149:    Default step computation uses the Newton direction, but a dogleg type update is also supported.
1150:    The 1- and infinity-norms are also supported, via `SNESNewtonTRSetNormType()`, when computing the trust region bounds.

1152: .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESSetObjective()`,
1153:           `SNESNewtonTRSetTolerances()`, `SNESNewtonTRSetUpdateParameters()`,
1154:           `SNESNewtonTRSetNormType()`, `SNESNewtonTRSetFallbackType()`, `SNESNewtonTRSetQNType()`,
1155:           `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRSetPreCheck()`
1156: M*/
1157: PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
1158: {
1159:   SNES_NEWTONTR *neP;

1161:   PetscFunctionBegin;
1162:   snes->ops->setup          = SNESSetUp_NEWTONTR;
1163:   snes->ops->solve          = SNESSolve_NEWTONTR;
1164:   snes->ops->reset          = SNESReset_NEWTONTR;
1165:   snes->ops->destroy        = SNESDestroy_NEWTONTR;
1166:   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
1167:   snes->ops->view           = SNESView_NEWTONTR;

1169:   PetscCall(SNESParametersInitialize(snes));
1170:   PetscObjectParameterSetDefault(snes, stol, 0.0);

1172:   snes->usesksp = PETSC_TRUE;
1173:   snes->npcside = PC_RIGHT;
1174:   snes->usesnpc = PETSC_TRUE;

1176:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

1178:   PetscCall(PetscNew(&neP));
1179:   snes->data = (void *)neP;

1181:   PetscObjectParameterSetDefault(neP, eta1, 0.001);
1182:   PetscObjectParameterSetDefault(neP, eta2, 0.25);
1183:   PetscObjectParameterSetDefault(neP, eta3, 0.75);
1184:   PetscObjectParameterSetDefault(neP, t1, 0.25);
1185:   PetscObjectParameterSetDefault(neP, t2, 2.0);
1186:   PetscObjectParameterSetDefault(neP, deltam, PetscDefined(USE_REAL_SINGLE) ? 1.e-6 : 1.e-12);
1187:   PetscObjectParameterSetDefault(neP, delta0, 0.2);
1188:   PetscObjectParameterSetDefault(neP, deltaM, 1.e10);

1190:   neP->norm     = NORM_2;
1191:   neP->fallback = SNES_TR_FALLBACK_NEWTON;
1192:   neP->kmdc     = 0.0; /* by default do not use sufficient decrease */

1194:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRSetTolerances_C", SNESNewtonTRSetTolerances_TR));
1195:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNewtonTRGetTolerances_C", SNESNewtonTRGetTolerances_TR));
1196:   PetscFunctionReturn(PETSC_SUCCESS);
1197: }