xref: /petsc/src/snes/impls/tr/tr.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 #include <../src/snes/impls/tr/trimpl.h> /*I   "petscsnes.h"   I*/
2 
3 typedef struct {
4   SNES snes;
5   /*  Information on the regular SNES convergence test; which may have been user provided */
6   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
7   PetscErrorCode (*convdestroy)(void *);
8   void *convctx;
9 } SNES_TR_KSPConverged_Ctx;
10 
11 static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) {
12   SNES_TR_KSPConverged_Ctx *ctx  = (SNES_TR_KSPConverged_Ctx *)cctx;
13   SNES                      snes = ctx->snes;
14   SNES_NEWTONTR            *neP  = (SNES_NEWTONTR *)snes->data;
15   Vec                       x;
16   PetscReal                 nrm;
17 
18   PetscFunctionBegin;
19   PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
20   if (*reason) { PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); }
21   /* Determine norm of solution */
22   PetscCall(KSPBuildSolution(ksp, NULL, &x));
23   PetscCall(VecNorm(x, NORM_2, &nrm));
24   if (nrm >= neP->delta) {
25     PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
26     *reason = KSP_CONVERGED_STEP_LENGTH;
27   }
28   PetscFunctionReturn(0);
29 }
30 
31 static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx) {
32   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx;
33 
34   PetscFunctionBegin;
35   PetscCall((*ctx->convdestroy)(ctx->convctx));
36   PetscCall(PetscFree(ctx));
37   PetscFunctionReturn(0);
38 }
39 
40 /* ---------------------------------------------------------------- */
41 /*
42    SNESTR_Converged_Private -test convergence JUST for
43    the trust region tolerance.
44 
45 */
46 static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) {
47   SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data;
48 
49   PetscFunctionBegin;
50   *reason = SNES_CONVERGED_ITERATING;
51   if (neP->delta < xnorm * snes->deltatol) {
52     PetscCall(PetscInfo(snes, "Converged due to trust region param %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol));
53     *reason = SNES_DIVERGED_TR_DELTA;
54   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
55     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
56     *reason = SNES_DIVERGED_FUNCTION_COUNT;
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 /*@C
62    SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined.
63        Allows the user a chance to change or override the decision of the line search routine.
64 
65    Logically Collective on snes
66 
67    Input Parameters:
68 +  snes - the nonlinear solver object
69 .  func - [optional] function evaluation routine, see SNESNewtonTRPreCheck()  for the calling sequence
70 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
71 
72    Level: intermediate
73 
74    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver.
75 
76 .seealso: `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`
77 @*/
78 PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) {
79   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
80 
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
83   if (func) tr->precheck = func;
84   if (ctx) tr->precheckctx = ctx;
85   PetscFunctionReturn(0);
86 }
87 
88 /*@C
89    SNESNewtonTRGetPreCheck - Gets the pre-check function
90 
91    Not collective
92 
93    Input Parameter:
94 .  snes - the nonlinear solver context
95 
96    Output Parameters:
97 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPreCheck()
98 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
99 
100    Level: intermediate
101 
102 .seealso: `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()`
103 @*/
104 PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) {
105   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
106 
107   PetscFunctionBegin;
108   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
109   if (func) *func = tr->precheck;
110   if (ctx) *ctx = tr->precheckctx;
111   PetscFunctionReturn(0);
112 }
113 
114 /*@C
115    SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
116        function evaluation. Allows the user a chance to change or override the decision of the line search routine
117 
118    Logically Collective on snes
119 
120    Input Parameters:
121 +  snes - the nonlinear solver object
122 .  func - [optional] function evaluation routine, see SNESNewtonTRPostCheck()  for the calling sequence
123 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
124 
125    Level: intermediate
126 
127    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTR solver while the function set in
128    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
129 
130 .seealso: `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`
131 @*/
132 PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) {
133   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
134 
135   PetscFunctionBegin;
136   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
137   if (func) tr->postcheck = func;
138   if (ctx) tr->postcheckctx = ctx;
139   PetscFunctionReturn(0);
140 }
141 
142 /*@C
143    SNESNewtonTRGetPostCheck - Gets the post-check function
144 
145    Not collective
146 
147    Input Parameter:
148 .  snes - the nonlinear solver context
149 
150    Output Parameters:
151 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck()
152 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
153 
154    Level: intermediate
155 
156 .seealso: `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()`
157 @*/
158 PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) {
159   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
160 
161   PetscFunctionBegin;
162   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
163   if (func) *func = tr->postcheck;
164   if (ctx) *ctx = tr->postcheckctx;
165   PetscFunctionReturn(0);
166 }
167 
168 /*@C
169    SNESNewtonTRPreCheck - Called before the step has been determined in SNESNEWTONTR
170 
171    Logically Collective on snes
172 
173    Input Parameters:
174 +  snes - the solver
175 .  X - The last solution
176 -  Y - The step direction
177 
178    Output Parameters:
179 .  changed_Y - Indicator that the step direction Y has been changed.
180 
181    Level: developer
182 
183 .seealso: `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`
184 @*/
185 static PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) {
186   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
187 
188   PetscFunctionBegin;
189   *changed_Y = PETSC_FALSE;
190   if (tr->precheck) {
191     PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
192     PetscValidLogicalCollectiveBool(snes, *changed_Y, 4);
193   }
194   PetscFunctionReturn(0);
195 }
196 
197 /*@C
198    SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation
199 
200    Logically Collective on snes
201 
202    Input Parameters:
203 +  snes - the solver
204 .  X - The last solution
205 .  Y - The full step direction
206 -  W - The updated solution, W = X - Y
207 
208    Output Parameters:
209 +  changed_Y - indicator if step has been changed
210 -  changed_W - Indicator if the new candidate solution W has been changed.
211 
212    Notes:
213      If Y is changed then W is recomputed as X - Y
214 
215    Level: developer
216 
217 .seealso: `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`
218 @*/
219 static PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) {
220   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
221 
222   PetscFunctionBegin;
223   *changed_Y = PETSC_FALSE;
224   *changed_W = PETSC_FALSE;
225   if (tr->postcheck) {
226     PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
227     PetscValidLogicalCollectiveBool(snes, *changed_Y, 5);
228     PetscValidLogicalCollectiveBool(snes, *changed_W, 6);
229   }
230   PetscFunctionReturn(0);
231 }
232 
233 /*
234    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
235    region approach for solving systems of nonlinear equations.
236 
237 */
238 static PetscErrorCode SNESSolve_NEWTONTR(SNES snes) {
239   SNES_NEWTONTR            *neP = (SNES_NEWTONTR *)snes->data;
240   Vec                       X, F, Y, G, Ytmp, W;
241   PetscInt                  maxits, i, lits;
242   PetscReal                 rho, fnorm, gnorm, gpnorm, xnorm = 0, delta, nrm, ynorm, norm1;
243   PetscScalar               cnorm;
244   KSP                       ksp;
245   SNESConvergedReason       reason   = SNES_CONVERGED_ITERATING;
246   PetscBool                 breakout = PETSC_FALSE;
247   SNES_TR_KSPConverged_Ctx *ctx;
248   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *);
249   void *convctx;
250 
251   PetscFunctionBegin;
252   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);
253 
254   maxits = snes->max_its;  /* maximum number of iterations */
255   X      = snes->vec_sol;  /* solution vector */
256   F      = snes->vec_func; /* residual vector */
257   Y      = snes->work[0];  /* work vectors */
258   G      = snes->work[1];
259   Ytmp   = snes->work[2];
260   W      = snes->work[3];
261 
262   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
263   snes->iter = 0;
264   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
265 
266   /* Set the linear stopping criteria to use the More' trick. */
267   PetscCall(SNESGetKSP(snes, &ksp));
268   PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
269   if (convtest != SNESTR_KSPConverged_Private) {
270     PetscCall(PetscNew(&ctx));
271     ctx->snes = snes;
272     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
273     PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy));
274     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n"));
275   }
276 
277   if (!snes->vec_func_init_set) {
278     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
279   } else snes->vec_func_init_set = PETSC_FALSE;
280 
281   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
282   SNESCheckFunctionNorm(snes, fnorm);
283   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
284   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
285   snes->norm = fnorm;
286   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
287   delta      = xnorm ? neP->delta0 * xnorm : neP->delta0;
288   neP->delta = delta;
289   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
290   PetscCall(SNESMonitor(snes, 0, fnorm));
291 
292   /* test convergence */
293   PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
294   if (snes->reason) PetscFunctionReturn(0);
295 
296   for (i = 0; i < maxits; i++) {
297     /* Call general purpose update function */
298     PetscTryTypeMethod(snes, update, snes->iter);
299 
300     /* Solve J Y = F, where J is Jacobian matrix */
301     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
302     SNESCheckJacobianDomainerror(snes);
303     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
304     PetscCall(KSPSolve(snes->ksp, F, Ytmp));
305     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
306     snes->linear_its += lits;
307 
308     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
309     PetscCall(VecNorm(Ytmp, NORM_2, &nrm));
310     norm1 = nrm;
311 
312     while (1) {
313       PetscBool changed_y;
314       PetscBool changed_w;
315       PetscCall(VecCopy(Ytmp, Y));
316       nrm = norm1;
317 
318       /* Scale Y if need be and predict new value of F norm */
319       if (nrm >= delta) {
320         nrm    = delta / nrm;
321         gpnorm = (1.0 - nrm) * fnorm;
322         cnorm  = nrm;
323         PetscCall(PetscInfo(snes, "Scaling direction by %g\n", (double)nrm));
324         PetscCall(VecScale(Y, cnorm));
325         nrm   = gpnorm;
326         ynorm = delta;
327       } else {
328         gpnorm = 0.0;
329         PetscCall(PetscInfo(snes, "Direction is in Trust Region\n"));
330         ynorm = nrm;
331       }
332       /* PreCheck() allows for updates to Y prior to W <- X - Y */
333       PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y));
334       PetscCall(VecWAXPY(W, -1.0, Y, X)); /* W <- X - Y */
335       PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w));
336       if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X));
337       PetscCall(VecCopy(Y, snes->vec_sol_update));
338       PetscCall(SNESComputeFunction(snes, W, G)); /*  F(X-Y) = G */
339       PetscCall(VecNorm(G, NORM_2, &gnorm));      /* gnorm <- || g || */
340       SNESCheckFunctionNorm(snes, gnorm);
341       if (fnorm == gpnorm) rho = 0.0;
342       else rho = (fnorm * fnorm - gnorm * gnorm) / (fnorm * fnorm - gpnorm * gpnorm);
343 
344       /* Update size of trust region */
345       if (rho < neP->mu) delta *= neP->delta1;
346       else if (rho < neP->eta) delta *= neP->delta2;
347       else delta *= neP->delta3;
348       PetscCall(PetscInfo(snes, "fnorm=%g, gnorm=%g, ynorm=%g\n", (double)fnorm, (double)gnorm, (double)ynorm));
349       PetscCall(PetscInfo(snes, "gpred=%g, rho=%g, delta=%g\n", (double)gpnorm, (double)rho, (double)delta));
350 
351       neP->delta = delta;
352       if (rho > neP->sigma) break;
353       PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
354 
355       /* check to see if progress is hopeless */
356       neP->itflag = PETSC_FALSE;
357       PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP));
358       if (!reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP);
359       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
360       if (reason) {
361         /* We're not progressing, so return with the current iterate */
362         PetscCall(SNESMonitor(snes, i + 1, fnorm));
363         breakout = PETSC_TRUE;
364         break;
365       }
366       snes->numFailures++;
367     }
368     if (!breakout) {
369       /* Update function and solution vectors */
370       fnorm = gnorm;
371       PetscCall(VecCopy(G, F));
372       PetscCall(VecCopy(W, X));
373       /* Monitor convergence */
374       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
375       snes->iter  = i + 1;
376       snes->norm  = fnorm;
377       snes->xnorm = xnorm;
378       snes->ynorm = ynorm;
379       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
380       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
381       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
382       /* Test for convergence, xnorm = || X || */
383       neP->itflag = PETSC_TRUE;
384       if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
385       PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP);
386       if (reason) break;
387     } else break;
388   }
389 
390   if (i == maxits) {
391     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
392     if (!reason) reason = SNES_DIVERGED_MAX_IT;
393   }
394   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
395   snes->reason = reason;
396   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
397   if (convtest != SNESTR_KSPConverged_Private) {
398     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
399     PetscCall(PetscFree(ctx));
400     PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
401   }
402   PetscFunctionReturn(0);
403 }
404 
405 /*------------------------------------------------------------*/
406 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes) {
407   PetscFunctionBegin;
408   PetscCall(SNESSetWorkVecs(snes, 4));
409   PetscCall(SNESSetUpMatrices(snes));
410   PetscFunctionReturn(0);
411 }
412 
413 PetscErrorCode SNESReset_NEWTONTR(SNES snes) {
414   PetscFunctionBegin;
415   PetscFunctionReturn(0);
416 }
417 
418 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes) {
419   PetscFunctionBegin;
420   PetscCall(SNESReset_NEWTONTR(snes));
421   PetscCall(PetscFree(snes->data));
422   PetscFunctionReturn(0);
423 }
424 /*------------------------------------------------------------*/
425 
426 static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject) {
427   SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data;
428 
429   PetscFunctionBegin;
430   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
431   PetscCall(PetscOptionsReal("-snes_trtol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL));
432   PetscCall(PetscOptionsReal("-snes_tr_mu", "mu", "None", ctx->mu, &ctx->mu, NULL));
433   PetscCall(PetscOptionsReal("-snes_tr_eta", "eta", "None", ctx->eta, &ctx->eta, NULL));
434   PetscCall(PetscOptionsReal("-snes_tr_sigma", "sigma", "None", ctx->sigma, &ctx->sigma, NULL));
435   PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
436   PetscCall(PetscOptionsReal("-snes_tr_delta1", "delta1", "None", ctx->delta1, &ctx->delta1, NULL));
437   PetscCall(PetscOptionsReal("-snes_tr_delta2", "delta2", "None", ctx->delta2, &ctx->delta2, NULL));
438   PetscCall(PetscOptionsReal("-snes_tr_delta3", "delta3", "None", ctx->delta3, &ctx->delta3, NULL));
439   PetscOptionsHeadEnd();
440   PetscFunctionReturn(0);
441 }
442 
443 static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer) {
444   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
445   PetscBool      iascii;
446 
447   PetscFunctionBegin;
448   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
449   if (iascii) {
450     PetscCall(PetscViewerASCIIPrintf(viewer, "  Trust region tolerance %g (-snes_trtol)\n", (double)snes->deltatol));
451     PetscCall(PetscViewerASCIIPrintf(viewer, "  mu=%g, eta=%g, sigma=%g\n", (double)tr->mu, (double)tr->eta, (double)tr->sigma));
452     PetscCall(PetscViewerASCIIPrintf(viewer, "  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n", (double)tr->delta0, (double)tr->delta1, (double)tr->delta2, (double)tr->delta3));
453   }
454   PetscFunctionReturn(0);
455 }
456 /* ------------------------------------------------------------ */
457 /*MC
458       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
459 
460    Options Database:
461 +    -snes_trtol <tol> - trust region tolerance
462 .    -snes_tr_mu <mu> - trust region parameter
463 .    -snes_tr_eta <eta> - trust region parameter
464 .    -snes_tr_sigma <sigma> - trust region parameter
465 .    -snes_tr_delta0 <delta0> -  initial size of the trust region is delta0*norm2(x)
466 .    -snes_tr_delta1 <delta1> - trust region parameter
467 .    -snes_tr_delta2 <delta2> - trust region parameter
468 -    -snes_tr_delta3 <delta3> - trust region parameter
469 
470    The basic algorithm is taken from "The Minpack Project", by More',
471    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
472    of Mathematical Software", Wayne Cowell, editor.
473 
474    Level: intermediate
475 
476 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`
477 
478 M*/
479 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes) {
480   SNES_NEWTONTR *neP;
481 
482   PetscFunctionBegin;
483   snes->ops->setup          = SNESSetUp_NEWTONTR;
484   snes->ops->solve          = SNESSolve_NEWTONTR;
485   snes->ops->destroy        = SNESDestroy_NEWTONTR;
486   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
487   snes->ops->view           = SNESView_NEWTONTR;
488   snes->ops->reset          = SNESReset_NEWTONTR;
489 
490   snes->usesksp = PETSC_TRUE;
491   snes->usesnpc = PETSC_FALSE;
492 
493   snes->alwayscomputesfinalresidual = PETSC_TRUE;
494 
495   PetscCall(PetscNewLog(snes, &neP));
496   snes->data  = (void *)neP;
497   neP->mu     = 0.25;
498   neP->eta    = 0.75;
499   neP->delta  = 0.0;
500   neP->delta0 = 0.2;
501   neP->delta1 = 0.3;
502   neP->delta2 = 0.75;
503   neP->delta3 = 2.0;
504   neP->sigma  = 0.0001;
505   neP->itflag = PETSC_FALSE;
506   neP->rnorm0 = 0.0;
507   neP->ttol   = 0.0;
508   PetscFunctionReturn(0);
509 }
510