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