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