xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision e91c04dfc8a52dee1965211bb1cc8e5bf775178f)
1 #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
2 #include <petscblaslapack.h>
3 #include <petscdm.h>
4 
5 const char *const SNESNGMRESRestartTypes[] = {"NONE", "PERIODIC", "DIFFERENCE", "SNESNGMRESRestartType", "SNES_NGMRES_RESTART_", NULL};
6 const char *const SNESNGMRESSelectTypes[]  = {"NONE", "DIFFERENCE", "LINESEARCH", "SNESNGMRESSelectType", "SNES_NGMRES_SELECT_", NULL};
7 
8 PetscErrorCode SNESReset_NGMRES(SNES snes)
9 {
10   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
11 
12   PetscFunctionBegin;
13   PetscCall(VecDestroyVecs(ngmres->msize, &ngmres->Fdot));
14   PetscCall(VecDestroyVecs(ngmres->msize, &ngmres->Xdot));
15   PetscCall(SNESLineSearchDestroy(&ngmres->additive_linesearch));
16   PetscFunctionReturn(PETSC_SUCCESS);
17 }
18 
19 PetscErrorCode SNESDestroy_NGMRES(SNES snes)
20 {
21   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
22 
23   PetscFunctionBegin;
24   PetscCall(SNESReset_NGMRES(snes));
25   PetscCall(PetscFree4(ngmres->h, ngmres->beta, ngmres->xi, ngmres->q));
26   PetscCall(PetscFree3(ngmres->xnorms, ngmres->fnorms, ngmres->s));
27 #if defined(PETSC_USE_COMPLEX)
28   PetscCall(PetscFree(ngmres->rwork));
29 #endif
30   PetscCall(PetscFree(ngmres->work));
31   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", NULL));
32   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", NULL));
33   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", NULL));
34   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", NULL));
35   PetscCall(PetscFree(snes->data));
36   PetscFunctionReturn(PETSC_SUCCESS);
37 }
38 
39 PetscErrorCode SNESSetUp_NGMRES(SNES snes)
40 {
41   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
42   PetscInt     msize, hsize;
43   DM           dm;
44 
45   PetscFunctionBegin;
46   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
47     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNGMRES does not support left preconditioning with unpreconditioned function");
48   }
49   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
50   PetscCall(SNESSetWorkVecs(snes, 5));
51 
52   if (!snes->vec_sol) {
53     PetscCall(SNESGetDM(snes, &dm));
54     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
55   }
56 
57   if (!ngmres->Xdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Xdot));
58   if (!ngmres->Fdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Fdot));
59   if (!ngmres->setup_called) {
60     msize = ngmres->msize; /* restart size */
61     hsize = msize * msize;
62 
63     /* explicit least squares minimization solve */
64     PetscCall(PetscCalloc4(hsize, &ngmres->h, msize, &ngmres->beta, msize, &ngmres->xi, hsize, &ngmres->q));
65     PetscCall(PetscMalloc3(msize, &ngmres->xnorms, msize, &ngmres->fnorms, msize, &ngmres->s));
66     ngmres->nrhs = 1;
67     PetscCall(PetscBLASIntCast(msize, &ngmres->lda));
68     PetscCall(PetscBLASIntCast(msize, &ngmres->ldb));
69     PetscCall(PetscBLASIntCast(12 * msize, &ngmres->lwork));
70 #if defined(PETSC_USE_COMPLEX)
71     PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->rwork));
72 #endif
73     PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->work));
74   }
75 
76   ngmres->setup_called = PETSC_TRUE;
77   PetscFunctionReturn(PETSC_SUCCESS);
78 }
79 
80 static PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes, PetscOptionItems *PetscOptionsObject)
81 {
82   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
83   PetscBool    debug  = PETSC_FALSE;
84 
85   PetscFunctionBegin;
86   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NGMRES options");
87   PetscCall(PetscOptionsEnum("-snes_ngmres_select_type", "Select type", "SNESNGMRESSetSelectType", SNESNGMRESSelectTypes, (PetscEnum)ngmres->select_type, (PetscEnum *)&ngmres->select_type, NULL));
88   PetscCall(PetscOptionsEnum("-snes_ngmres_restart_type", "Restart type", "SNESNGMRESSetRestartType", SNESNGMRESRestartTypes, (PetscEnum)ngmres->restart_type, (PetscEnum *)&ngmres->restart_type, NULL));
89   PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES", ngmres->candidate, &ngmres->candidate, NULL));
90   PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc", "Linearly approximate the function", "SNES", ngmres->approxfunc, &ngmres->approxfunc, NULL));
91   PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, NULL));
92   PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES", ngmres->restart_periodic, &ngmres->restart_periodic, NULL));
93   PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart", "SNES", ngmres->restart_it, &ngmres->restart_it, NULL));
94   PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
95   if (debug) ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
96   PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, NULL));
97   PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, NULL));
98   PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, NULL));
99   PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, NULL));
100   PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise", ngmres->restart_fm_rise, &ngmres->restart_fm_rise, NULL));
101   PetscOptionsHeadEnd();
102   if (ngmres->gammaA > ngmres->gammaC && ngmres->gammaC > 2.) ngmres->gammaC = ngmres->gammaA;
103   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
104     PetscCall(SNESNGMRESGetAdditiveLineSearch_Private(snes, &ngmres->additive_linesearch));
105     PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch));
106   }
107   PetscFunctionReturn(PETSC_SUCCESS);
108 }
109 
110 PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
111 {
112   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
113   PetscBool    iascii;
114 
115   PetscFunctionBegin;
116   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
117   if (iascii) {
118     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of stored past updates: %" PetscInt_FMT "\n", ngmres->msize));
119     if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
120       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", (double)ngmres->gammaA, (double)ngmres->gammaC));
121       PetscCall(PetscViewerASCIIPrintf(viewer, "  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", (double)ngmres->epsilonB, (double)ngmres->deltaB));
122       PetscCall(PetscViewerASCIIPrintf(viewer, "  Restart on F_M residual increase: %s\n", PetscBools[ngmres->restart_fm_rise]));
123     }
124     if (ngmres->additive_linesearch) {
125       PetscCall(PetscViewerASCIIPrintf(viewer, "  Additive line-search details:\n"));
126       PetscCall(SNESLineSearchView(ngmres->additive_linesearch, viewer));
127     }
128   }
129   PetscFunctionReturn(PETSC_SUCCESS);
130 }
131 
132 static PetscErrorCode SNESSolve_NGMRES(SNES snes)
133 {
134   SNES_NGMRES         *ngmres = (SNES_NGMRES *)snes->data;
135   Vec                  X, F, B, D, Y;         /* present solution, residual, and preconditioned residual */
136   Vec                  XA, FA, XM, FM;        /* candidate linear combination answers */
137   PetscReal            fnorm, fMnorm, fAnorm; /* coefficients and RHS to the minimization problem */
138   PetscReal            xnorm, xMnorm, xAnorm;
139   PetscReal            ynorm, yMnorm, yAnorm;
140   PetscInt             k, k_restart, l, ivec, restart_count = 0;
141   PetscReal            objmin, objM, objA, obj; /* support for objective functions minimization */
142   PetscBool            selectRestart;           /* solution selection data */
143   SNESConvergedReason  reason;
144   SNESLineSearchReason lssucceed;
145   SNESObjectiveFn     *objective;
146   /*
147       These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed
148       to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case
149       so the code is correct as written.
150   */
151   PetscReal dnorm = 0.0, dminnorm = 0.0;
152 
153   PetscFunctionBegin;
154   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);
155 
156   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
157   /* variable initialization */
158   snes->reason = SNES_CONVERGED_ITERATING;
159   X            = snes->vec_sol;
160   F            = snes->vec_func;
161   B            = snes->vec_rhs;
162   XA           = snes->work[2];
163   FA           = snes->work[0];
164   D            = snes->work[1];
165 
166   /* work for the line search */
167   Y  = snes->vec_sol_update;
168   XM = snes->work[3];
169   FM = snes->work[4];
170 
171   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
172   snes->iter = 0;
173   snes->norm = 0.;
174   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
175 
176   /* initialization */
177 
178   if (snes->npc && snes->npcside == PC_LEFT) {
179     PetscCall(SNESApplyNPC(snes, X, NULL, F));
180     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
181     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
182       snes->reason = SNES_DIVERGED_INNER;
183       PetscFunctionReturn(PETSC_SUCCESS);
184     }
185     PetscCall(VecNorm(F, NORM_2, &fnorm));
186   } else {
187     if (!snes->vec_func_init_set) {
188       PetscCall(SNESComputeFunction(snes, X, F));
189     } else snes->vec_func_init_set = PETSC_FALSE;
190 
191     PetscCall(VecNorm(F, NORM_2, &fnorm));
192     SNESCheckFunctionNorm(snes, fnorm);
193   }
194   PetscCall(SNESGetObjective(snes, &objective, NULL));
195   objmin = fnorm;
196   if (objective) PetscCall(SNESComputeObjective(snes, X, &objmin));
197   obj = objmin;
198 
199   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
200   snes->norm = fnorm;
201   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
202   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
203   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
204   PetscCall(SNESMonitor(snes, 0, fnorm));
205   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
206   PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, F, fnorm, X));
207 
208   k_restart = 1;
209   l         = 1;
210   ivec      = 0;
211   for (k = 1; k < snes->max_its + 1; k++) {
212     /* Call general purpose update function */
213     PetscTryTypeMethod(snes, update, snes->iter);
214 
215     /* Computation of x^M */
216     if (snes->npc && snes->npcside == PC_RIGHT) {
217       PetscCall(VecCopy(X, XM));
218       PetscCall(SNESSetInitialFunction(snes->npc, F));
219 
220       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0));
221       PetscCall(SNESSolve(snes->npc, B, XM));
222       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0));
223 
224       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
225       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
226         snes->reason = SNES_DIVERGED_INNER;
227         PetscFunctionReturn(PETSC_SUCCESS);
228       }
229       PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm));
230     } else {
231       /* no preconditioner -- just take gradient descent with line search */
232       PetscCall(VecCopy(F, Y));
233       PetscCall(VecCopy(F, FM));
234       PetscCall(VecCopy(X, XM));
235 
236       fMnorm = fnorm;
237 
238       PetscCall(SNESLineSearchApply(snes->linesearch, XM, FM, &fMnorm, Y));
239       PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
240       if (lssucceed) {
241         if (++snes->numFailures >= snes->maxFailures) {
242           snes->reason = SNES_DIVERGED_LINE_SEARCH;
243           PetscFunctionReturn(PETSC_SUCCESS);
244         }
245       }
246     }
247     if (objective) PetscCall(SNESComputeObjective(snes, XM, &objM));
248     else objM = fMnorm;
249     objmin = PetscMin(objmin, objM);
250 
251     PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA));
252 
253     /* differences for selection and restart */
254     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
255       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
256     } else {
257       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
258     }
259     if (objective) PetscCall(SNESComputeObjective(snes, XA, &objA));
260     else objA = fAnorm;
261     SNESCheckFunctionNorm(snes, fnorm);
262 
263     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
264     PetscCall(SNESNGMRESSelect_Private(snes, k_restart, XM, FM, xMnorm, fMnorm, yMnorm, objM, XA, FA, xAnorm, fAnorm, yAnorm, objA, dnorm, objmin, dminnorm, X, F, Y, &xnorm, &fnorm, &ynorm));
265     if (objective) PetscCall(SNESComputeObjective(snes, X, &obj));
266     else obj = fnorm;
267     selectRestart = PETSC_FALSE;
268 
269     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
270       PetscCall(SNESNGMRESSelectRestart_Private(snes, l, obj, objM, objA, dnorm, objmin, dminnorm, &selectRestart));
271 
272       /* if the restart conditions persist for more than restart_it iterations, restart. */
273       if (selectRestart) restart_count++;
274       else restart_count = 0;
275     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
276       if (k_restart > ngmres->restart_periodic) {
277         if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart));
278         restart_count = ngmres->restart_it;
279       }
280     }
281 
282     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
283 
284     /* restart after restart conditions have persisted for a fixed number of iterations */
285     if (restart_count >= ngmres->restart_it) {
286       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart));
287       restart_count = 0;
288       k_restart     = 1;
289       l             = 1;
290       ivec          = 0;
291       /* q_{00} = nu */
292       PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, FM, fMnorm, XM));
293     } else {
294       /* select the current size of the subspace */
295       if (l < ngmres->msize) l++;
296       k_restart++;
297       /* place the current entry in the list of previous entries */
298       if (ngmres->candidate) {
299         objmin = PetscMin(objmin, objM);
300         PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fMnorm, XM));
301       } else {
302         objmin = PetscMin(objmin, obj);
303         PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, F, fnorm, X));
304       }
305     }
306 
307     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
308     snes->iter  = k;
309     snes->norm  = fnorm;
310     snes->ynorm = ynorm;
311     snes->xnorm = xnorm;
312     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
313     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
314     PetscCall(SNESConverged(snes, snes->iter, 0, 0, fnorm));
315     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
316     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
317   }
318   snes->reason = SNES_DIVERGED_MAX_IT;
319   PetscFunctionReturn(PETSC_SUCCESS);
320 }
321 
322 /*@
323   SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M inside a `SNESNGMRES` solve
324 
325   Input Parameters:
326 + snes - the `SNES` context.
327 - flg  - boolean value deciding whether to use the option or not, default is `PETSC_FALSE`
328 
329   Options Database Key:
330 . -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M
331 
332   Level: advanced
333 
334   Notes:
335   If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area.
336   To help the solver do that, remove the current stored solutions and residuals whenever F_M increases.
337 
338   This option must be used with the `SNESNGMRES` `SNESNGMRESRestartType` of `SNES_NGMRES_RESTART_DIFFERENCE`
339 
340 .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartType()`
341   @*/
342 PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes, PetscBool flg)
343 {
344   PetscErrorCode (*f)(SNES, PetscBool);
345 
346   PetscFunctionBegin;
347   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", &f));
348   if (f) PetscCall((f)(snes, flg));
349   PetscFunctionReturn(PETSC_SUCCESS);
350 }
351 
352 static PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes, PetscBool flg)
353 {
354   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
355 
356   PetscFunctionBegin;
357   ngmres->restart_fm_rise = flg;
358   PetscFunctionReturn(PETSC_SUCCESS);
359 }
360 
361 PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes, PetscBool *flg)
362 {
363   PetscErrorCode (*f)(SNES, PetscBool *);
364 
365   PetscFunctionBegin;
366   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", &f));
367   if (f) PetscCall((f)(snes, flg));
368   PetscFunctionReturn(PETSC_SUCCESS);
369 }
370 
371 static PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes, PetscBool *flg)
372 {
373   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
374 
375   PetscFunctionBegin;
376   *flg = ngmres->restart_fm_rise;
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379 
380 /*@
381   SNESNGMRESSetRestartType - Sets the restart type for `SNESNGMRES`.
382 
383   Logically Collective
384 
385   Input Parameters:
386 + snes  - the iterative context
387 - rtype - restart type, see `SNESNGMRESRestartType`
388 
389   Options Database Keys:
390 + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
391 - -snes_ngmres_restart <30>                           - sets the number of iterations before restart for periodic
392 
393   Level: intermediate
394 
395 .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartFmRise()`,
396           `SNESNGMRESSetSelectType()`
397 @*/
398 PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype)
399 {
400   PetscFunctionBegin;
401   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
402   PetscTryMethod(snes, "SNESNGMRESSetRestartType_C", (SNES, SNESNGMRESRestartType), (snes, rtype));
403   PetscFunctionReturn(PETSC_SUCCESS);
404 }
405 
406 /*@
407   SNESNGMRESSetSelectType - Sets the selection type for `SNESNGMRES`.  This determines how the candidate solution and
408   combined solution are used to create the next iterate.
409 
410   Logically Collective
411 
412   Input Parameters:
413 + snes  - the iterative context
414 - stype - selection type, see `SNESNGMRESSelectType`
415 
416   Options Database Key:
417 . -snes_ngmres_select_type<difference,none,linesearch> - select type
418 
419   Level: intermediate
420 
421   Note:
422   The default line search used is the `SNESLINESEARCHL2` line search and it requires two additional function evaluations.
423 
424 .seealso: [](ch_snes), `SNES`, `SNESNGMRES`, `SNESNGMRESSelectType`, `SNES_NGMRES_SELECT_NONE`, `SNES_NGMRES_SELECT_DIFFERENCE`, `SNES_NGMRES_SELECT_LINESEARCH`,
425           `SNESNGMRESSetRestartType()`
426 @*/
427 PetscErrorCode SNESNGMRESSetSelectType(SNES snes, SNESNGMRESSelectType stype)
428 {
429   PetscFunctionBegin;
430   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
431   PetscTryMethod(snes, "SNESNGMRESSetSelectType_C", (SNES, SNESNGMRESSelectType), (snes, stype));
432   PetscFunctionReturn(PETSC_SUCCESS);
433 }
434 
435 static PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype)
436 {
437   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
438 
439   PetscFunctionBegin;
440   ngmres->select_type = stype;
441   PetscFunctionReturn(PETSC_SUCCESS);
442 }
443 
444 static PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype)
445 {
446   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
447 
448   PetscFunctionBegin;
449   ngmres->restart_type = rtype;
450   PetscFunctionReturn(PETSC_SUCCESS);
451 }
452 
453 /*MC
454   SNESNGMRES - An implementation of the Nonlinear Generalized Minimum Residual method, Nonlinear GMRES, or N-GMRES {cite}`ow1`, {cite}`bruneknepleysmithtu15` for solving
455                nonlinear systems with `SNES`.
456 
457    Level: beginner
458 
459    Options Database Keys:
460 +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
461 .  -snes_ngmres_restart_type<difference,none,periodic>  - choose the restart conditions
462 .  -snes_ngmres_candidate                               - Use `SNESNGMRES` variant which combines candidate solutions instead of actual solutions
463 .  -snes_ngmres_m                                       - Number of stored previous solutions and residuals
464 .  -snes_ngmres_restart_it                              - Number of iterations the restart conditions hold before restart
465 .  -snes_ngmres_gammaA                                  - Residual tolerance for solution select between the candidate and combination
466 .  -snes_ngmres_gammaC                                  - Residual tolerance for restart
467 .  -snes_ngmres_epsilonB                                - Difference tolerance between subsequent solutions triggering restart
468 .  -snes_ngmres_deltaB                                  - Difference tolerance between residuals triggering restart
469 .  -snes_ngmres_restart_fm_rise                         - Restart on residual rise from $x_M$ step
470 .  -snes_ngmres_monitor                                 - Prints relevant information about the nonlinear GNMRES iterations
471 .  -snes_linesearch_type <basic,l2,cp>                  - Line search type used for the default smoother
472 -  -snes_ngmres_additive_snes_linesearch_type           - line search type used to select between the candidate and combined solution with additive select type
473 
474    Notes:
475    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
476    optimization problem at each iteration.
477 
478    Very similar to the `SNESANDERSON` algorithm.
479 
480    Unlike the linear GMRES algorithm this algorithm does not compute a Krylov subspace using the Arnoldi process. Instead it stores a
481    collection of previous solutions and the residuals $ F(x) - b $ at those solutions.
482 
483    This algorithm ignores any Jacobian provided with `SNESSetJacobian()`
484 
485 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESANDERSON`, `SNESNGMRESSetSelectType()`, `SNESNGMRESSetRestartType()`,
486           `SNESNGMRESSetRestartFmRise()`, `SNESNGMRESSelectType`, `SNESNGMRESRestartType`
487 M*/
488 
489 PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
490 {
491   SNES_NGMRES   *ngmres;
492   SNESLineSearch linesearch;
493 
494   PetscFunctionBegin;
495   snes->ops->destroy        = SNESDestroy_NGMRES;
496   snes->ops->setup          = SNESSetUp_NGMRES;
497   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
498   snes->ops->view           = SNESView_NGMRES;
499   snes->ops->solve          = SNESSolve_NGMRES;
500   snes->ops->reset          = SNESReset_NGMRES;
501 
502   snes->usesnpc = PETSC_TRUE;
503   snes->usesksp = PETSC_FALSE;
504   snes->npcside = PC_RIGHT;
505 
506   snes->alwayscomputesfinalresidual = PETSC_TRUE;
507 
508   PetscCall(PetscNew(&ngmres));
509   snes->data    = (void *)ngmres;
510   ngmres->msize = 30;
511 
512   PetscCall(SNESParametersInitialize(snes));
513   PetscObjectParameterSetDefault(snes, max_funcs, 30000);
514   PetscObjectParameterSetDefault(snes, max_its, 10000);
515 
516   ngmres->candidate = PETSC_FALSE;
517 
518   PetscCall(SNESGetLineSearch(snes, &linesearch));
519   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
520 
521   ngmres->additive_linesearch = NULL;
522   ngmres->approxfunc          = PETSC_FALSE;
523   ngmres->restart_it          = 2;
524   ngmres->restart_periodic    = 30;
525   ngmres->gammaA              = 2.0;
526   ngmres->gammaC              = 2.0;
527   ngmres->deltaB              = 0.9;
528   ngmres->epsilonB            = 0.1;
529   ngmres->restart_fm_rise     = PETSC_FALSE;
530 
531   ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
532   ngmres->select_type  = SNES_NGMRES_SELECT_DIFFERENCE;
533 
534   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", SNESNGMRESSetSelectType_NGMRES));
535   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", SNESNGMRESSetRestartType_NGMRES));
536   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", SNESNGMRESSetRestartFmRise_NGMRES));
537   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", SNESNGMRESGetRestartFmRise_NGMRES));
538   PetscFunctionReturn(PETSC_SUCCESS);
539 }
540