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