xref: /petsc/src/snes/impls/tr/tr.c (revision 4a221d5978a6d86d172efe6a78940aca257e68d3)
1c6db04a5SJed Brown #include <../src/snes/impls/tr/trimpl.h> /*I   "petscsnes.h"   I*/
24800dd8cSBarry Smith 
3971273eeSBarry Smith typedef struct {
4971273eeSBarry Smith   SNES snes;
5df8705c3SBarry Smith   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *);
6df8705c3SBarry Smith   PetscErrorCode (*convdestroy)(void *);
7df8705c3SBarry Smith   void *convctx;
8971273eeSBarry Smith } SNES_TR_KSPConverged_Ctx;
9971273eeSBarry Smith 
10*4a221d59SStefano Zampini const char *const SNESNewtonTRFallbackTypes[] = {"NEWTON", "CAUCHY", "DOGLEG", "SNESNewtonTRFallbackType", "SNES_TR_FALLBACK_", NULL};
11*4a221d59SStefano Zampini 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp, PetscInt n, PetscReal rnorm, KSPConvergedReason *reason, void *cctx)
13d71ae5a4SJacob Faibussowitsch {
14971273eeSBarry Smith   SNES_TR_KSPConverged_Ctx *ctx  = (SNES_TR_KSPConverged_Ctx *)cctx;
15971273eeSBarry Smith   SNES                      snes = ctx->snes;
1604d7464bSBarry Smith   SNES_NEWTONTR            *neP  = (SNES_NEWTONTR *)snes->data;
17df60cc22SBarry Smith   Vec                       x;
18064f8208SBarry Smith   PetscReal                 nrm;
19df60cc22SBarry Smith 
203a40ed3dSBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   PetscCall((*ctx->convtest)(ksp, n, rnorm, reason, ctx->convctx));
2248a46eb9SPierre Jolivet   if (*reason) PetscCall(PetscInfo(snes, "Default or user provided convergence test KSP iterations=%" PetscInt_FMT ", rnorm=%g\n", n, (double)rnorm));
23a935fc98SLois Curfman McInnes   /* Determine norm of solution */
249566063dSJacob Faibussowitsch   PetscCall(KSPBuildSolution(ksp, NULL, &x));
259566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &nrm));
26064f8208SBarry Smith   if (nrm >= neP->delta) {
279566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Ending linear iteration early, delta=%g, length=%g\n", (double)neP->delta, (double)nrm));
28329f5518SBarry Smith     *reason = KSP_CONVERGED_STEP_LENGTH;
29df60cc22SBarry Smith   }
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31df60cc22SBarry Smith }
3282bf6240SBarry Smith 
33d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx)
34d71ae5a4SJacob Faibussowitsch {
35971273eeSBarry Smith   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx *)cctx;
36971273eeSBarry Smith 
37971273eeSBarry Smith   PetscFunctionBegin;
389566063dSJacob Faibussowitsch   PetscCall((*ctx->convdestroy)(ctx->convctx));
399566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41971273eeSBarry Smith }
42971273eeSBarry Smith 
43d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESTR_Converged_Private(SNES snes, PetscInt it, PetscReal xnorm, PetscReal pnorm, PetscReal fnorm, SNESConvergedReason *reason, void *dummy)
44d71ae5a4SJacob Faibussowitsch {
4504d7464bSBarry Smith   SNES_NEWTONTR *neP = (SNES_NEWTONTR *)snes->data;
4685385478SLisandro Dalcin 
4785385478SLisandro Dalcin   PetscFunctionBegin;
4885385478SLisandro Dalcin   *reason = SNES_CONVERGED_ITERATING;
49*4a221d59SStefano Zampini   if (neP->delta < snes->deltatol) {
50*4a221d59SStefano Zampini     PetscCall(PetscInfo(snes, "Diverged due to too small a trust region %g<%g\n", (double)neP->delta, (double)snes->deltatol));
511c6b2ff8SBarry Smith     *reason = SNES_DIVERGED_TR_DELTA;
52e71169deSBarry Smith   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
5363a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Exceeded maximum number of function evaluations: %" PetscInt_FMT "\n", snes->max_funcs));
5485385478SLisandro Dalcin     *reason = SNES_DIVERGED_FUNCTION_COUNT;
5585385478SLisandro Dalcin   }
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5785385478SLisandro Dalcin }
5885385478SLisandro Dalcin 
59*4a221d59SStefano Zampini /*@
60*4a221d59SStefano Zampini   SNESNewtonTRSetFallbackType - Set the type of fallback if the solution of the trust region subproblem is outside the radius
61*4a221d59SStefano Zampini 
62*4a221d59SStefano Zampini   Input Parameters:
63*4a221d59SStefano Zampini + snes - the nonlinear solver object
64*4a221d59SStefano Zampini - ftype - the fallback type, see `SNESNewtonTRFallbackType`
65*4a221d59SStefano Zampini 
66*4a221d59SStefano Zampini   Level: intermediate
67*4a221d59SStefano Zampini 
68*4a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPreCheck()`,
69*4a221d59SStefano Zampini           `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`
70*4a221d59SStefano Zampini @*/
71*4a221d59SStefano Zampini PetscErrorCode SNESNewtonTRSetFallbackType(SNES snes, SNESNewtonTRFallbackType ftype)
72*4a221d59SStefano Zampini {
73*4a221d59SStefano Zampini   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
74*4a221d59SStefano Zampini   PetscBool      flg;
75*4a221d59SStefano Zampini 
76*4a221d59SStefano Zampini   PetscFunctionBegin;
77*4a221d59SStefano Zampini   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
78*4a221d59SStefano Zampini   PetscValidLogicalCollectiveEnum(snes, ftype, 2);
79*4a221d59SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
80*4a221d59SStefano Zampini   if (flg) tr->fallback = ftype;
81*4a221d59SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
82*4a221d59SStefano Zampini }
83*4a221d59SStefano Zampini 
847cb011f5SBarry Smith /*@C
85c9368356SGlenn Hammond    SNESNewtonTRSetPreCheck - Sets a user function that is called before the search step has been determined.
86*4a221d59SStefano Zampini        Allows the user a chance to change or override the trust region decision.
87f6dfbefdSBarry Smith 
88c3339decSBarry Smith    Logically Collective
89c9368356SGlenn Hammond 
90c9368356SGlenn Hammond    Input Parameters:
91c9368356SGlenn Hammond +  snes - the nonlinear solver object
92f6dfbefdSBarry Smith .  func - [optional] function evaluation routine, see `SNESNewtonTRPreCheck()`  for the calling sequence
93c9368356SGlenn Hammond -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
94c9368356SGlenn Hammond 
95c9368356SGlenn Hammond    Level: intermediate
96c9368356SGlenn Hammond 
97f6dfbefdSBarry Smith    Note:
98*4a221d59SStefano Zampini    This function is called BEFORE the function evaluation within the solver.
99c9368356SGlenn Hammond 
100*4a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`,
101c9368356SGlenn Hammond @*/
102d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPreCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, PetscBool *, void *), void *ctx)
103d71ae5a4SJacob Faibussowitsch {
104c9368356SGlenn Hammond   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
105*4a221d59SStefano Zampini   PetscBool      flg;
106c9368356SGlenn Hammond 
107c9368356SGlenn Hammond   PetscFunctionBegin;
108c9368356SGlenn Hammond   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
109*4a221d59SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
110*4a221d59SStefano Zampini   if (flg) {
111c9368356SGlenn Hammond     if (func) tr->precheck = func;
112c9368356SGlenn Hammond     if (ctx) tr->precheckctx = ctx;
113*4a221d59SStefano Zampini   }
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115c9368356SGlenn Hammond }
116c9368356SGlenn Hammond 
117c9368356SGlenn Hammond /*@C
118c9368356SGlenn Hammond    SNESNewtonTRGetPreCheck - Gets the pre-check function
119c9368356SGlenn Hammond 
120c9368356SGlenn Hammond    Not collective
121c9368356SGlenn Hammond 
122c9368356SGlenn Hammond    Input Parameter:
123c9368356SGlenn Hammond .  snes - the nonlinear solver context
124c9368356SGlenn Hammond 
125c9368356SGlenn Hammond    Output Parameters:
126f6dfbefdSBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPreCheck()`
127c9368356SGlenn Hammond -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
128c9368356SGlenn Hammond 
129c9368356SGlenn Hammond    Level: intermediate
130c9368356SGlenn Hammond 
131*4a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRPreCheck()`
132c9368356SGlenn Hammond @*/
133d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPreCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, PetscBool *, void *), void **ctx)
134d71ae5a4SJacob Faibussowitsch {
135c9368356SGlenn Hammond   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
136*4a221d59SStefano Zampini   PetscBool      flg;
137c9368356SGlenn Hammond 
138c9368356SGlenn Hammond   PetscFunctionBegin;
139c9368356SGlenn Hammond   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
140*4a221d59SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
141*4a221d59SStefano Zampini   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
142c9368356SGlenn Hammond   if (func) *func = tr->precheck;
143c9368356SGlenn Hammond   if (ctx) *ctx = tr->precheckctx;
1443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
145c9368356SGlenn Hammond }
146c9368356SGlenn Hammond 
147c9368356SGlenn Hammond /*@C
1487cb011f5SBarry Smith    SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
149*4a221d59SStefano Zampini        function evaluation. Allows the user a chance to change or override the internal decision of the solver
150f6dfbefdSBarry Smith 
151c3339decSBarry Smith    Logically Collective
1527cb011f5SBarry Smith 
1537cb011f5SBarry Smith    Input Parameters:
1547cb011f5SBarry Smith +  snes - the nonlinear solver object
155f6dfbefdSBarry Smith .  func - [optional] function evaluation routine, see `SNESNewtonTRPostCheck()`  for the calling sequence
1567cb011f5SBarry Smith -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
1577cb011f5SBarry Smith 
1587cb011f5SBarry Smith    Level: intermediate
1597cb011f5SBarry Smith 
160f6dfbefdSBarry Smith    Note:
161*4a221d59SStefano Zampini    This function is called BEFORE the function evaluation within the solver while the function set in
162f6dfbefdSBarry Smith    `SNESLineSearchSetPostCheck()` is called AFTER the function evaluation.
1637cb011f5SBarry Smith 
164*4a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`
1657cb011f5SBarry Smith @*/
166d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void *ctx)
167d71ae5a4SJacob Faibussowitsch {
1687cb011f5SBarry Smith   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
169*4a221d59SStefano Zampini   PetscBool      flg;
1707cb011f5SBarry Smith 
1717cb011f5SBarry Smith   PetscFunctionBegin;
1727cb011f5SBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
173*4a221d59SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
174*4a221d59SStefano Zampini   if (flg) {
1757cb011f5SBarry Smith     if (func) tr->postcheck = func;
1767cb011f5SBarry Smith     if (ctx) tr->postcheckctx = ctx;
177*4a221d59SStefano Zampini   }
1783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1797cb011f5SBarry Smith }
1807cb011f5SBarry Smith 
1817cb011f5SBarry Smith /*@C
1827cb011f5SBarry Smith    SNESNewtonTRGetPostCheck - Gets the post-check function
1837cb011f5SBarry Smith 
1847cb011f5SBarry Smith    Not collective
1857cb011f5SBarry Smith 
1867cb011f5SBarry Smith    Input Parameter:
1877cb011f5SBarry Smith .  snes - the nonlinear solver context
1887cb011f5SBarry Smith 
1897cb011f5SBarry Smith    Output Parameters:
190f6dfbefdSBarry Smith +  func - [optional] function evaluation routine, see for the calling sequence `SNESNewtonTRPostCheck()`
1917cb011f5SBarry Smith -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
1927cb011f5SBarry Smith 
1937cb011f5SBarry Smith    Level: intermediate
1947cb011f5SBarry Smith 
195*4a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRPostCheck()`
1967cb011f5SBarry Smith @*/
197d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES, Vec, Vec, Vec, PetscBool *, PetscBool *, void *), void **ctx)
198d71ae5a4SJacob Faibussowitsch {
1997cb011f5SBarry Smith   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
200*4a221d59SStefano Zampini   PetscBool      flg;
2017cb011f5SBarry Smith 
2027cb011f5SBarry Smith   PetscFunctionBegin;
2037cb011f5SBarry Smith   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
204*4a221d59SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
205*4a221d59SStefano Zampini   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
2067cb011f5SBarry Smith   if (func) *func = tr->postcheck;
2077cb011f5SBarry Smith   if (ctx) *ctx = tr->postcheckctx;
2083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2097cb011f5SBarry Smith }
2107cb011f5SBarry Smith 
2117cb011f5SBarry Smith /*@C
212*4a221d59SStefano Zampini    SNESNewtonTRPreCheck - Runs the precheck routine
213c9368356SGlenn Hammond 
214c3339decSBarry Smith    Logically Collective
215c9368356SGlenn Hammond 
216c9368356SGlenn Hammond    Input Parameters:
217c9368356SGlenn Hammond +  snes - the solver
218c9368356SGlenn Hammond .  X - The last solution
219c9368356SGlenn Hammond -  Y - The step direction
220c9368356SGlenn Hammond 
221c9368356SGlenn Hammond    Output Parameters:
222c9368356SGlenn Hammond .  changed_Y - Indicator that the step direction Y has been changed.
223c9368356SGlenn Hammond 
224*4a221d59SStefano Zampini    Level: intermediate
225c9368356SGlenn Hammond 
226*4a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRPostCheck()`
227c9368356SGlenn Hammond @*/
228*4a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPreCheck(SNES snes, Vec X, Vec Y, PetscBool *changed_Y)
229d71ae5a4SJacob Faibussowitsch {
230c9368356SGlenn Hammond   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
231*4a221d59SStefano Zampini   PetscBool      flg;
232c9368356SGlenn Hammond 
233c9368356SGlenn Hammond   PetscFunctionBegin;
234*4a221d59SStefano Zampini   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
235*4a221d59SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
236*4a221d59SStefano Zampini   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
237c9368356SGlenn Hammond   *changed_Y = PETSC_FALSE;
238c9368356SGlenn Hammond   if (tr->precheck) {
2399566063dSJacob Faibussowitsch     PetscCall((*tr->precheck)(snes, X, Y, changed_Y, tr->precheckctx));
240c9368356SGlenn Hammond     PetscValidLogicalCollectiveBool(snes, *changed_Y, 4);
241c9368356SGlenn Hammond   }
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
243c9368356SGlenn Hammond }
244c9368356SGlenn Hammond 
245c9368356SGlenn Hammond /*@C
246*4a221d59SStefano Zampini    SNESNewtonTRPostCheck - Runs the postcheck routine
2477cb011f5SBarry Smith 
248c3339decSBarry Smith    Logically Collective
2497cb011f5SBarry Smith 
2507cb011f5SBarry Smith    Input Parameters:
2516b867d5aSJose E. Roman +  snes - the solver
2526b867d5aSJose E. Roman .  X - The last solution
2537cb011f5SBarry Smith .  Y - The full step direction
2543312a946SBarry Smith -  W - The updated solution, W = X - Y
2557cb011f5SBarry Smith 
2567cb011f5SBarry Smith    Output Parameters:
2573312a946SBarry Smith +  changed_Y - indicator if step has been changed
2583312a946SBarry Smith -  changed_W - Indicator if the new candidate solution W has been changed.
2597cb011f5SBarry Smith 
260f6dfbefdSBarry Smith    Note:
2613312a946SBarry Smith      If Y is changed then W is recomputed as X - Y
2627cb011f5SBarry Smith 
263*4a221d59SStefano Zampini    Level: intermediate
2647cb011f5SBarry Smith 
265*4a221d59SStefano Zampini .seealso: `SNESNEWTONTR`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`, `SNESNewtonTRPreCheck()
2667cb011f5SBarry Smith @*/
267*4a221d59SStefano Zampini PetscErrorCode SNESNewtonTRPostCheck(SNES snes, Vec X, Vec Y, Vec W, PetscBool *changed_Y, PetscBool *changed_W)
268d71ae5a4SJacob Faibussowitsch {
2697cb011f5SBarry Smith   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
270*4a221d59SStefano Zampini   PetscBool      flg;
2717cb011f5SBarry Smith 
2727cb011f5SBarry Smith   PetscFunctionBegin;
273*4a221d59SStefano Zampini   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
274*4a221d59SStefano Zampini   PetscCall(PetscObjectTypeCompare((PetscObject)snes, SNESNEWTONTR, &flg));
275*4a221d59SStefano Zampini   PetscAssert(flg, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONG, "Not for type %s", ((PetscObject)snes)->type_name);
276c9368356SGlenn Hammond   *changed_Y = PETSC_FALSE;
2777cb011f5SBarry Smith   *changed_W = PETSC_FALSE;
2787cb011f5SBarry Smith   if (tr->postcheck) {
2799566063dSJacob Faibussowitsch     PetscCall((*tr->postcheck)(snes, X, Y, W, changed_Y, changed_W, tr->postcheckctx));
280c9368356SGlenn Hammond     PetscValidLogicalCollectiveBool(snes, *changed_Y, 5);
2817cb011f5SBarry Smith     PetscValidLogicalCollectiveBool(snes, *changed_W, 6);
2827cb011f5SBarry Smith   }
2833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2847cb011f5SBarry Smith }
28585385478SLisandro Dalcin 
286*4a221d59SStefano Zampini static inline void PetscQuadraticRoots(PetscReal a, PetscReal b, PetscReal c, PetscReal *xm, PetscReal *xp)
287*4a221d59SStefano Zampini {
288*4a221d59SStefano Zampini   PetscReal temp = -0.5 * (b + PetscCopysignReal(1.0, b) * PetscSqrtReal(b * b - 4 * a * c));
289*4a221d59SStefano Zampini   PetscReal x1   = temp / a;
290*4a221d59SStefano Zampini   PetscReal x2   = c / temp;
291*4a221d59SStefano Zampini   *xm            = PetscMin(x1, x2);
292*4a221d59SStefano Zampini   *xp            = PetscMax(x1, x2);
293*4a221d59SStefano Zampini }
294*4a221d59SStefano Zampini 
295df60cc22SBarry Smith /*
296*4a221d59SStefano Zampini    SNESSolve_NEWTONTR - Implements Newton's Method with trust-region subproblem and adds dogleg Cauchy
297*4a221d59SStefano Zampini    (Steepest Descent direction) step and direction if the trust region is not satisfied for solving system of
298*4a221d59SStefano Zampini    nonlinear equations
2994800dd8cSBarry Smith 
3004800dd8cSBarry Smith */
301d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
302d71ae5a4SJacob Faibussowitsch {
30304d7464bSBarry Smith   SNES_NEWTONTR            *neP = (SNES_NEWTONTR *)snes->data;
304*4a221d59SStefano Zampini   Vec                       X, F, Y, G, W, GradF, YU;
305*4a221d59SStefano Zampini   PetscInt                  maxits, lits;
306*4a221d59SStefano Zampini   PetscReal                 rho, fnorm, gnorm, xnorm = 0, delta, ynorm;
307*4a221d59SStefano Zampini   PetscReal                 deltaM, fk, fkp1, deltaqm, gTy, yTHy;
308*4a221d59SStefano Zampini   PetscReal                 auk, gfnorm, ycnorm, gTBg;
309df60cc22SBarry Smith   KSP                       ksp;
310*4a221d59SStefano Zampini   PetscBool                 already_done = PETSC_FALSE;
311*4a221d59SStefano Zampini   PetscBool                 clear_converged_test, rho_satisfied;
312*4a221d59SStefano Zampini   PetscVoidFunction         ksp_has_radius;
313df8705c3SBarry Smith   SNES_TR_KSPConverged_Ctx *ctx;
3145e28bcb6Sprj-   void                     *convctx;
315*4a221d59SStefano Zampini   PetscErrorCode (*convtest)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), (*convdestroy)(void *);
316*4a221d59SStefano Zampini   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
3174800dd8cSBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
319*4a221d59SStefano Zampini   PetscCall(SNESGetObjective(snes, &objective, NULL));
320c579b300SPatrick Farrell 
321fbe28522SBarry Smith   maxits = snes->max_its;                               /* maximum number of iterations */
322fbe28522SBarry Smith   X      = snes->vec_sol;                               /* solution vector */
32339e2f89bSBarry Smith   F      = snes->vec_func;                              /* residual vector */
324*4a221d59SStefano Zampini   Y      = snes->vec_sol_update;                        /* update vector */
325*4a221d59SStefano Zampini   G      = snes->work[0];                               /* updated residual */
326*4a221d59SStefano Zampini   W      = snes->work[1];                               /* temporary vector */
327*4a221d59SStefano Zampini   GradF  = !objective ? snes->work[2] : snes->vec_func; /* grad f = J^T F */
328*4a221d59SStefano Zampini   YU     = snes->work[3];                               /* work vector for dogleg method */
329*4a221d59SStefano Zampini 
330*4a221d59SStefano Zampini   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);
3314800dd8cSBarry Smith 
3329566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
3334c49b128SBarry Smith   snes->iter = 0;
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3354800dd8cSBarry Smith 
336*4a221d59SStefano Zampini   /* Set the linear stopping criteria to use the More' trick if needed */
337*4a221d59SStefano Zampini   clear_converged_test = PETSC_FALSE;
3389566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
3399566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergenceTest(ksp, &convtest, &convctx, &convdestroy));
340*4a221d59SStefano Zampini   PetscCall(PetscObjectQueryFunction((PetscObject)ksp, "KSPCGSetRadius_C", &ksp_has_radius));
341*4a221d59SStefano Zampini   if (convtest != SNESTR_KSPConverged_Private && !ksp_has_radius) {
342*4a221d59SStefano Zampini     clear_converged_test = PETSC_TRUE;
3439566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
344df8705c3SBarry Smith     ctx->snes = snes;
3459566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
3469566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, SNESTR_KSPConverged_Private, ctx, SNESTR_KSPConverged_Destroy));
3479566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Using Krylov convergence test SNESTR_KSPConverged_Private\n"));
348df8705c3SBarry Smith   }
349df8705c3SBarry Smith 
350e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
3519566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F)); /* F(X) */
3521aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
353e4ed7901SPeter Brune 
3549566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- || F || */
355422a814eSBarry Smith   SNESCheckFunctionNorm(snes, fnorm);
3569566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &xnorm)); /* xnorm <- || X || */
357*4a221d59SStefano Zampini 
3589566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
359fbe28522SBarry Smith   snes->norm = fnorm;
3609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
361*4a221d59SStefano Zampini   delta      = neP->delta0;
362*4a221d59SStefano Zampini   deltaM     = neP->deltaM;
3634800dd8cSBarry Smith   neP->delta = delta;
3649566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
3659566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
366b37302c6SLois Curfman McInnes 
36785385478SLisandro Dalcin   /* test convergence */
368*4a221d59SStefano Zampini   rho_satisfied = PETSC_FALSE;
369dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
3703ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
3713f149594SLisandro Dalcin 
372*4a221d59SStefano Zampini   if (objective) PetscCall(SNESComputeObjective(snes, X, &fk));
373*4a221d59SStefano Zampini   else fk = 0.5 * PetscSqr(fnorm); /* obj(x) = 0.5 * ||F(x)||^2 */
37499a96b7cSMatthew Knepley 
375*4a221d59SStefano Zampini   while (snes->iter < maxits) {
376c9368356SGlenn Hammond     PetscBool changed_y;
3777cb011f5SBarry Smith     PetscBool changed_w;
3786b5873e3SBarry Smith 
379*4a221d59SStefano Zampini     /* solve trust-region subproblem */
380*4a221d59SStefano Zampini     if (!already_done) {
381*4a221d59SStefano Zampini       PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
382*4a221d59SStefano Zampini       SNESCheckJacobianDomainerror(snes);
383fbe28522SBarry Smith     }
384*4a221d59SStefano Zampini     PetscCall(KSPCGSetRadius(snes->ksp, delta));
385*4a221d59SStefano Zampini     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
386*4a221d59SStefano Zampini     PetscCall(KSPSolve(snes->ksp, F, Y));
387*4a221d59SStefano Zampini     SNESCheckKSPSolve(snes);
388*4a221d59SStefano Zampini     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
3894800dd8cSBarry Smith 
390*4a221d59SStefano Zampini     /* calculating GradF of minimization function only once */
391*4a221d59SStefano Zampini     if (!already_done) {
392*4a221d59SStefano Zampini       if (objective) gfnorm = fnorm;
393*4a221d59SStefano Zampini       else {
394*4a221d59SStefano Zampini         PetscCall(MatMultTranspose(snes->jacobian, F, GradF)); /* grad f = J^T F */
395*4a221d59SStefano Zampini         PetscCall(VecNorm(GradF, NORM_2, &gfnorm));
396*4a221d59SStefano Zampini       }
397*4a221d59SStefano Zampini       already_done = PETSC_TRUE;
398*4a221d59SStefano Zampini     }
3991aa26658SKarl Rupp 
400*4a221d59SStefano Zampini     /* decide what to do when the update is outside of trust region */
401*4a221d59SStefano Zampini     PetscCall(VecNorm(Y, NORM_2, &ynorm));
402*4a221d59SStefano Zampini     if (ynorm > delta) {
403*4a221d59SStefano Zampini       switch (neP->fallback) {
404*4a221d59SStefano Zampini       case SNES_TR_FALLBACK_NEWTON:
405*4a221d59SStefano Zampini         auk = delta / ynorm;
406*4a221d59SStefano Zampini         PetscCall(VecScale(Y, auk));
407*4a221d59SStefano Zampini         break;
408*4a221d59SStefano Zampini       case SNES_TR_FALLBACK_CAUCHY:
409*4a221d59SStefano Zampini       case SNES_TR_FALLBACK_DOGLEG:
410*4a221d59SStefano Zampini         PetscCall(MatMult(snes->jacobian, GradF, W));
411*4a221d59SStefano Zampini         if (objective) PetscCall(VecDotRealPart(GradF, W, &gTBg));
412*4a221d59SStefano Zampini         else PetscCall(VecDotRealPart(W, W, &gTBg)); /* B = J^t * J */
413*4a221d59SStefano Zampini         /* Eqs 4.7 and 4.8 in Nocedal and Wright */
414*4a221d59SStefano Zampini         auk = delta / gfnorm;
415*4a221d59SStefano Zampini         if (gTBg > 0.0) auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1);
416*4a221d59SStefano Zampini         ycnorm = auk * gfnorm;
417*4a221d59SStefano Zampini         if (neP->fallback == SNES_TR_FALLBACK_CAUCHY || gTBg <= 0.0) {
418*4a221d59SStefano Zampini           /* Cauchy solution */
419*4a221d59SStefano Zampini           PetscCall(VecAXPBY(Y, auk, 0.0, GradF));
420*4a221d59SStefano Zampini           PetscCall(PetscInfo(snes, "CP evaluated. delta: %g, ynorm: %g, ycnorm: %g, gTBg: %g\n", (double)delta, (double)ynorm, (double)ycnorm, (double)gTBg));
421*4a221d59SStefano Zampini         } else { /* take linear combination of Cauchy and Newton direction and step */
422*4a221d59SStefano Zampini           PetscReal c0, c1, c2, tau = 0.0, tpos, tneg;
423*4a221d59SStefano Zampini           PetscBool noroots;
424284fb49fSHeeho Park 
425*4a221d59SStefano Zampini           auk = gfnorm * gfnorm / gTBg;
426*4a221d59SStefano Zampini           PetscCall(VecAXPBY(YU, auk, 0.0, GradF));
427*4a221d59SStefano Zampini           PetscCall(VecAXPY(Y, -1.0, YU));
428*4a221d59SStefano Zampini           PetscCall(VecNorm(Y, NORM_2, &c0));
429*4a221d59SStefano Zampini           PetscCall(VecDotRealPart(YU, Y, &c1));
430*4a221d59SStefano Zampini           c0 = PetscSqr(c0);
431*4a221d59SStefano Zampini           c2 = PetscSqr(ycnorm) - PetscSqr(delta);
432*4a221d59SStefano Zampini           PetscQuadraticRoots(c0, c1, c2, &tneg, &tpos);
433*4a221d59SStefano Zampini 
434*4a221d59SStefano Zampini           noroots = PetscIsInfOrNanReal(tneg);
435*4a221d59SStefano Zampini           if (noroots) { /*  No roots, select Cauchy point */
436*4a221d59SStefano Zampini             auk = delta / gfnorm;
437*4a221d59SStefano Zampini             auk *= PetscMin(gfnorm * gfnorm * gfnorm / (delta * gTBg), 1);
438*4a221d59SStefano Zampini             PetscCall(VecAXPBY(Y, auk, 0.0, GradF));
439*4a221d59SStefano Zampini           } else { /* Here roots corresponds to tau-1 in Nocedal and Wright */
440*4a221d59SStefano Zampini             tpos += 1.0;
441*4a221d59SStefano Zampini             tneg += 1.0;
442*4a221d59SStefano Zampini             tau = PetscClipInterval(tpos, 0.0, 2.0); /* clip to tau [0,2] */
443*4a221d59SStefano Zampini             if (tau < 1.0) {
444*4a221d59SStefano Zampini               PetscCall(VecAXPBY(Y, tau, 0.0, YU));
445*4a221d59SStefano Zampini             } else {
446*4a221d59SStefano Zampini               PetscCall(VecAXPBY(Y, 1.0, tau - 1, YU));
447*4a221d59SStefano Zampini             }
448*4a221d59SStefano Zampini           }
449*4a221d59SStefano Zampini           PetscCall(VecNorm(Y, NORM_2, &c0)); /* this norm will be cached and reused later */
450*4a221d59SStefano Zampini           PetscCall(PetscInfo(snes, "%s evaluated. roots: (%g, %g), tau %g, ynorm: %g, ycnorm: %g, ydlnorm %g, gTBg: %g\n", noroots ? "CP" : "DL", (double)tneg, (double)tpos, (double)tau, (double)ynorm, (double)ycnorm, (double)c0, (double)gTBg));
451*4a221d59SStefano Zampini         }
452*4a221d59SStefano Zampini         break;
453*4a221d59SStefano Zampini       default:
454*4a221d59SStefano Zampini         SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_SUP, "Unknown fallback mode");
455454a90a3SBarry Smith         break;
45652392280SLois Curfman McInnes       }
4574800dd8cSBarry Smith     }
458*4a221d59SStefano Zampini 
459*4a221d59SStefano Zampini     /* Evaluate the solution to meet the improvement ratio criteria */
460*4a221d59SStefano Zampini 
461*4a221d59SStefano Zampini     /* compute the final ynorm */
462*4a221d59SStefano Zampini     PetscCall(SNESNewtonTRPreCheck(snes, X, Y, &changed_y));
463*4a221d59SStefano Zampini     PetscCall(VecNorm(Y, NORM_2, &ynorm));
464*4a221d59SStefano Zampini 
465*4a221d59SStefano Zampini     /* the quadratic model difference */
466*4a221d59SStefano Zampini     PetscCall(MatMult(snes->jacobian, Y, W));
467*4a221d59SStefano Zampini     if (objective) PetscCall(VecDotRealPart(Y, W, &yTHy));
468*4a221d59SStefano Zampini     else PetscCall(VecDotRealPart(W, W, &yTHy)); /* Gauss-Newton approximation J^t * J */
469*4a221d59SStefano Zampini     PetscCall(VecDotRealPart(GradF, Y, &gTy));
470*4a221d59SStefano Zampini     deltaqm = -(-gTy + 0.5 * yTHy); /* difference in quadratic model, -gTy because SNES solves it this way */
471*4a221d59SStefano Zampini 
472*4a221d59SStefano Zampini     /* update */
473*4a221d59SStefano Zampini     PetscCall(VecWAXPY(W, -1.0, Y, X)); /* Xkp1 */
474*4a221d59SStefano Zampini     PetscCall(SNESNewtonTRPostCheck(snes, X, Y, W, &changed_y, &changed_w));
475*4a221d59SStefano Zampini     if (changed_y) {
476*4a221d59SStefano Zampini       /* Need to recompute the quadratic model difference */
477*4a221d59SStefano Zampini       PetscCall(MatMult(snes->jacobian, Y, W));
478*4a221d59SStefano Zampini       if (objective) PetscCall(VecDotRealPart(Y, W, &yTHy));
479*4a221d59SStefano Zampini       else PetscCall(VecDotRealPart(W, W, &yTHy));
480*4a221d59SStefano Zampini       PetscCall(VecDotRealPart(GradF, Y, &gTy));
481*4a221d59SStefano Zampini       deltaqm = -(-gTy + 0.5 * yTHy);
482*4a221d59SStefano Zampini       /* User changed Y but not W */
483*4a221d59SStefano Zampini       if (!changed_w) PetscCall(VecWAXPY(W, -1.0, Y, X));
484*4a221d59SStefano Zampini     }
485*4a221d59SStefano Zampini 
486*4a221d59SStefano Zampini     /* Compute new objective function */
487*4a221d59SStefano Zampini     PetscCall(SNESComputeFunction(snes, W, G)); /*  F(Xkp1) = G */
488*4a221d59SStefano Zampini     PetscCall(VecNorm(G, NORM_2, &gnorm));
489*4a221d59SStefano Zampini     if (objective) PetscCall(SNESComputeObjective(snes, W, &fkp1));
490*4a221d59SStefano Zampini     else fkp1 = 0.5 * PetscSqr(gnorm);
491*4a221d59SStefano Zampini     SNESCheckFunctionNorm(snes, fkp1);
492*4a221d59SStefano Zampini 
493*4a221d59SStefano Zampini     if (deltaqm > 0.0) rho = (fk - fkp1) / deltaqm; /* actual improvement over predicted improvement */
494*4a221d59SStefano Zampini     else rho = -1.0;                                /*  no reduction in quadratic model, step must be rejected */
495*4a221d59SStefano Zampini     PetscCall(PetscInfo(snes, "rho=%g, delta=%g, fk=%g, fkp1=%g, deltaqm=%g, gTy=%g, yTHy=%g\n", (double)rho, (double)delta, (double)fk, (double)fkp1, (double)deltaqm, (double)gTy, (double)yTHy));
496*4a221d59SStefano Zampini 
497*4a221d59SStefano Zampini     if (rho < neP->eta2) delta *= neP->t1;      /* shrink the region */
498*4a221d59SStefano Zampini     else if (rho > neP->eta3) delta *= neP->t2; /* expand the region */
499*4a221d59SStefano Zampini     delta = PetscMin(delta, deltaM);            /* but not greater than deltaM */
500*4a221d59SStefano Zampini 
501*4a221d59SStefano Zampini     neP->delta = delta;
502*4a221d59SStefano Zampini     if (rho >= neP->eta1) {
503*4a221d59SStefano Zampini       rho_satisfied = PETSC_TRUE;
504*4a221d59SStefano Zampini     } else {
505*4a221d59SStefano Zampini       rho_satisfied = PETSC_FALSE;
506*4a221d59SStefano Zampini       PetscCall(PetscInfo(snes, "Trying again in smaller region\n"));
507*4a221d59SStefano Zampini       /* check to see if progress is hopeless */
508*4a221d59SStefano Zampini       PetscCall(SNESTR_Converged_Private(snes, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP));
509*4a221d59SStefano Zampini       if (!snes->reason) PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
510*4a221d59SStefano Zampini       if (snes->reason == SNES_CONVERGED_SNORM_RELATIVE) snes->reason = SNES_DIVERGED_INNER;
511*4a221d59SStefano Zampini       snes->numFailures++;
512*4a221d59SStefano Zampini       /* We're not progressing, so return with the current iterate */
513*4a221d59SStefano Zampini       if (snes->reason) break;
514*4a221d59SStefano Zampini     }
515*4a221d59SStefano Zampini     if (rho_satisfied) {
516*4a221d59SStefano Zampini       /* Update function values */
517*4a221d59SStefano Zampini       already_done = PETSC_FALSE;
5184800dd8cSBarry Smith       fnorm        = gnorm;
519*4a221d59SStefano Zampini       fk           = fkp1;
520*4a221d59SStefano Zampini 
521*4a221d59SStefano Zampini       /* New residual and linearization point */
5229566063dSJacob Faibussowitsch       PetscCall(VecCopy(G, F));
5239566063dSJacob Faibussowitsch       PetscCall(VecCopy(W, X));
524*4a221d59SStefano Zampini 
52585385478SLisandro Dalcin       /* Monitor convergence */
5269566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
527*4a221d59SStefano Zampini       snes->iter++;
528fbe28522SBarry Smith       snes->norm  = fnorm;
529c1e67a49SFande Kong       snes->xnorm = xnorm;
530c1e67a49SFande Kong       snes->ynorm = ynorm;
5319566063dSJacob Faibussowitsch       PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
5329566063dSJacob Faibussowitsch       PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
5339566063dSJacob Faibussowitsch       PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
534*4a221d59SStefano Zampini 
53585385478SLisandro Dalcin       /* Test for convergence, xnorm = || X || */
536*4a221d59SStefano Zampini       PetscCall(VecNorm(X, NORM_2, &xnorm));
537*4a221d59SStefano Zampini       PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
538*4a221d59SStefano Zampini       if (snes->reason) break;
539*4a221d59SStefano Zampini     }
54038442cffSBarry Smith   }
541284fb49fSHeeho Park 
542*4a221d59SStefano Zampini   if (snes->iter == maxits) {
5439566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
544*4a221d59SStefano Zampini     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
54552392280SLois Curfman McInnes   }
546*4a221d59SStefano Zampini   if (clear_converged_test) {
5479566063dSJacob Faibussowitsch     PetscCall(KSPGetAndClearConvergenceTest(ksp, &ctx->convtest, &ctx->convctx, &ctx->convdestroy));
5489566063dSJacob Faibussowitsch     PetscCall(PetscFree(ctx));
5499566063dSJacob Faibussowitsch     PetscCall(KSPSetConvergenceTest(ksp, convtest, convctx, convdestroy));
5505e28bcb6Sprj-   }
5513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5524800dd8cSBarry Smith }
553284fb49fSHeeho Park 
554d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
555d71ae5a4SJacob Faibussowitsch {
5563a40ed3dSBarry Smith   PetscFunctionBegin;
5579566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 4));
5589566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
5593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5604800dd8cSBarry Smith }
5616b8b9a38SLisandro Dalcin 
562d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
563d71ae5a4SJacob Faibussowitsch {
5643a40ed3dSBarry Smith   PetscFunctionBegin;
5659566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
5663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5674800dd8cSBarry Smith }
5684800dd8cSBarry Smith 
569d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONTR(SNES snes, PetscOptionItems *PetscOptionsObject)
570d71ae5a4SJacob Faibussowitsch {
57104d7464bSBarry Smith   SNES_NEWTONTR *ctx = (SNES_NEWTONTR *)snes->data;
5724800dd8cSBarry Smith 
5733a40ed3dSBarry Smith   PetscFunctionBegin;
574d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES trust region options for nonlinear equations");
575*4a221d59SStefano Zampini   PetscCall(PetscOptionsReal("-snes_tr_tol", "Trust region tolerance", "SNESSetTrustRegionTolerance", snes->deltatol, &snes->deltatol, NULL));
576*4a221d59SStefano Zampini   PetscCall(PetscOptionsReal("-snes_tr_eta1", "eta1", "None", ctx->eta1, &ctx->eta1, NULL));
577*4a221d59SStefano Zampini   PetscCall(PetscOptionsReal("-snes_tr_eta2", "eta2", "None", ctx->eta2, &ctx->eta2, NULL));
578*4a221d59SStefano Zampini   PetscCall(PetscOptionsReal("-snes_tr_eta3", "eta3", "None", ctx->eta3, &ctx->eta3, NULL));
579*4a221d59SStefano Zampini   PetscCall(PetscOptionsReal("-snes_tr_t1", "t1", "None", ctx->t1, &ctx->t1, NULL));
580*4a221d59SStefano Zampini   PetscCall(PetscOptionsReal("-snes_tr_t2", "t2", "None", ctx->t2, &ctx->t2, NULL));
581*4a221d59SStefano Zampini   PetscCall(PetscOptionsReal("-snes_tr_deltaM", "deltaM", "None", ctx->deltaM, &ctx->deltaM, NULL));
5829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_tr_delta0", "delta0", "None", ctx->delta0, &ctx->delta0, NULL));
583*4a221d59SStefano Zampini   PetscCall(PetscOptionsEnum("-snes_tr_fallback_type", "Type of fallback if subproblem solution is outside of the trust region", "SNESNewtonTRSetFallbackType", SNESNewtonTRFallbackTypes, (PetscEnum)ctx->fallback, (PetscEnum *)&ctx->fallback, NULL));
584d0609cedSBarry Smith   PetscOptionsHeadEnd();
5853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5864800dd8cSBarry Smith }
5874800dd8cSBarry Smith 
588d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONTR(SNES snes, PetscViewer viewer)
589d71ae5a4SJacob Faibussowitsch {
59004d7464bSBarry Smith   SNES_NEWTONTR *tr = (SNES_NEWTONTR *)snes->data;
591ace3abfcSBarry Smith   PetscBool      iascii;
592a935fc98SLois Curfman McInnes 
5933a40ed3dSBarry Smith   PetscFunctionBegin;
5949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
59532077d6dSBarry Smith   if (iascii) {
596*4a221d59SStefano Zampini     PetscCall(PetscViewerASCIIPrintf(viewer, "  Trust region tolerance %g\n", (double)snes->deltatol));
597*4a221d59SStefano Zampini     PetscCall(PetscViewerASCIIPrintf(viewer, "  eta1=%g, eta2=%g, eta3=%g\n", (double)tr->eta1, (double)tr->eta2, (double)tr->eta3));
598*4a221d59SStefano Zampini     PetscCall(PetscViewerASCIIPrintf(viewer, "  delta0=%g, t1=%g, t2=%g, deltaM=%g\n", (double)tr->delta0, (double)tr->t1, (double)tr->t2, (double)tr->deltaM));
599*4a221d59SStefano Zampini     PetscCall(PetscViewerASCIIPrintf(viewer, "  fallback=%s\n", SNESNewtonTRFallbackTypes[tr->fallback]));
60019bcc07fSBarry Smith   }
6013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
602a935fc98SLois Curfman McInnes }
603f6dfbefdSBarry Smith 
604ebe3b25bSBarry Smith /*MC
605*4a221d59SStefano Zampini       SNESNEWTONTR - Newton based nonlinear solver that uses trust-region dogleg method with Cauchy direction
606f6dfbefdSBarry Smith 
607f6dfbefdSBarry Smith    Options Database Keys:
608*4a221d59SStefano Zampini +   -snes_tr_tol <tol> - trust region tolerance
609*4a221d59SStefano Zampini .   -snes_tr_eta1 <eta1> - trust region parameter 0.0 <= eta1 <= eta2, rho >= eta1 breaks out of the inner iteration (default: eta1=0.001)
610*4a221d59SStefano Zampini .   -snes_tr_eta2 <eta2> - trust region parameter 0.0 <= eta1 <= eta2, rho <= eta2 shrinks the trust region (default: eta2=0.25)
611*4a221d59SStefano Zampini .   -snes_tr_eta3 <eta3> - trust region parameter eta3 > eta2, rho >= eta3 expands the trust region (default: eta3=0.75)
612*4a221d59SStefano Zampini .   -snes_tr_t1 <t1> - trust region parameter, shrinking factor of trust region (default: 0.25)
613*4a221d59SStefano Zampini .   -snes_tr_t2 <t2> - trust region parameter, expanding factor of trust region (default: 2.0)
614*4a221d59SStefano Zampini .   -snes_tr_deltaM <deltaM> - trust region parameter, max size of trust region (default: MAX_REAL)
615*4a221d59SStefano Zampini .   -snes_tr_delta0 <delta0> - trust region parameter, initial size of trust region (default: 0.2)
616*4a221d59SStefano Zampini -   -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.
617ebe3b25bSBarry Smith 
618f6dfbefdSBarry Smith     Reference:
619*4a221d59SStefano Zampini .   * - "Numerical Optimization" by Nocedal and Wright, chapter 4.
620ebe3b25bSBarry Smith 
621ee3001cbSBarry Smith    Level: intermediate
622ee3001cbSBarry Smith 
623*4a221d59SStefano Zampini .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESSetTrustRegionTolerance()`,
624*4a221d59SStefano Zampini           `SNESNewtonTRPreCheck()`, `SNESNewtonTRGetPreCheck()`, `SNESNewtonTRSetPostCheck()`, `SNESNewtonTRGetPostCheck()`,
625*4a221d59SStefano Zampini           `SNESNewtonTRSetPreCheck()`, `SNESNewtonTRSetFallbackType()`
626ebe3b25bSBarry Smith M*/
627d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
628d71ae5a4SJacob Faibussowitsch {
62904d7464bSBarry Smith   SNES_NEWTONTR *neP;
6304800dd8cSBarry Smith 
6313a40ed3dSBarry Smith   PetscFunctionBegin;
63204d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONTR;
63304d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONTR;
63404d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONTR;
63504d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
63604d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONTR;
637fbe28522SBarry Smith 
63842f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
639efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
64042f4f86dSBarry Smith 
6414fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
6424fc747eaSLawrence Mitchell 
6434dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
644fbe28522SBarry Smith   snes->data    = (void *)neP;
645fbe28522SBarry Smith   neP->delta    = 0.0;
646fbe28522SBarry Smith   neP->delta0   = 0.2;
647*4a221d59SStefano Zampini   neP->eta1     = 0.001;
648*4a221d59SStefano Zampini   neP->eta2     = 0.25;
649*4a221d59SStefano Zampini   neP->eta3     = 0.75;
650*4a221d59SStefano Zampini   neP->t1       = 0.25;
651*4a221d59SStefano Zampini   neP->t2       = 2.0;
652*4a221d59SStefano Zampini   neP->deltaM   = 1.e10;
653*4a221d59SStefano Zampini   neP->fallback = SNES_TR_FALLBACK_NEWTON;
6543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6554800dd8cSBarry Smith }
656