xref: /petsc/src/snes/impls/ntrdc/ntrdc.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 
2 #include <../src/snes/impls/ntrdc/ntrdcimpl.h> /*I   "petscsnes.h"   I*/
3 
4 typedef struct {
5   SNES snes;
6   /*  Information on the regular SNES convergence test; which may have been user provided
7       Copied from tr.c (maybe able to disposed, but this is a private function) - Heeho
8       Same with SNESTR_KSPConverged_Private, SNESTR_KSPConverged_Destroy, and SNESTR_Converged_Private
9  */
10 
11   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
12   PetscErrorCode (*convdestroy)(void *);
13   void *convctx;
14 } SNES_TRDC_KSPConverged_Ctx;
15 
16 static PetscErrorCode SNESTRDC_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx) {
17   SNES_TRDC_KSPConverged_Ctx *ctx  = (SNES_TRDC_KSPConverged_Ctx *)cctx;
18   SNES                        snes = ctx->snes;
19   SNES_NEWTONTRDC            *neP  = (SNES_NEWTONTRDC *)snes->data;
20   Vec                         x;
21   PetscReal                   nrm;
22 
23   PetscFunctionBegin;
24   PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
25   if (*reason) { PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm)); }
26   /* Determine norm of solution */
27   PetscCall(KSPBuildSolution(ksp, NULL, &x));
28   PetscCall(VecNorm(x, NORM_2, &nrm));
29   if (nrm >= neP->delta) {
30     PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
31     *reason = KSP_CONVERGED_STEP_LENGTH;
32   }
33   PetscFunctionReturn(0);
34 }
35 
36 static PetscErrorCode SNESTRDC_KSPConverged_Destroy(void *cctx) {
37   SNES_TRDC_KSPConverged_Ctx *ctx = (SNES_TRDC_KSPConverged_Ctx *)cctx;
38 
39   PetscFunctionBegin;
40   PetscCall((*ctx->convdestroy)(ctx->convctx));
41   PetscCall(PetscFree(ctx));
42 
43   PetscFunctionReturn(0);
44 }
45 
46 /* ---------------------------------------------------------------- */
47 /*
48    SNESTRDC_Converged_Private -test convergence JUST for
49    the trust region tolerance.
50 
51 */
52 static PetscErrorCode SNESTRDC_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy) {
53   SNES_NEWTONTRDC *neP = (SNES_NEWTONTRDC *)snes->data;
54 
55   PetscFunctionBegin;
56   *reason = SNES_CONVERGED_ITERATING;
57   if (neP->delta < xnorm * snes->deltatol) {
58     PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g*%g\n", (double)neP->delta, (double)xnorm, (double)snes->deltatol));
59     *reason = SNES_DIVERGED_TR_DELTA;
60   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
61     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
62     *reason = SNES_DIVERGED_FUNCTION_COUNT;
63   }
64   PetscFunctionReturn(0);
65 }
66 
67 /*@
68   SNESNewtonTRDCGetRhoFlag - Get whether the solution update is within the trust-region.
69 
70   Input Parameters:
71 . snes - the nonlinear solver object
72 
73   Output Parameters:
74 . rho_flag: PETSC_TRUE if the solution update is in the trust-region; otherwise, PETSC_FALSE
75 
76   Level: developer
77 
78 @*/
79 PetscErrorCode SNESNewtonTRDCGetRhoFlag(SNES snes, PetscBool *rho_flag) {
80   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
81 
82   PetscFunctionBegin;
83   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
84   PetscValidBoolPointer(rho_flag, 2);
85   *rho_flag = tr->rho_satisfied;
86   PetscFunctionReturn(0);
87 }
88 
89 /*@C
90    SNESNewtonTRDCSetPreCheck - Sets a user function that is called before the search step has been determined.
91        Allows the user a chance to change or override the trust region decision.
92 
93    Logically Collective on snes
94 
95    Input Parameters:
96 +  snes - the nonlinear solver object
97 .  func - [optional] function evaluation routine, see SNESNewtonTRDCPreCheck()  for the calling sequence
98 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
99 
100    Level: intermediate
101 
102    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver.
103 
104 .seealso: `SNESNewtonTRDCPreCheck()`, `SNESNewtonTRDCGetPreCheck()`, `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
105 @*/
106 PetscErrorCode SNESNewtonTRDCSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx) {
107   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
108 
109   PetscFunctionBegin;
110   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
111   if (func) tr->precheck = func;
112   if (ctx) tr->precheckctx = ctx;
113   PetscFunctionReturn(0);
114 }
115 
116 /*@C
117    SNESNewtonTRDCGetPreCheck - Gets the pre-check function
118 
119    Not collective
120 
121    Input Parameter:
122 .  snes - the nonlinear solver context
123 
124    Output Parameters:
125 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPreCheck()
126 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
127 
128    Level: intermediate
129 
130 .seealso: `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCPreCheck()`
131 @*/
132 PetscErrorCode SNESNewtonTRDCGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx) {
133   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
134 
135   PetscFunctionBegin;
136   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
137   if (func) *func = tr->precheck;
138   if (ctx) *ctx = tr->precheckctx;
139   PetscFunctionReturn(0);
140 }
141 
142 /*@C
143    SNESNewtonTRDCSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
144        function evaluation. Allows the user a chance to change or override the decision of the line search routine
145 
146    Logically Collective on snes
147 
148    Input Parameters:
149 +  snes - the nonlinear solver object
150 .  func - [optional] function evaluation routine, see SNESNewtonTRDCPostCheck()  for the calling sequence
151 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
152 
153    Level: intermediate
154 
155    Note: This function is called BEFORE the function evaluation within the SNESNEWTONTRDC solver while the function set in
156    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
157 
158 .seealso: `SNESNewtonTRDCPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
159 @*/
160 PetscErrorCode SNESNewtonTRDCSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx) {
161   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
162 
163   PetscFunctionBegin;
164   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
165   if (func) tr->postcheck = func;
166   if (ctx) tr->postcheckctx = ctx;
167   PetscFunctionReturn(0);
168 }
169 
170 /*@C
171    SNESNewtonTRDCGetPostCheck - Gets the post-check function
172 
173    Not collective
174 
175    Input Parameter:
176 .  snes - the nonlinear solver context
177 
178    Output Parameters:
179 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRDCPostCheck()
180 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
181 
182    Level: intermediate
183 
184 .seealso: `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCPostCheck()`
185 @*/
186 PetscErrorCode SNESNewtonTRDCGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx) {
187   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
188 
189   PetscFunctionBegin;
190   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
191   if (func) *func = tr->postcheck;
192   if (ctx) *ctx = tr->postcheckctx;
193   PetscFunctionReturn(0);
194 }
195 
196 /*@C
197    SNESNewtonTRDCPreCheck - Called before the step has been determined in SNESNEWTONTRDC
198 
199    Logically Collective on snes
200 
201    Input Parameters:
202 +  snes - the solver
203 .  X - The last solution
204 -  Y - The step direction
205 
206    Output Parameters:
207 .  changed_Y - Indicator that the step direction Y has been changed.
208 
209    Level: developer
210 
211 .seealso: `SNESNewtonTRDCSetPreCheck()`, `SNESNewtonTRDCGetPreCheck()`
212 @*/
213 static PetscErrorCode SNESNewtonTRDCPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y) {
214   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
215 
216   PetscFunctionBegin;
217   *changed_Y = PETSC_FALSE;
218   if (tr->precheck) {
219     PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
220     PetscValidLogicalCollectiveBool(snes, *changed_Y, 4);
221   }
222   PetscFunctionReturn(0);
223 }
224 
225 /*@C
226    SNESNewtonTRDCPostCheck - Called after the step has been determined in SNESNEWTONTRDC but before the function evaluation
227 
228    Logically Collective on snes
229 
230    Input Parameters:
231 +  snes - the solver
232 .  X - The last solution
233 .  Y - The full step direction
234 -  W - The updated solution, W = X - Y
235 
236    Output Parameters:
237 +  changed_Y - indicator if step has been changed
238 -  changed_W - Indicator if the new candidate solution W has been changed.
239 
240    Notes:
241      If Y is changed then W is recomputed as X - Y
242 
243    Level: developer
244 
245 .seealso: `SNESNewtonTRDCSetPostCheck()`, `SNESNewtonTRDCGetPostCheck()`
246 @*/
247 static PetscErrorCode SNESNewtonTRDCPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W) {
248   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)snes->data;
249 
250   PetscFunctionBegin;
251   *changed_Y = PETSC_FALSE;
252   *changed_W = PETSC_FALSE;
253   if (tr->postcheck) {
254     PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
255     PetscValidLogicalCollectiveBool(snes, *changed_Y, 5);
256     PetscValidLogicalCollectiveBool(snes, *changed_W, 6);
257   }
258   PetscFunctionReturn(0);
259 }
260 
261 /*
262    SNESSolve_NEWTONTRDC - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
263    (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
264    nonlinear equations
265 
266 */
267 static PetscErrorCode SNESSolve_NEWTONTRDC(SNES snes) {
268   SNES_NEWTONTRDC            *neP = (SNES_NEWTONTRDC *)snes->data;
269   Vec                         X, F, Y, G, W, GradF, YNtmp;
270   Vec                         YCtmp;
271   Mat                         jac;
272   PetscInt                    maxits, i, j, lits, inner_count, bs;
273   PetscReal                   rho, fnorm, gnorm, xnorm = 0, delta, ynorm, temp_xnorm, temp_ynorm; /* TRDC inner iteration */
274   PetscReal                   inorms[99];                                                         /* need to make it dynamic eventually, fixed max block size of 99 for now */
275   PetscReal                   deltaM, ynnorm, f0, mp, gTy, g, yTHy;                               /* rho calculation */
276   PetscReal                   auk, gfnorm, ycnorm, c0, c1, c2, tau, tau_pos, tau_neg, gTBg;       /* Cauchy Point */
277   KSP                         ksp;
278   SNESConvergedReason         reason   = SNES_CONVERGED_ITERATING;
279   PetscBool                   breakout = PETSC_FALSE;
280   SNES_TRDC_KSPConverged_Ctx *ctx;
281   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *);
282   void *convctx;
283 
284   PetscFunctionBegin;
285   maxits = snes->max_its;  /* maximum number of iterations */
286   X      = snes->vec_sol;  /* solution vector */
287   F      = snes->vec_func; /* residual vector */
288   Y      = snes->work[0];  /* update vector */
289   G      = snes->work[1];  /* updated residual */
290   W      = snes->work[2];  /* temporary vector */
291   GradF  = snes->work[3];  /* grad f = J^T F */
292   YNtmp  = snes->work[4];  /* Newton solution */
293   YCtmp  = snes->work[5];  /* Cauchy solution */
294 
295   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);
296 
297   PetscCall(VecGetBlockSize(YNtmp, &bs));
298 
299   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
300   snes->iter = 0;
301   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
302 
303   /* Set the linear stopping criteria to use the More' trick. From tr.c */
304   PetscCall(SNESGetKSP(snes, &ksp));
305   PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
306   if (convtest != SNESTRDC_KSPConverged_Private) {
307     PetscCall(PetscNew(&ctx));
308     ctx->snes = snes;
309     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
310     PetscCall(KSPSetConvergenceTest(ksp, SNESTRDC_KSPConverged_Private, ctx, SNESTRDC_KSPConverged_Destroy));
311     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTRDC_KSPConverged_Private\n"));
312   }
313 
314   if (!snes->vec_func_init_set) {
315     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
316   } else snes->vec_func_init_set = PETSC_FALSE;
317 
318   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
319   SNESCheckFunctionNorm(snes, fnorm);
320   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
321 
322   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
323   snes->norm = fnorm;
324   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
325   delta      = xnorm ? neP->delta0 * xnorm : neP->delta0; /* initial trust region size scaled by xnorm */
326   deltaM     = xnorm ? neP->deltaM * xnorm : neP->deltaM; /* maximum trust region size scaled by xnorm */
327   neP->delta = delta;
328   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
329   PetscCall(SNESMonitor(snes, 0, fnorm));
330 
331   neP->rho_satisfied = PETSC_FALSE;
332 
333   /* test convergence */
334   PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
335   if (snes->reason) PetscFunctionReturn(0);
336 
337   for (i = 0; i < maxits; i++) {
338     PetscBool changed_y;
339     PetscBool changed_w;
340 
341     /* dogleg method */
342     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
343     SNESCheckJacobianDomainerror(snes);
344     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian));
345     PetscCall(KSPSolve(snes->ksp, F, YNtmp)); /* Quasi Newton Solution */
346     SNESCheckKSPSolve(snes);                  /* this is necessary but old tr.c did not have it*/
347     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
348     PetscCall(SNESGetJacobian(snes, &jac, NULL, NULL, NULL));
349 
350     /* rescale Jacobian, Newton solution update, and re-calculate delta for multiphase (multivariable)
351        for inner iteration and Cauchy direction calculation
352     */
353     if (bs > 1 && neP->auto_scale_multiphase) {
354       PetscCall(VecStrideNormAll(YNtmp, NORM_INFINITY, inorms));
355       for (j = 0; j < bs; j++) {
356         if (neP->auto_scale_max > 1.0) {
357           if (inorms[j] < 1.0 / neP->auto_scale_max) { inorms[j] = 1.0 / neP->auto_scale_max; }
358         }
359         PetscCall(VecStrideSet(W, j, inorms[j]));
360         PetscCall(VecStrideScale(YNtmp, j, 1.0 / inorms[j]));
361         PetscCall(VecStrideScale(X, j, 1.0 / inorms[j]));
362       }
363       PetscCall(VecNorm(X, NORM_2, &xnorm));
364       if (i == 0) {
365         delta = neP->delta0 * xnorm;
366       } else {
367         delta = neP->delta * xnorm;
368       }
369       deltaM = neP->deltaM * xnorm;
370       PetscCall(MatDiagonalScale(jac, PETSC_NULL, W));
371     }
372 
373     /* calculating GradF of minimization function */
374     PetscCall(MatMultTranspose(jac, F, GradF)); /* grad f = J^T F */
375     PetscCall(VecNorm(YNtmp, NORM_2, &ynnorm)); /* ynnorm <- || Y_newton || */
376 
377     inner_count        = 0;
378     neP->rho_satisfied = PETSC_FALSE;
379     while (1) {
380       if (ynnorm <= delta) { /* see if the Newton solution is within the trust region */
381         PetscCall(VecCopy(YNtmp, Y));
382       } else if (neP->use_cauchy) { /* use Cauchy direction if enabled */
383         PetscCall(MatMult(jac, GradF, W));
384         PetscCall(VecDotRealPart(W, W, &gTBg));     /* completes GradF^T J^T J GradF */
385         PetscCall(VecNorm(GradF, NORM_2, &gfnorm)); /* grad f norm <- || grad f || */
386         if (gTBg <= 0.0) {
387           auk = PETSC_MAX_REAL;
388         } else {
389           auk = PetscSqr(gfnorm) / gTBg;
390         }
391         auk = PetscMin(delta / gfnorm, auk);
392         PetscCall(VecCopy(GradF, YCtmp));           /* this could be improved */
393         PetscCall(VecScale(YCtmp, auk));            /* YCtmp, Cauchy solution*/
394         PetscCall(VecNorm(YCtmp, NORM_2, &ycnorm)); /* ycnorm <- || Y_cauchy || */
395         if (ycnorm >= delta) {                      /* see if the Cauchy solution meets the criteria */
396           PetscCall(VecCopy(YCtmp, Y));
397           PetscCall(PetscInfo(snes, "DL evaluated. delta: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)delta, (double)ynnorm, (double)ycnorm));
398         } else {                                  /* take ratio, tau, of Cauchy and Newton direction and step */
399           PetscCall(VecAXPY(YNtmp, -1.0, YCtmp)); /* YCtmp = A, YNtmp = B */
400           PetscCall(VecNorm(YNtmp, NORM_2, &c0)); /* this could be improved */
401           c0 = PetscSqr(c0);
402           PetscCall(VecDotRealPart(YCtmp, YNtmp, &c1));
403           c1 = 2.0 * c1;
404           PetscCall(VecNorm(YCtmp, NORM_2, &c2)); /* this could be improved */
405           c2      = PetscSqr(c2) - PetscSqr(delta);
406           tau_pos = (c1 + PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0); /* quadratic formula */
407           tau_neg = (c1 - PetscSqrtReal(PetscSqr(c1) - 4. * c0 * c2)) / (2. * c0);
408           tau     = PetscMax(tau_pos, tau_neg); /* can tau_neg > tau_pos? I don't think so, but just in case. */
409           PetscCall(PetscInfo(snes, "DL evaluated. tau: %8.4e, ynnorm: %8.4e, ycnorm: %8.4e\n", (double)tau, (double)ynnorm, (double)ycnorm));
410           PetscCall(VecWAXPY(W, tau, YNtmp, YCtmp));
411           PetscCall(VecAXPY(W, -tau, YCtmp));
412           PetscCall(VecCopy(W, Y)); /* this could be improved */
413         }
414       } else {
415         /* if Cauchy is disabled, only use Newton direction */
416         auk = delta / ynnorm;
417         PetscCall(VecScale(YNtmp, auk));
418         PetscCall(VecCopy(YNtmp, Y)); /* this could be improved (many VecCopy, VecNorm)*/
419       }
420 
421       PetscCall(VecNorm(Y, NORM_2, &ynorm)); /* compute the final ynorm  */
422       f0 = 0.5 * PetscSqr(fnorm);            /* minimizing function f(X) */
423       PetscCall(MatMult(jac, Y, W));
424       PetscCall(VecDotRealPart(W, W, &yTHy)); /* completes GradY^T J^T J GradY */
425       PetscCall(VecDotRealPart(GradF, Y, &gTy));
426       mp = f0 - gTy + 0.5 * yTHy; /* quadratic model to satisfy, -gTy because our update is X-Y*/
427 
428       /* scale back solution update */
429       if (bs > 1 && neP->auto_scale_multiphase) {
430         for (j = 0; j < bs; j++) {
431           PetscCall(VecStrideScale(Y, j, inorms[j]));
432           if (inner_count == 0) {
433             /* TRDC inner algorithm does not need scaled X after calculating delta in the outer iteration */
434             /* need to scale back X to match Y and provide proper update to the external code */
435             PetscCall(VecStrideScale(X, j, inorms[j]));
436           }
437         }
438         if (inner_count == 0) PetscCall(VecNorm(X, NORM_2, &temp_xnorm)); /* only in the first iteration */
439         PetscCall(VecNorm(Y, NORM_2, &temp_ynorm));
440       } else {
441         temp_xnorm = xnorm;
442         temp_ynorm = ynorm;
443       }
444       inner_count++;
445 
446       /* Evaluate the solution to meet the improvement ratio criteria */
447       PetscCall(SNESNewtonTRDCPreCheck(snes, X, Y, &changed_y));
448       PetscCall(VecWAXPY(W, -1.0, Y, X));
449       PetscCall(SNESNewtonTRDCPostCheck(snes, X, Y, W, &changed_y, &changed_w));
450       if (changed_y) PetscCall(VecWAXPY(W, -1.0, Y, X));
451       PetscCall(VecCopy(Y, snes->vec_sol_update));
452       PetscCall(SNESComputeFunction(snes, W, G)); /*  F(X-Y) = G */
453       PetscCall(VecNorm(G, NORM_2, &gnorm));      /* gnorm <- || g || */
454       SNESCheckFunctionNorm(snes, gnorm);
455       g = 0.5 * PetscSqr(gnorm); /* minimizing function g(W) */
456       if (f0 == mp) rho = 0.0;
457       else rho = (f0 - g) / (f0 - mp); /* actual improvement over predicted improvement */
458 
459       if (rho < neP->eta2) {
460         delta *= neP->t1; /* shrink the region */
461       } else if (rho > neP->eta3) {
462         delta = PetscMin(neP->t2 * delta, deltaM); /* expand the region, but not greater than deltaM */
463       }
464 
465       neP->delta = delta;
466       if (rho >= neP->eta1) {
467         /* unscale delta and xnorm before going to the next outer iteration */
468         if (bs > 1 && neP->auto_scale_multiphase) {
469           neP->delta = delta / xnorm;
470           xnorm      = temp_xnorm;
471           ynorm      = temp_ynorm;
472         }
473         neP->rho_satisfied = PETSC_TRUE;
474         break; /* the improvement ratio is satisfactory */
475       }
476       PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
477 
478       /* check to see if progress is hopeless */
479       neP->itflag = PETSC_FALSE;
480       /* both delta, ynorm, and xnorm are either scaled or unscaled */
481       PetscCall(SNESTRDC_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP));
482       if (!reason) {
483         /* temp_xnorm, temp_ynorm is always unscaled */
484         /* also the inner iteration already calculated the Jacobian and solved the matrix */
485         /* therefore, it should be passing iteration number of iter+1 instead of iter+0 in the first iteration and after */
486         PetscCall((*snes->ops->converged)(snes, snes->iter + 1, temp_xnorm, temp_ynorm, fnorm, &reason, snes->cnvP));
487       }
488       /* if multiphase state changes, break out inner iteration */
489       if (reason == SNES_BREAKOUT_INNER_ITER) {
490         if (bs > 1 && neP->auto_scale_multiphase) {
491           /* unscale delta and xnorm before going to the next outer iteration */
492           neP->delta = delta / xnorm;
493           xnorm      = temp_xnorm;
494           ynorm      = temp_ynorm;
495         }
496         reason = SNES_CONVERGED_ITERATING;
497         break;
498       }
499       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
500       if (reason) {
501         if (reason < 0) {
502           /* We're not progressing, so return with the current iterate */
503           PetscCall(SNESMonitor(snes, i + 1, fnorm));
504           breakout = PETSC_TRUE;
505           break;
506         } else if (reason > 0) {
507           /* We're converged, so return with the current iterate and update solution */
508           PetscCall(SNESMonitor(snes, i + 1, fnorm));
509           breakout = PETSC_FALSE;
510           break;
511         }
512       }
513       snes->numFailures++;
514     }
515     if (!breakout) {
516       /* Update function and solution vectors */
517       fnorm = gnorm;
518       PetscCall(VecCopy(G, F));
519       PetscCall(VecCopy(W, X));
520       /* Monitor convergence */
521       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
522       snes->iter  = i + 1;
523       snes->norm  = fnorm;
524       snes->xnorm = xnorm;
525       snes->ynorm = ynorm;
526       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
527       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
528       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
529       /* Test for convergence, xnorm = || X || */
530       neP->itflag = PETSC_TRUE;
531       if (snes->ops->converged != SNESConvergedSkip) PetscCall(VecNorm(X, NORM_2, &xnorm));
532       PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &reason, snes->cnvP);
533       if (reason) break;
534     } else break;
535   }
536 
537   /* PetscCall(PetscFree(inorms)); */
538   if (i == maxits) {
539     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
540     if (!reason) reason = SNES_DIVERGED_MAX_IT;
541   }
542   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
543   snes->reason = reason;
544   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
545   if (convtest != SNESTRDC_KSPConverged_Private) {
546     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
547     PetscCall(PetscFree(ctx));
548     PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
549   }
550   PetscFunctionReturn(0);
551 }
552 
553 /*------------------------------------------------------------*/
554 static PetscErrorCode SNESSetUp_NEWTONTRDC(SNES snes) {
555   PetscFunctionBegin;
556   PetscCall(SNESSetWorkVecs(snes, 6));
557   PetscCall(SNESSetUpMatrices(snes));
558   PetscFunctionReturn(0);
559 }
560 
561 PetscErrorCode SNESReset_NEWTONTRDC(SNES snes) {
562   PetscFunctionBegin;
563   PetscFunctionReturn(0);
564 }
565 
566 static PetscErrorCode SNESDestroy_NEWTONTRDC(SNES snes) {
567   PetscFunctionBegin;
568   PetscCall(SNESReset_NEWTONTRDC(snes));
569   PetscCall(PetscFree(snes->data));
570   PetscFunctionReturn(0);
571 }
572 /*------------------------------------------------------------*/
573 
574 static PetscErrorCode SNESSetFromOptions_NEWTONTRDC(SNES snes, PetscOptionItems *PetscOptionsObject) {
575   SNES_NEWTONTRDC *ctx = (SNES_NEWTONTRDC *)snes->data;
576 
577   PetscFunctionBegin;
578   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
579   PetscCall(PetscOptionsReal("-snes_trdc_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL));
580   PetscCall(PetscOptionsReal("-snes_trdc_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL));
581   PetscCall(PetscOptionsReal("-snes_trdc_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL));
582   PetscCall(PetscOptionsReal("-snes_trdc_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL));
583   PetscCall(PetscOptionsReal("-snes_trdc_t1", "t1", "None", ctx->t1, &ctx->t1, NULL));
584   PetscCall(PetscOptionsReal("-snes_trdc_t2", "t2", "None", ctx->t2, &ctx->t2, NULL));
585   PetscCall(PetscOptionsReal("-snes_trdc_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL));
586   PetscCall(PetscOptionsReal("-snes_trdc_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
587   PetscCall(PetscOptionsReal("-snes_trdc_auto_scale_max", "auto_scale_max", "None", ctx->auto_scale_max, &ctx->auto_scale_max, NULL));
588   PetscCall(PetscOptionsBool("-snes_trdc_use_cauchy", "use_cauchy", "use Cauchy step and direction", ctx->use_cauchy, &ctx->use_cauchy, NULL));
589   PetscCall(PetscOptionsBool("-snes_trdc_auto_scale_multiphase", "auto_scale_multiphase", "Auto scaling for proper cauchy direction", ctx->auto_scale_multiphase, &ctx->auto_scale_multiphase, NULL));
590   PetscOptionsHeadEnd();
591   PetscFunctionReturn(0);
592 }
593 
594 static PetscErrorCode SNESView_NEWTONTRDC(SNES snes, PetscViewer viewer) {
595   SNES_NEWTONTRDC *tr = (SNES_NEWTONTRDC *)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 (-snes_trtol)\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   }
605   PetscFunctionReturn(0);
606 }
607 /* ------------------------------------------------------------ */
608 /*MC
609       SNESNEWTONTRDC - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
610 
611    Options Database:
612 +   -snes_trdc_tol <tol> - trust region tolerance
613 .   -snes_trdc_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
614 .   -snes_trdc_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
615 .   -snes_trdc_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
616 .   -snes_trdc_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
617 .   -snes_trdc_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
618 .   -snes_trdc_deltaM <deltaM> - trust region parameter, max size of trust region, deltaM*norm2(x) (default: 0.5)
619 .   -snes_trdc_delta0 <delta0> - trust region parameter, initial size of trust region, delta0*norm2(x) (default: 0.1)
620 .   -snes_trdc_auto_scale_max <auto_scale_max> - used with auto_scale_multiphase, caps the maximum auto-scaling factor
621 .   -snes_trdc_use_cauchy <use_cauchy> - True uses dogleg Cauchy (Steepest Descent direction) step & direction in the trust region algorithm
622 -   -snes_trdc_auto_scale_multiphase <auto_scale_multiphase> - True turns on auto-scaling for multivariable block matrix for Cauchy and trust region
623 
624     Notes:
625     The algorithm is taken from "Linear and Nonlinear Solvers for Simulating Multiphase Flow
626     within Large-Scale Engineered Subsurface Systems" by Heeho D. Park, Glenn E. Hammond,
627     Albert J. Valocchi, Tara LaForce.
628 
629    Level: intermediate
630 
631 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`, `SNESNEWTONTRDC`
632 
633 M*/
634 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTRDC(SNES snes) {
635   SNES_NEWTONTRDC *neP;
636 
637   PetscFunctionBegin;
638   snes->ops->setup          = SNESSetUp_NEWTONTRDC;
639   snes->ops->solve          = SNESSolve_NEWTONTRDC;
640   snes->ops->destroy        = SNESDestroy_NEWTONTRDC;
641   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTRDC;
642   snes->ops->view           = SNESView_NEWTONTRDC;
643   snes->ops->reset          = SNESReset_NEWTONTRDC;
644 
645   snes->usesksp = PETSC_TRUE;
646   snes->usesnpc = PETSC_FALSE;
647 
648   snes->alwayscomputesfinalresidual = PETSC_TRUE;
649 
650   PetscCall(PetscNewLog(snes, &neP));
651   snes->data                 = (void *)neP;
652   neP->delta                 = 0.0;
653   neP->delta0                = 0.1;
654   neP->eta1                  = 0.001;
655   neP->eta2                  = 0.25;
656   neP->eta3                  = 0.75;
657   neP->t1                    = 0.25;
658   neP->t2                    = 2.0;
659   neP->deltaM                = 0.5;
660   neP->sigma                 = 0.0001;
661   neP->itflag                = PETSC_FALSE;
662   neP->rnorm0                = 0.0;
663   neP->ttol                  = 0.0;
664   neP->use_cauchy            = PETSC_TRUE;
665   neP->auto_scale_multiphase = PETSC_FALSE;
666   neP->auto_scale_max        = -1.0;
667   neP->rho_satisfied         = PETSC_FALSE;
668   snes->deltatol             = 1.e-12;
669 
670   /* for multiphase (multivariable) scaling */
671   /* may be used for dynamic allocation of inorms, but it fails snes_tutorials-ex3_13
672      on test forced DIVERGED_JACOBIAN_DOMAIN test. I will use static array for now.
673   PetscCall(VecGetBlockSize(snes->work[0],&neP->bs));
674   PetscCall(PetscCalloc1(neP->bs,&neP->inorms));
675   */
676 
677   PetscFunctionReturn(0);
678 }
679