xref: /petsc/src/snes/impls/tr/tr.c (revision 6d5076fa512d0d8e4e12a3685978afdbfee2b1a4)
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 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: `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: `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: `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: `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: `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 Parameters:
224 .  changed_Y - Indicator that the step direction Y has been changed.
225 
226    Level: intermediate
227 
228 .seealso: `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: `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, xnorm = 0, delta, ynorm;
321   PetscReal                 deltaM, fk, fkp1, deltaqm, gTy, yTHy;
322   PetscReal                 auk, gfnorm, ycnorm, gTBg;
323   KSP                       ksp;
324   PetscBool                 already_done = PETSC_FALSE;
325   PetscBool                 clear_converged_test, rho_satisfied, has_objective;
326   PetscVoidFunction         ksp_has_radius;
327   SNES_TR_KSPConverged_Ctx *ctx;
328   void                     *convctx;
329   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *);
330   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
331 
332   PetscFunctionBegin;
333   PetscCall(SNESGetObjective(snes, &objective, NULL));
334   has_objective = objective ? PETSC_TRUE : PETSC_FALSE;
335 
336   maxits = snes->max_its;                                   /* maximum number of iterations */
337   X      = snes->vec_sol;                                   /* solution vector */
338   F      = snes->vec_func;                                  /* residual vector */
339   Y      = snes->vec_sol_update;                            /* update vector */
340   G      = snes->work[0];                                   /* updated residual */
341   W      = snes->work[1];                                   /* temporary vector */
342   GradF  = !has_objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */
343   YU     = snes->work[3];                                   /* work vector for dogleg method */
344 
345   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);
346 
347   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
348   snes->iter = 0;
349   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
350 
351   /* Set the linear stopping criteria to use the More' trick if needed */
352   clear_converged_test = PETSC_FALSE;
353   PetscCall(SNESGetKSP(snes, &ksp));
354   PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
355   PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &ksp_has_radius));
356   if (convtest != SNESTR_KSPConverged_Private && !ksp_has_radius) {
357     clear_converged_test = PETSC_TRUE;
358     PetscCall(PetscNew(&ctx));
359     ctx->snes = snes;
360     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
361     PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy));
362     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n"));
363   }
364 
365   if (!snes->vec_func_init_set) {
366     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
367   } else snes->vec_func_init_set = PETSC_FALSE;
368 
369   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
370   SNESCheckFunctionNorm(snes, fnorm);
371   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
372 
373   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
374   snes->norm = fnorm;
375   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
376   delta      = neP->delta0;
377   deltaM     = neP->deltaM;
378   neP->delta = delta;
379   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
380   PetscCall(SNESMonitor(snes, 0, fnorm));
381 
382   /* test convergence */
383   rho_satisfied = PETSC_FALSE;
384   PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
385   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
386 
387   if (has_objective) PetscCall(SNESComputeObjective(snes, X, &fk));
388   else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */
389 
390   while (snes->iter < maxits) {
391     PetscBool changed_y;
392     PetscBool changed_w;
393 
394     /* calculating GradF of minimization function only once */
395     if (!already_done) {
396       PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
397       SNESCheckJacobianDomainerror(snes);
398       if (has_objective) gfnorm = fnorm;
399       else {
400         PetscCall(MatMultTranspose(snes->jacobian, F, GradF)); /* grad f = J^T F */
401         PetscCall(VecNorm(GradF, NORM_2, &gfnorm));
402       }
403     }
404     already_done = PETSC_TRUE;
405 
406     /* solve trust-region subproblem */
407     PetscCall(KSPCGSetRadius(snes->ksp, delta));
408     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
409     PetscCall(KSPSolve(snes->ksp, F, Y));
410     SNESCheckKSPSolve(snes);
411     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
412 
413     /* decide what to do when the update is outside of trust region */
414     PetscCall(VecNorm(Y, NORM_2, &ynorm));
415     if (ynorm > delta) {
416       switch (neP->fallback) {
417       case SNES_TR_FALLBACK_NEWTON:
418         auk = delta / ynorm;
419         PetscCall(VecScale(Y, auk));
420         break;
421       case SNES_TR_FALLBACK_CAUCHY:
422       case SNES_TR_FALLBACK_DOGLEG:
423         PetscCall(MatMult(snes->jacobian, GradF, W));
424         if (has_objective) PetscCall(VecDotRealPart(GradF, W, &gTBg));
425         else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */
426         /* Eqs 4.7 and 4.8 in Nocedal and Wright */
427         auk = delta / gfnorm;
428         if (gTBg > 0.0) auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1);
429         ycnorm = auk * gfnorm;
430         if (neP->fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) {
431           /* Cauchy solution */
432           PetscCall(VecAXPBY(Y, auk, 0.0, GradF));
433           PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg));
434         } else { /* take linear combination of Cauchy and Newton direction and step */
435           PetscReal c0, c1, c2, tau = 0.0, tpos, tneg;
436           PetscBool noroots;
437 
438           auk = gfnorm * gfnorm / gTBg;
439           PetscCall(VecAXPBY(YU, auk, 0.0, GradF));
440           PetscCall(VecAXPY(Y, -1.0, YU));
441           PetscCall(VecNorm(Y, NORM_2, &c0));
442           PetscCall(VecDotRealPart(YU, Y, &c1));
443           c0 = PetscSqr(c0);
444           c2 = PetscSqr(ycnorm) - PetscSqr(delta);
445           PetscQuadraticRoots(c0, c1, c2, &tneg, &tpos);
446 
447           noroots = PetscIsInfOrNanReal(tneg);
448           if (noroots) { /*  No roots, select Cauchy point */
449             auk = delta / gfnorm;
450             auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1);
451             PetscCall(VecAXPBY(Y, auk, 0.0, GradF));
452           } else { /* Here roots corresponds to tau-1 in Nocedal and Wright */
453             tpos += 1.0;
454             tneg += 1.0;
455             tau = PetscClipInterval(tpos, 0.0, 2.0); /* clip to tau [0,2] */
456             if (tau < 1.0) {
457               PetscCall(VecAXPBY(Y, tau, 0.0, YU));
458             } else {
459               PetscCall(VecAXPBY(Y, 1.0, tau - 1, YU));
460             }
461           }
462           PetscCall(VecNorm(Y, NORM_2, &c0)); /* this norm will be cached and reused later */
463           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));
464         }
465         break;
466       default:
467         SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode");
468         break;
469       }
470     }
471 
472     /* Evaluate the solution to meet the improvement ratio criteria */
473 
474     /* compute the final ynorm */
475     PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y));
476     PetscCall(VecNorm(Y, NORM_2, &ynorm));
477 
478     /* compute the quadratic model difference */
479     PetscCall(SNESNewtonTRQuadraticDelta(snes, has_objective, Y, GradF, W, &yTHy, &gTy, &deltaqm));
480 
481     /* update */
482     PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */
483     PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w));
484     if (changed_y) {
485       /* Need to recompute the quadratic model difference */
486       PetscCall(SNESNewtonTRQuadraticDelta(snes, has_objective, Y, GradF, YU, &yTHy, &gTy, &deltaqm));
487       /* User changed Y but not W */
488       if (!changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X));
489     }
490 
491     /* Compute new objective function */
492     PetscCall(SNESComputeFunction(snes, W, G)); /*  F(Xkp1) = G */
493     PetscCall(VecNorm(G, NORM_2, &gnorm));
494     if (has_objective) PetscCall(SNESComputeObjective(snes, W, &fkp1));
495     else fkp1 = 0.5 * PetscSqr(gnorm);
496     SNESCheckFunctionNorm(snes, fkp1);
497 
498     if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */
499     else rho = -1.0;                                /*  no reduction in quadratic model, step must be rejected */
500     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));
501 
502     if (rho < neP->eta2) delta *= neP->t1;      /* shrink the region */
503     else if (rho > neP->eta3) delta *= neP->t2; /* expand the region */
504     delta = PetscMin(delta, deltaM);            /* but not greater than deltaM */
505 
506     neP->delta = delta;
507     if (rho >= neP->eta1) {
508       rho_satisfied = PETSC_TRUE;
509     } else {
510       rho_satisfied = PETSC_FALSE;
511       PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
512       /* check to see if progress is hopeless */
513       PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP));
514       if (!snes->reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
515       if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_INNER;
516       snes->numFailures++;
517       /* We're not progressing, so return with the current iterate */
518       if (snes->reason) break;
519     }
520     if (rho_satisfied) {
521       /* Update function values */
522       already_done = PETSC_FALSE;
523       fnorm        = gnorm;
524       fk           = fkp1;
525 
526       /* New residual and linearization point */
527       PetscCall(VecCopy(G, F));
528       PetscCall(VecCopy(W, X));
529 
530       /* Monitor convergence */
531       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
532       snes->iter++;
533       snes->norm  = fnorm;
534       snes->xnorm = xnorm;
535       snes->ynorm = ynorm;
536       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
537       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
538       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
539 
540       /* Test for convergence, xnorm = || X || */
541       PetscCall(VecNorm(X, NORM_2, &xnorm));
542       PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
543       if (snes->reason) break;
544     }
545   }
546 
547   if (snes->iter == maxits) {
548     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
549     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
550   }
551   if (clear_converged_test) {
552     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
553     PetscCall(PetscFree(ctx));
554     PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
555   }
556   PetscFunctionReturn(PETSC_SUCCESS);
557 }
558 
559 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
560 {
561   PetscFunctionBegin;
562   PetscCall(SNESSetWorkVecs(snes, 4));
563   PetscCall(SNESSetUpMatrices(snes));
564   PetscFunctionReturn(PETSC_SUCCESS);
565 }
566 
567 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
568 {
569   PetscFunctionBegin;
570   PetscCall(PetscFree(snes->data));
571   PetscFunctionReturn(PETSC_SUCCESS);
572 }
573 
574 static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject)
575 {
576   SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data;
577 
578   PetscFunctionBegin;
579   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
580   PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL));
581   PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL));
582   PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL));
583   PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL));
584   PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL));
585   PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL));
586   PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL));
587   PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
588   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));
589   PetscOptionsHeadEnd();
590   PetscFunctionReturn(PETSC_SUCCESS);
591 }
592 
593 static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer)
594 {
595   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
596   PetscBool      iascii;
597 
598   PetscFunctionBegin;
599   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
600   if (iascii) {
601     PetscCall(PetscViewerASCIIPrintf(viewer, "  Trust region tolerance %g\n", (double)snes->deltatol));
602     PetscCall(PetscViewerASCIIPrintf(viewer, "  eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
603     PetscCall(PetscViewerASCIIPrintf(viewer, "  delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM));
604     PetscCall(PetscViewerASCIIPrintf(viewer, "  fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback]));
605   }
606   PetscFunctionReturn(PETSC_SUCCESS);
607 }
608 
609 /*MC
610       SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
611 
612    Options Database Keys:
613 +   -snes_tr_tol <tol> - trust region tolerance
614 .   -snes_tr_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
615 .   -snes_tr_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
616 .   -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
617 .   -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
618 .   -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
619 .   -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL)
620 .   -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2)
621 -   -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.
622 
623     Reference:
624 .   * - "Numerical Optimization" by Nocedal and Wright, chapter 4.
625 
626    Level: deprecated (since 3.19)
627 
628 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`,
629           `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`,
630           `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()`
631 M*/
632 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
633 {
634   SNES_NEWTONTR *neP;
635 
636   PetscFunctionBegin;
637   snes->ops->setup          = SNESSetUp_NEWTONTR;
638   snes->ops->solve          = SNESSolve_NEWTONTR;
639   snes->ops->destroy        = SNESDestroy_NEWTONTR;
640   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
641   snes->ops->view           = SNESView_NEWTONTR;
642 
643   snes->usesksp = PETSC_TRUE;
644   snes->usesnpc = PETSC_FALSE;
645 
646   snes->alwayscomputesfinalresidual = PETSC_TRUE;
647 
648   PetscCall(PetscNew(&neP));
649   snes->data    = (void *)neP;
650   neP->delta    = 0.0;
651   neP->delta0   = 0.2;
652   neP->eta1     = 0.001;
653   neP->eta2     = 0.25;
654   neP->eta3     = 0.75;
655   neP->t1       = 0.25;
656   neP->t2       = 2.0;
657   neP->deltaM   = 1.e10;
658   neP->fallback = SNES_TR_FALLBACK_NEWTON;
659   PetscFunctionReturn(PETSC_SUCCESS);
660 }
661