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