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