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