xref: /petsc/src/snes/impls/ngmres/snesngmres.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
113a62661SPeter Brune #include <../src/snes/impls/ngmres/snesngmres.h> /*I "petscsnes.h" I*/
219653cdaSPeter Brune #include <petscblaslapack.h>
3fced5a79SAsbjørn Nilsen Riseth #include <petscdm.h>
4a312c225SMatthew G Knepley 
59e5d0892SLisandro Dalcin const char *const SNESNGMRESRestartTypes[] = {"NONE", "PERIODIC", "DIFFERENCE", "SNESNGMRESRestartType", "SNES_NGMRES_RESTART_", NULL};
69e5d0892SLisandro Dalcin const char *const SNESNGMRESSelectTypes[]  = {"NONE", "DIFFERENCE", "LINESEARCH", "SNESNGMRESSelectType", "SNES_NGMRES_SELECT_", NULL};
713a62661SPeter Brune 
8d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NGMRES(SNES snes)
9d71ae5a4SJacob Faibussowitsch {
10a312c225SMatthew G Knepley   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
11a312c225SMatthew G Knepley 
12a312c225SMatthew G Knepley   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ngmres->msize, &ngmres->Fdot));
149566063dSJacob Faibussowitsch   PetscCall(VecDestroyVecs(ngmres->msize, &ngmres->Xdot));
159566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchDestroy(&ngmres->additive_linesearch));
163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17a312c225SMatthew G Knepley }
18a312c225SMatthew G Knepley 
19d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NGMRES(SNES snes)
20d71ae5a4SJacob Faibussowitsch {
2178440776SJed Brown   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
22a312c225SMatthew G Knepley 
23a312c225SMatthew G Knepley   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(SNESReset_NGMRES(snes));
259566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ngmres->h, ngmres->beta, ngmres->xi, ngmres->q));
269566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ngmres->xnorms, ngmres->fnorms, ngmres->s));
27dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX)
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(ngmres->rwork));
2919653cdaSPeter Brune #endif
309566063dSJacob Faibussowitsch   PetscCall(PetscFree(ngmres->work));
312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", NULL));
322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", NULL));
332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", NULL));
342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", NULL));
359566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37a312c225SMatthew G Knepley }
38a312c225SMatthew G Knepley 
39d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NGMRES(SNES snes)
40d71ae5a4SJacob Faibussowitsch {
41a312c225SMatthew G Knepley   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
42e7058c64SPeter Brune   const char  *optionsprefix;
4319653cdaSPeter Brune   PetscInt     msize, hsize;
44fced5a79SAsbjørn Nilsen Riseth   DM           dm;
45a312c225SMatthew G Knepley 
46a312c225SMatthew G Knepley   PetscFunctionBegin;
47efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
4846159c86SPeter Brune     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNGMRES does not support left preconditioning with unpreconditioned function");
4946159c86SPeter Brune   }
50efd4aadfSBarry Smith   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
519566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 5));
52fced5a79SAsbjørn Nilsen Riseth 
53fced5a79SAsbjørn Nilsen Riseth   if (!snes->vec_sol) {
549566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
559566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
56fced5a79SAsbjørn Nilsen Riseth   }
57fced5a79SAsbjørn Nilsen Riseth 
589566063dSJacob Faibussowitsch   if (!ngmres->Xdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Xdot));
599566063dSJacob Faibussowitsch   if (!ngmres->Fdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Fdot));
6078440776SJed Brown   if (!ngmres->setup_called) {
61087dfb9eSxuemin     msize = ngmres->msize; /* restart size */
6219653cdaSPeter Brune     hsize = msize * msize;
63087dfb9eSxuemin 
6498b3e84cSPeter Brune     /* explicit least squares minimization solve */
659566063dSJacob Faibussowitsch     PetscCall(PetscCalloc4(hsize, &ngmres->h, msize, &ngmres->beta, msize, &ngmres->xi, hsize, &ngmres->q));
669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(msize, &ngmres->xnorms, msize, &ngmres->fnorms, msize, &ngmres->s));
6719653cdaSPeter Brune     ngmres->nrhs  = 1;
6819653cdaSPeter Brune     ngmres->lda   = msize;
6919653cdaSPeter Brune     ngmres->ldb   = msize;
7019653cdaSPeter Brune     ngmres->lwork = 12 * msize;
71dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX)
729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->rwork));
7319653cdaSPeter Brune #endif
749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->work));
7578440776SJed Brown   }
76e7058c64SPeter Brune 
77e7058c64SPeter Brune   /* linesearch setup */
789566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
79e7058c64SPeter Brune 
8013a62661SPeter Brune   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
819566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &ngmres->additive_linesearch));
829566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetSNES(ngmres->additive_linesearch, snes));
8348a46eb9SPierre Jolivet     if (!((PetscObject)ngmres->additive_linesearch)->type_name) PetscCall(SNESLineSearchSetType(ngmres->additive_linesearch, SNESLINESEARCHL2));
849566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, "additive_"));
859566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, optionsprefix));
869566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch));
87e7058c64SPeter Brune   }
88e7058c64SPeter Brune 
8978440776SJed Brown   ngmres->setup_called = PETSC_TRUE;
903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91a312c225SMatthew G Knepley }
92a312c225SMatthew G Knepley 
9366976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes, PetscOptionItems *PetscOptionsObject)
94d71ae5a4SJacob Faibussowitsch {
95a312c225SMatthew G Knepley   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
9694ae4db5SBarry Smith   PetscBool    debug  = PETSC_FALSE;
970adebc6cSBarry Smith 
98a312c225SMatthew G Knepley   PetscFunctionBegin;
99d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NGMRES options");
100d0609cedSBarry Smith   PetscCall(PetscOptionsEnum("-snes_ngmres_select_type", "Select type", "SNESNGMRESSetSelectType", SNESNGMRESSelectTypes, (PetscEnum)ngmres->select_type, (PetscEnum *)&ngmres->select_type, NULL));
101d0609cedSBarry Smith   PetscCall(PetscOptionsEnum("-snes_ngmres_restart_type", "Restart type", "SNESNGMRESSetRestartType", SNESNGMRESRestartTypes, (PetscEnum)ngmres->restart_type, (PetscEnum *)&ngmres->restart_type, NULL));
1029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES", ngmres->candidate, &ngmres->candidate, NULL));
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc", "Linearly approximate the function", "SNES", ngmres->approxfunc, &ngmres->approxfunc, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES", ngmres->restart_periodic, &ngmres->restart_periodic, NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart", "SNES", ngmres->restart_it, &ngmres->restart_it, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
108ad540459SPierre Jolivet   if (debug) ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
1099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, NULL));
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, NULL));
1119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, NULL));
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, NULL));
1139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES", ngmres->singlereduction, &ngmres->singlereduction, NULL));
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise", ngmres->restart_fm_rise, &ngmres->restart_fm_rise, NULL));
115d0609cedSBarry Smith   PetscOptionsHeadEnd();
1166a7cf640SPeter Brune   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118a312c225SMatthew G Knepley }
119a312c225SMatthew G Knepley 
120d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
121d71ae5a4SJacob Faibussowitsch {
122a312c225SMatthew G Knepley   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
123a312c225SMatthew G Knepley   PetscBool    iascii;
124a312c225SMatthew G Knepley 
125a312c225SMatthew G Knepley   PetscFunctionBegin;
1269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
127a312c225SMatthew G Knepley   if (iascii) {
12863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of stored past updates: %" PetscInt_FMT "\n", ngmres->msize));
12963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", (double)ngmres->gammaA, (double)ngmres->gammaC));
13063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", (double)ngmres->epsilonB, (double)ngmres->deltaB));
13163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Restart on F_M residual increase: %s\n", PetscBools[ngmres->restart_fm_rise]));
132a312c225SMatthew G Knepley   }
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
134a312c225SMatthew G Knepley }
135a312c225SMatthew G Knepley 
13666976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NGMRES(SNES snes)
137d71ae5a4SJacob Faibussowitsch {
138087dfb9eSxuemin   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
13998b3e84cSPeter Brune   /* present solution, residual, and preconditioned residual */
1409f425c49SPeter Brune   Vec X, F, B, D, Y;
141f109b39eSPeter Brune 
142f109b39eSPeter Brune   /* candidate linear combination answers */
143ddd40ce5SPeter Brune   Vec XA, FA, XM, FM;
14419653cdaSPeter Brune 
14598b3e84cSPeter Brune   /* coefficients and RHS to the minimization problem */
14618aa0c0cSPeter Brune   PetscReal fnorm, fMnorm, fAnorm;
147b3c6a99cSPeter Brune   PetscReal xnorm, xMnorm, xAnorm;
148b3c6a99cSPeter Brune   PetscReal ynorm, yMnorm, yAnorm;
14938774f0aSPeter Brune   PetscInt  k, k_restart, l, ivec, restart_count = 0;
15019653cdaSPeter Brune 
15198b3e84cSPeter Brune   /* solution selection data */
15238774f0aSPeter Brune   PetscBool selectRestart;
15361ba4676SBarry Smith   /*
15461ba4676SBarry Smith       These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed
15561ba4676SBarry Smith       to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case
15661ba4676SBarry Smith       so the code is correct as written.
15761ba4676SBarry Smith   */
15861ba4676SBarry Smith   PetscReal dnorm = 0.0, dminnorm = 0.0;
159b3c6a99cSPeter Brune   PetscReal fminnorm;
16019653cdaSPeter Brune 
1611e633543SBarry Smith   SNESConvergedReason  reason;
162422a814eSBarry Smith   SNESLineSearchReason lssucceed;
163a312c225SMatthew G Knepley 
164a312c225SMatthew G Knepley   PetscFunctionBegin;
1650b121fc5SBarry Smith   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);
166c579b300SPatrick Farrell 
1679566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
16898b3e84cSPeter Brune   /* variable initialization */
169a312c225SMatthew G Knepley   snes->reason = SNES_CONVERGED_ITERATING;
170f109b39eSPeter Brune   X            = snes->vec_sol;
171f109b39eSPeter Brune   F            = snes->vec_func;
172f109b39eSPeter Brune   B            = snes->vec_rhs;
173f109b39eSPeter Brune   XA           = snes->vec_sol_update;
174f109b39eSPeter Brune   FA           = snes->work[0];
175f109b39eSPeter Brune   D            = snes->work[1];
176f109b39eSPeter Brune 
177f109b39eSPeter Brune   /* work for the line search */
178f109b39eSPeter Brune   Y  = snes->work[2];
1799f425c49SPeter Brune   XM = snes->work[3];
1809f425c49SPeter Brune   FM = snes->work[4];
181a312c225SMatthew G Knepley 
1829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
183a312c225SMatthew G Knepley   snes->iter = 0;
184a312c225SMatthew G Knepley   snes->norm = 0.;
1859566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
18619653cdaSPeter Brune 
18798b3e84cSPeter Brune   /* initialization */
18819653cdaSPeter Brune 
189efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT) {
1909566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1919566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
1923a2ae377SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1933a2ae377SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1943ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
1953a2ae377SPeter Brune     }
1969566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
1973a2ae377SPeter Brune   } else {
198e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1999566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
2001aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
201c1c75074SPeter Brune 
2029566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
203422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2043a2ae377SPeter Brune   }
205e4ed7901SPeter Brune   fminnorm = fnorm;
20619653cdaSPeter Brune 
2079566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
208f109b39eSPeter Brune   snes->norm = fnorm;
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2109566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
2112d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
2129566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
2133ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
2143ba16761SJacob Faibussowitsch   PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, F, fnorm, X));
215a312c225SMatthew G Knepley 
21619653cdaSPeter Brune   k_restart = 1;
21719653cdaSPeter Brune   l         = 1;
218b3c6a99cSPeter Brune   ivec      = 0;
21909c08436SPeter Brune   for (k = 1; k < snes->max_its + 1; k++) {
220d52c4f7dSRené Chenard     /* Call general purpose update function */
221d52c4f7dSRené Chenard     PetscTryTypeMethod(snes, update, snes->iter);
222d52c4f7dSRené Chenard 
22398b3e84cSPeter Brune     /* Computation of x^M */
224efd4aadfSBarry Smith     if (snes->npc && snes->npcside == PC_RIGHT) {
2259566063dSJacob Faibussowitsch       PetscCall(VecCopy(X, XM));
2269566063dSJacob Faibussowitsch       PetscCall(SNESSetInitialFunction(snes->npc, F));
22763e7833aSPeter Brune 
2289566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0));
2299566063dSJacob Faibussowitsch       PetscCall(SNESSolve(snes->npc, B, XM));
2309566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0));
23163e7833aSPeter Brune 
2329566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
2338cc86e31SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2348cc86e31SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
2353ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
2368cc86e31SPeter Brune       }
2379566063dSJacob Faibussowitsch       PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm));
2388cc86e31SPeter Brune     } else {
239f109b39eSPeter Brune       /* no preconditioner -- just take gradient descent with line search */
2409566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, Y));
2419566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, FM));
2429566063dSJacob Faibussowitsch       PetscCall(VecCopy(X, XM));
2431aa26658SKarl Rupp 
244e7058c64SPeter Brune       fMnorm = fnorm;
2451aa26658SKarl Rupp 
2469566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchApply(snes->linesearch, XM, FM, &fMnorm, Y));
2479566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
248422a814eSBarry Smith       if (lssucceed) {
249f109b39eSPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
250f109b39eSPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
2513ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
252f109b39eSPeter Brune         }
253f109b39eSPeter Brune       }
2546634f59bSPeter Brune     }
25523b3e82cSAsbjørn Nilsen Riseth 
2569566063dSJacob Faibussowitsch     PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA));
25798b3e84cSPeter Brune     /* r = F(x) */
2589f425c49SPeter Brune     if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */
25919653cdaSPeter Brune 
2609f425c49SPeter Brune     /* differences for selection and restart */
26113a62661SPeter Brune     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
2629566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
26313a62661SPeter Brune     } else {
2649566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
26513a62661SPeter Brune     }
266422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
2671aa26658SKarl Rupp 
2689f425c49SPeter Brune     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
2699566063dSJacob Faibussowitsch     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));
27019653cdaSPeter Brune     selectRestart = PETSC_FALSE;
27123b3e82cSAsbjørn Nilsen Riseth 
27213a62661SPeter Brune     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
2739566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESSelectRestart_Private(snes, l, fMnorm, fAnorm, dnorm, fminnorm, dminnorm, &selectRestart));
27423b3e82cSAsbjørn Nilsen Riseth 
27528ed4a04SPeter Brune       /* if the restart conditions persist for more than restart_it iterations, restart. */
2761aa26658SKarl Rupp       if (selectRestart) restart_count++;
2771aa26658SKarl Rupp       else restart_count = 0;
27813a62661SPeter Brune     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
27913a62661SPeter Brune       if (k_restart > ngmres->restart_periodic) {
28063a3b9bcSJacob Faibussowitsch         if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart));
28113a62661SPeter Brune         restart_count = ngmres->restart_it;
28213a62661SPeter Brune       }
28313a62661SPeter Brune     }
28423b3e82cSAsbjørn Nilsen Riseth 
285b3c6a99cSPeter Brune     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
28623b3e82cSAsbjørn Nilsen Riseth 
28728ed4a04SPeter Brune     /* restart after restart conditions have persisted for a fixed number of iterations */
28828ed4a04SPeter Brune     if (restart_count >= ngmres->restart_it) {
28948a46eb9SPierre Jolivet       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart));
29028ed4a04SPeter Brune       restart_count = 0;
29119653cdaSPeter Brune       k_restart     = 1;
29219653cdaSPeter Brune       l             = 1;
293b3c6a99cSPeter Brune       ivec          = 0;
29498b3e84cSPeter Brune       /* q_{00} = nu */
2959566063dSJacob Faibussowitsch       PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, FM, fMnorm, XM));
296d2e16ddcSPeter Brune     } else {
29798b3e84cSPeter Brune       /* select the current size of the subspace */
2981e633543SBarry Smith       if (l < ngmres->msize) l++;
29919653cdaSPeter Brune       k_restart++;
30098b3e84cSPeter Brune       /* place the current entry in the list of previous entries */
30138774f0aSPeter Brune       if (ngmres->candidate) {
30238774f0aSPeter Brune         if (fminnorm > fMnorm) fminnorm = fMnorm;
3039566063dSJacob Faibussowitsch         PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fMnorm, XM));
304d2e16ddcSPeter Brune       } else {
30538774f0aSPeter Brune         if (fminnorm > fnorm) fminnorm = fnorm;
3069566063dSJacob Faibussowitsch         PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, F, fnorm, X));
30719653cdaSPeter Brune       }
308d2e16ddcSPeter Brune     }
30919653cdaSPeter Brune 
3109566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
311087dfb9eSxuemin     snes->iter  = k;
312f109b39eSPeter Brune     snes->norm  = fnorm;
31338606ef4SRené Chenard     snes->ynorm = ynorm;
31438606ef4SRené Chenard     snes->xnorm = xnorm;
3159566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3169566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
3172d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, 0, 0, fnorm));
3189566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
3193ba16761SJacob Faibussowitsch     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
320a312c225SMatthew G Knepley   }
321a312c225SMatthew G Knepley   snes->reason = SNES_DIVERGED_MAX_IT;
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323a312c225SMatthew G Knepley }
324a312c225SMatthew G Knepley 
32523b3e82cSAsbjørn Nilsen Riseth /*@
32623b3e82cSAsbjørn Nilsen Riseth   SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M
32723b3e82cSAsbjørn Nilsen Riseth 
32823b3e82cSAsbjørn Nilsen Riseth   Input Parameters:
329f6dfbefdSBarry Smith + snes - the `SNES` context.
330f6dfbefdSBarry Smith - flg  - boolean value deciding whether to use the option or not, default is `PETSC_FALSE`
33123b3e82cSAsbjørn Nilsen Riseth 
332f6dfbefdSBarry Smith   Options Database Key:
333f6dfbefdSBarry Smith . -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M
33423b3e82cSAsbjørn Nilsen Riseth 
33523b3e82cSAsbjørn Nilsen Riseth   Level: intermediate
33623b3e82cSAsbjørn Nilsen Riseth 
33723b3e82cSAsbjørn Nilsen Riseth   Notes:
33823b3e82cSAsbjørn Nilsen Riseth   If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area.
33923b3e82cSAsbjørn Nilsen Riseth   To help the solver do that, reset the Krylov subspace whenever F_M increases.
34023b3e82cSAsbjørn Nilsen Riseth 
341f6dfbefdSBarry Smith   This option must be used with the `SNESNGMRES` `SNESNGMRESRestartType` of `SNES_NGMRES_RESTART_DIFFERENCE`
34223b3e82cSAsbjørn Nilsen Riseth 
343420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartType()`
34423b3e82cSAsbjørn Nilsen Riseth   @*/
345d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes, PetscBool flg)
346d71ae5a4SJacob Faibussowitsch {
34723b3e82cSAsbjørn Nilsen Riseth   PetscErrorCode (*f)(SNES, PetscBool);
34823b3e82cSAsbjørn Nilsen Riseth 
34923b3e82cSAsbjørn Nilsen Riseth   PetscFunctionBegin;
3509566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", &f));
3519566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, flg));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35323b3e82cSAsbjørn Nilsen Riseth }
35423b3e82cSAsbjørn Nilsen Riseth 
35566976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes, PetscBool flg)
356d71ae5a4SJacob Faibussowitsch {
35723b3e82cSAsbjørn Nilsen Riseth   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
35823b3e82cSAsbjørn Nilsen Riseth 
35923b3e82cSAsbjørn Nilsen Riseth   PetscFunctionBegin;
36023b3e82cSAsbjørn Nilsen Riseth   ngmres->restart_fm_rise = flg;
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36223b3e82cSAsbjørn Nilsen Riseth }
36323b3e82cSAsbjørn Nilsen Riseth 
364d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes, PetscBool *flg)
365d71ae5a4SJacob Faibussowitsch {
36623b3e82cSAsbjørn Nilsen Riseth   PetscErrorCode (*f)(SNES, PetscBool *);
36723b3e82cSAsbjørn Nilsen Riseth 
36823b3e82cSAsbjørn Nilsen Riseth   PetscFunctionBegin;
3699566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", &f));
3709566063dSJacob Faibussowitsch   if (f) PetscCall((f)(snes, flg));
3713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37223b3e82cSAsbjørn Nilsen Riseth }
37323b3e82cSAsbjørn Nilsen Riseth 
37466976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes, PetscBool *flg)
375d71ae5a4SJacob Faibussowitsch {
37623b3e82cSAsbjørn Nilsen Riseth   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
37723b3e82cSAsbjørn Nilsen Riseth 
37823b3e82cSAsbjørn Nilsen Riseth   PetscFunctionBegin;
37923b3e82cSAsbjørn Nilsen Riseth   *flg = ngmres->restart_fm_rise;
3803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38123b3e82cSAsbjørn Nilsen Riseth }
38223b3e82cSAsbjørn Nilsen Riseth 
38313a62661SPeter Brune /*@
384f6dfbefdSBarry Smith   SNESNGMRESSetRestartType - Sets the restart type for `SNESNGMRES`.
38513a62661SPeter Brune 
386c3339decSBarry Smith   Logically Collective
38713a62661SPeter Brune 
38813a62661SPeter Brune   Input Parameters:
38913a62661SPeter Brune + snes  - the iterative context
390ceaaa498SBarry Smith - rtype - restart type, see `SNESNGMRESRestartType`
39113a62661SPeter Brune 
392da81f932SPierre Jolivet   Options Database Keys:
39313a62661SPeter Brune + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
394420bcc1bSBarry Smith - -snes_ngmres_restart <30>                           - sets the number of iterations before restart for periodic
39513a62661SPeter Brune 
39613a62661SPeter Brune   Level: intermediate
39713a62661SPeter Brune 
398420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartFmRise()`,
399420bcc1bSBarry Smith           `SNESNGMRESSetSelectType()`
40013a62661SPeter Brune @*/
401d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype)
402d71ae5a4SJacob Faibussowitsch {
40313a62661SPeter Brune   PetscFunctionBegin;
40413a62661SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
405cac4c232SBarry Smith   PetscTryMethod(snes, "SNESNGMRESSetRestartType_C", (SNES, SNESNGMRESRestartType), (snes, rtype));
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40713a62661SPeter Brune }
40813a62661SPeter Brune 
40913a62661SPeter Brune /*@
410f6dfbefdSBarry Smith   SNESNGMRESSetSelectType - Sets the selection type for `SNESNGMRES`.  This determines how the candidate solution and
41113a62661SPeter Brune   combined solution are used to create the next iterate.
41213a62661SPeter Brune 
413c3339decSBarry Smith   Logically Collective
41413a62661SPeter Brune 
41513a62661SPeter Brune   Input Parameters:
41613a62661SPeter Brune + snes  - the iterative context
417ceaaa498SBarry Smith - stype - selection type, see `SNESNGMRESSelectType`
41813a62661SPeter Brune 
419f6dfbefdSBarry Smith   Options Database Key:
42067b8a455SSatish Balay . -snes_ngmres_select_type<difference,none,linesearch> - select type
42113a62661SPeter Brune 
42213a62661SPeter Brune   Level: intermediate
42313a62661SPeter Brune 
424f6dfbefdSBarry Smith   Note:
425f6dfbefdSBarry Smith   The default line search used is the `SNESLINESEARCHL2` line search and it requires two additional function evaluations.
42613a62661SPeter Brune 
427420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNGMRES`, `SNESNGMRESSelectType`, `SNES_NGMRES_SELECT_NONE`, `SNES_NGMRES_SELECT_DIFFERENCE`, `SNES_NGMRES_SELECT_LINESEARCH`,
428420bcc1bSBarry Smith           `SNESNGMRESSetRestartType()`
42913a62661SPeter Brune @*/
430d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetSelectType(SNES snes, SNESNGMRESSelectType stype)
431d71ae5a4SJacob Faibussowitsch {
43213a62661SPeter Brune   PetscFunctionBegin;
43313a62661SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
434cac4c232SBarry Smith   PetscTryMethod(snes, "SNESNGMRESSetSelectType_C", (SNES, SNESNGMRESSelectType), (snes, stype));
4353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43613a62661SPeter Brune }
43713a62661SPeter Brune 
43866976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype)
439d71ae5a4SJacob Faibussowitsch {
44013a62661SPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
4415fd66863SKarl Rupp 
44213a62661SPeter Brune   PetscFunctionBegin;
44313a62661SPeter Brune   ngmres->select_type = stype;
4443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44513a62661SPeter Brune }
44613a62661SPeter Brune 
44766976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype)
448d71ae5a4SJacob Faibussowitsch {
44913a62661SPeter Brune   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
4505fd66863SKarl Rupp 
45113a62661SPeter Brune   PetscFunctionBegin;
45213a62661SPeter Brune   ngmres->restart_type = rtype;
4533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45413a62661SPeter Brune }
45513a62661SPeter Brune 
456dfbf837cSBarry Smith /*MC
457*1d27aa22SBarry Smith   SNESNGMRES - The Nonlinear Generalized Minimum Residual method {cite}`ow1`, {cite}`bruneknepleysmithtu15`
458a312c225SMatthew G Knepley 
459dfbf837cSBarry Smith    Level: beginner
460dfbf837cSBarry Smith 
461f6dfbefdSBarry Smith    Options Database Keys:
46213a62661SPeter Brune +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
46338774f0aSPeter Brune .  -snes_ngmres_restart_type<difference,none,periodic>  - choose the restart conditions
464f6dfbefdSBarry Smith .  -snes_ngmres_candidate                               - Use `SNESNGMRES` variant which combines candidate solutions instead of actual solutions
46513a62661SPeter Brune .  -snes_ngmres_m                                       - Number of stored previous solutions and residuals
46613a62661SPeter Brune .  -snes_ngmres_restart_it                              - Number of iterations the restart conditions hold before restart
46713a62661SPeter Brune .  -snes_ngmres_gammaA                                  - Residual tolerance for solution select between the candidate and combination
46813a62661SPeter Brune .  -snes_ngmres_gammaC                                  - Residual tolerance for restart
46913a62661SPeter Brune .  -snes_ngmres_epsilonB                                - Difference tolerance between subsequent solutions triggering restart
47013a62661SPeter Brune .  -snes_ngmres_deltaB                                  - Difference tolerance between residuals triggering restart
47123b3e82cSAsbjørn Nilsen Riseth .  -snes_ngmres_restart_fm_rise                         - Restart on residual rise from x_M step
47213a62661SPeter Brune .  -snes_ngmres_monitor                                 - Prints relevant information about the ngmres iteration
4735c3e6ab7SPeter Brune .  -snes_linesearch_type <basic,l2,cp>                  - Line search type used for the default smoother
47413a62661SPeter Brune -  -additive_snes_linesearch_type                       - line search type used to select between the candidate and combined solution with additive select type
4751867fe5bSPeter Brune 
4761867fe5bSPeter Brune    Notes:
4771867fe5bSPeter Brune    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
4781867fe5bSPeter Brune    optimization problem at each iteration.
4791867fe5bSPeter Brune 
480f6dfbefdSBarry Smith    Very similar to the `SNESANDERSON` algorithm.
4814f02bc6aSBarry Smith 
482420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESANDERSON`, `SNESNGMRESSetSelectType()`, `SNESNGMRESSetRestartType()`,
483ceaaa498SBarry Smith           `SNESNGMRESSetRestartFmRise()`, `SNESNGMRESSelectType`, ``SNESNGMRESRestartType`
484dfbf837cSBarry Smith M*/
485a312c225SMatthew G Knepley 
486d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
487d71ae5a4SJacob Faibussowitsch {
488a312c225SMatthew G Knepley   SNES_NGMRES   *ngmres;
489d8d34be6SBarry Smith   SNESLineSearch linesearch;
490a312c225SMatthew G Knepley 
491a312c225SMatthew G Knepley   PetscFunctionBegin;
492a312c225SMatthew G Knepley   snes->ops->destroy        = SNESDestroy_NGMRES;
493a312c225SMatthew G Knepley   snes->ops->setup          = SNESSetUp_NGMRES;
494a312c225SMatthew G Knepley   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
495a312c225SMatthew G Knepley   snes->ops->view           = SNESView_NGMRES;
496a312c225SMatthew G Knepley   snes->ops->solve          = SNESSolve_NGMRES;
497a312c225SMatthew G Knepley   snes->ops->reset          = SNESReset_NGMRES;
498a312c225SMatthew G Knepley 
499efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
5002c155ee1SBarry Smith   snes->usesksp = PETSC_FALSE;
501efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
5022c155ee1SBarry Smith 
5034fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
5044fc747eaSLawrence Mitchell 
5054dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ngmres));
506a312c225SMatthew G Knepley   snes->data    = (void *)ngmres;
507d2e16ddcSPeter Brune   ngmres->msize = 30;
50819653cdaSPeter Brune 
50988976e71SPeter Brune   if (!snes->tolerancesset) {
5100e444f03SPeter Brune     snes->max_funcs = 30000;
5110e444f03SPeter Brune     snes->max_its   = 10000;
51288976e71SPeter Brune   }
5130e444f03SPeter Brune 
51438774f0aSPeter Brune   ngmres->candidate = PETSC_FALSE;
515d2e16ddcSPeter Brune 
5169566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
51748a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
518d8d34be6SBarry Smith 
5190298fd71SBarry Smith   ngmres->additive_linesearch = NULL;
520077c4231SPeter Brune   ngmres->approxfunc          = PETSC_FALSE;
52128ed4a04SPeter Brune   ngmres->restart_it          = 2;
52213a62661SPeter Brune   ngmres->restart_periodic    = 30;
523f109b39eSPeter Brune   ngmres->gammaA              = 2.0;
524f109b39eSPeter Brune   ngmres->gammaC              = 2.0;
525cac108bcSPeter Brune   ngmres->deltaB              = 0.9;
526cac108bcSPeter Brune   ngmres->epsilonB            = 0.1;
52723b3e82cSAsbjørn Nilsen Riseth   ngmres->restart_fm_rise     = PETSC_FALSE;
528e7058c64SPeter Brune 
52913a62661SPeter Brune   ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
53013a62661SPeter Brune   ngmres->select_type  = SNES_NGMRES_SELECT_DIFFERENCE;
53113a62661SPeter Brune 
5329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", SNESNGMRESSetSelectType_NGMRES));
5339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", SNESNGMRESSetRestartType_NGMRES));
5349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", SNESNGMRESSetRestartFmRise_NGMRES));
5359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", SNESNGMRESGetRestartFmRise_NGMRES));
5363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
537a312c225SMatthew G Knepley }
538