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
SNESReset_NGMRES(SNES snes)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
SNESDestroy_NGMRES(SNES snes)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
SNESSetUp_NGMRES(SNES snes)39d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NGMRES(SNES snes)
40d71ae5a4SJacob Faibussowitsch {
41a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
4219653cdaSPeter Brune PetscInt msize, hsize;
43fced5a79SAsbjørn Nilsen Riseth DM dm;
44a312c225SMatthew G Knepley
45a312c225SMatthew G Knepley PetscFunctionBegin;
46966bd95aSPierre Jolivet PetscCheck(!snes->npc || snes->npcside != PC_LEFT || snes->functype != SNES_FUNCTION_UNPRECONDITIONED, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNGMRES does not support left preconditioning with unpreconditioned function");
47efd4aadfSBarry Smith if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
489566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 5));
49fced5a79SAsbjørn Nilsen Riseth
50fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) {
519566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
529566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
53fced5a79SAsbjørn Nilsen Riseth }
54fced5a79SAsbjørn Nilsen Riseth
559566063dSJacob Faibussowitsch if (!ngmres->Xdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Xdot));
569566063dSJacob Faibussowitsch if (!ngmres->Fdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Fdot));
5778440776SJed Brown if (!ngmres->setup_called) {
58087dfb9eSxuemin msize = ngmres->msize; /* restart size */
5919653cdaSPeter Brune hsize = msize * msize;
60087dfb9eSxuemin
6198b3e84cSPeter Brune /* explicit least squares minimization solve */
629566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(hsize, &ngmres->h, msize, &ngmres->beta, msize, &ngmres->xi, hsize, &ngmres->q));
639566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(msize, &ngmres->xnorms, msize, &ngmres->fnorms, msize, &ngmres->s));
6419653cdaSPeter Brune ngmres->nrhs = 1;
656497c311SBarry Smith PetscCall(PetscBLASIntCast(msize, &ngmres->lda));
666497c311SBarry Smith PetscCall(PetscBLASIntCast(msize, &ngmres->ldb));
676497c311SBarry Smith PetscCall(PetscBLASIntCast(12 * msize, &ngmres->lwork));
68dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX)
699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->rwork));
7019653cdaSPeter Brune #endif
719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->work));
7278440776SJed Brown }
73e7058c64SPeter Brune
7478440776SJed Brown ngmres->setup_called = PETSC_TRUE;
753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
76a312c225SMatthew G Knepley }
77a312c225SMatthew G Knepley
SNESSetFromOptions_NGMRES(SNES snes,PetscOptionItems PetscOptionsObject)78ce78bad3SBarry Smith static PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes, PetscOptionItems PetscOptionsObject)
79d71ae5a4SJacob Faibussowitsch {
80a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
8194ae4db5SBarry Smith PetscBool debug = PETSC_FALSE;
820adebc6cSBarry Smith
83a312c225SMatthew G Knepley PetscFunctionBegin;
84d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES NGMRES options");
85d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-snes_ngmres_select_type", "Select type", "SNESNGMRESSetSelectType", SNESNGMRESSelectTypes, (PetscEnum)ngmres->select_type, (PetscEnum *)&ngmres->select_type, NULL));
86d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-snes_ngmres_restart_type", "Restart type", "SNESNGMRESSetRestartType", SNESNGMRESRestartTypes, (PetscEnum)ngmres->restart_type, (PetscEnum *)&ngmres->restart_type, NULL));
879566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES", ngmres->candidate, &ngmres->candidate, NULL));
889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc", "Linearly approximate the function", "SNES", ngmres->approxfunc, &ngmres->approxfunc, NULL));
899566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, NULL));
909566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES", ngmres->restart_periodic, &ngmres->restart_periodic, NULL));
919566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart", "SNES", ngmres->restart_it, &ngmres->restart_it, NULL));
929566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
93ad540459SPierre Jolivet if (debug) ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, NULL));
959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, NULL));
969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, NULL));
979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, NULL));
989566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise", ngmres->restart_fm_rise, &ngmres->restart_fm_rise, NULL));
99d0609cedSBarry Smith PetscOptionsHeadEnd();
1004936d809SStefano Zampini if (ngmres->gammaA > ngmres->gammaC && ngmres->gammaC > 2.) ngmres->gammaC = ngmres->gammaA;
1014936d809SStefano Zampini if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
1024936d809SStefano Zampini PetscCall(SNESNGMRESGetAdditiveLineSearch_Private(snes, &ngmres->additive_linesearch));
1034936d809SStefano Zampini PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch));
1044936d809SStefano Zampini }
1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
106a312c225SMatthew G Knepley }
107a312c225SMatthew G Knepley
SNESView_NGMRES(SNES snes,PetscViewer viewer)108d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
109d71ae5a4SJacob Faibussowitsch {
110a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
1119f196a02SMartin Diehl PetscBool isascii;
112a312c225SMatthew G Knepley
113a312c225SMatthew G Knepley PetscFunctionBegin;
1149f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1159f196a02SMartin Diehl if (isascii) {
11663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of stored past updates: %" PetscInt_FMT "\n", ngmres->msize));
1174936d809SStefano Zampini if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
11863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", (double)ngmres->gammaA, (double)ngmres->gammaC));
11963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", (double)ngmres->epsilonB, (double)ngmres->deltaB));
12063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Restart on F_M residual increase: %s\n", PetscBools[ngmres->restart_fm_rise]));
121a312c225SMatthew G Knepley }
1224936d809SStefano Zampini if (ngmres->additive_linesearch) {
1234936d809SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Additive line-search details:\n"));
1244936d809SStefano Zampini PetscCall(SNESLineSearchView(ngmres->additive_linesearch, viewer));
1254936d809SStefano Zampini }
1264936d809SStefano Zampini }
1273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
128a312c225SMatthew G Knepley }
129a312c225SMatthew G Knepley
SNESSolve_NGMRES(SNES snes)13066976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NGMRES(SNES snes)
131d71ae5a4SJacob Faibussowitsch {
132087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
1336b72add0SBarry Smith Vec X, F, B, D, Y; /* present solution, residual, and preconditioned residual */
1346b72add0SBarry Smith Vec XA, FA, XM, FM; /* candidate linear combination answers */
1356b72add0SBarry Smith PetscReal fnorm, fMnorm, fAnorm; /* coefficients and RHS to the minimization problem */
136b3c6a99cSPeter Brune PetscReal xnorm, xMnorm, xAnorm;
137b3c6a99cSPeter Brune PetscReal ynorm, yMnorm, yAnorm;
13838774f0aSPeter Brune PetscInt k, k_restart, l, ivec, restart_count = 0;
1396b72add0SBarry Smith PetscReal objmin, objM, objA, obj; /* support for objective functions minimization */
1406b72add0SBarry Smith PetscBool selectRestart; /* solution selection data */
1416b72add0SBarry Smith SNESConvergedReason reason;
1426b72add0SBarry Smith SNESObjectiveFn *objective;
14361ba4676SBarry Smith /*
14461ba4676SBarry Smith These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed
14561ba4676SBarry Smith to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case
14661ba4676SBarry Smith so the code is correct as written.
14761ba4676SBarry Smith */
14861ba4676SBarry Smith PetscReal dnorm = 0.0, dminnorm = 0.0;
14919653cdaSPeter Brune
150a312c225SMatthew G Knepley PetscFunctionBegin;
1510b121fc5SBarry 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);
152c579b300SPatrick Farrell
1539566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
15498b3e84cSPeter Brune /* variable initialization */
155a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING;
156f109b39eSPeter Brune X = snes->vec_sol;
157f109b39eSPeter Brune F = snes->vec_func;
158f109b39eSPeter Brune B = snes->vec_rhs;
1595ef867c2SRené Chenard XA = snes->work[2];
160f109b39eSPeter Brune FA = snes->work[0];
161f109b39eSPeter Brune D = snes->work[1];
162f109b39eSPeter Brune
163f109b39eSPeter Brune /* work for the line search */
1645ef867c2SRené Chenard Y = snes->vec_sol_update;
1659f425c49SPeter Brune XM = snes->work[3];
1669f425c49SPeter Brune FM = snes->work[4];
167a312c225SMatthew G Knepley
1689566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
169a312c225SMatthew G Knepley snes->iter = 0;
170a312c225SMatthew G Knepley snes->norm = 0.;
1719566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
17219653cdaSPeter Brune
17398b3e84cSPeter Brune /* initialization */
17419653cdaSPeter Brune
175efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT) {
1769566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F));
1779566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
1783a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1793a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER;
1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1813a2ae377SPeter Brune }
1823a2ae377SPeter Brune } else {
183ac530a7eSPierre Jolivet if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
184ac530a7eSPierre Jolivet else snes->vec_func_init_set = PETSC_FALSE;
1853a2ae377SPeter Brune }
186*76c63389SBarry Smith PetscCall(VecNorm(F, NORM_2, &fnorm));
187*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
1884936d809SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL));
1894936d809SStefano Zampini objmin = fnorm;
190*76c63389SBarry Smith if (objective) {
191*76c63389SBarry Smith PetscCall(SNESComputeObjective(snes, X, &objmin));
192*76c63389SBarry Smith SNESCheckObjectiveDomainError(snes, objmin);
193*76c63389SBarry Smith }
1944936d809SStefano Zampini obj = objmin;
19519653cdaSPeter Brune
1969566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
197f109b39eSPeter Brune snes->norm = fnorm;
1989566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1999566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
2002d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
2019566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm));
2023ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
2033ba16761SJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, F, fnorm, X));
204a312c225SMatthew G Knepley
20519653cdaSPeter Brune k_restart = 1;
20619653cdaSPeter Brune l = 1;
207b3c6a99cSPeter Brune ivec = 0;
20809c08436SPeter Brune for (k = 1; k < snes->max_its + 1; k++) {
209d52c4f7dSRené Chenard /* Call general purpose update function */
210d52c4f7dSRené Chenard PetscTryTypeMethod(snes, update, snes->iter);
211d52c4f7dSRené Chenard
21298b3e84cSPeter Brune /* Computation of x^M */
213efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_RIGHT) {
2149566063dSJacob Faibussowitsch PetscCall(VecCopy(X, XM));
2159566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc, F));
21663e7833aSPeter Brune
2179566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0));
2189566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, B, XM));
2199566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0));
22063e7833aSPeter Brune
2219566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason));
2228cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
2238cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER;
2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2258cc86e31SPeter Brune }
2269566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm));
2278cc86e31SPeter Brune } else {
228f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */
2299566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Y));
2309566063dSJacob Faibussowitsch PetscCall(VecCopy(F, FM));
2319566063dSJacob Faibussowitsch PetscCall(VecCopy(X, XM));
2321aa26658SKarl Rupp
233e7058c64SPeter Brune fMnorm = fnorm;
2341aa26658SKarl Rupp
2359566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, XM, FM, &fMnorm, Y));
236*76c63389SBarry Smith if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
237f109b39eSPeter Brune }
238*76c63389SBarry Smith if (objective) {
239*76c63389SBarry Smith PetscCall(SNESComputeObjective(snes, XM, &objM));
240*76c63389SBarry Smith SNESCheckObjectiveDomainError(snes, objM);
241*76c63389SBarry Smith } else objM = fMnorm;
2424936d809SStefano Zampini objmin = PetscMin(objmin, objM);
24323b3e82cSAsbjørn Nilsen Riseth
2449566063dSJacob Faibussowitsch PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA));
24519653cdaSPeter Brune
2469f425c49SPeter Brune /* differences for selection and restart */
24713a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
2489566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
24913a62661SPeter Brune } else {
2509566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
25113a62661SPeter Brune }
252*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
253*76c63389SBarry Smith if (objective) {
254*76c63389SBarry Smith PetscCall(SNESComputeObjective(snes, XA, &objA));
255*76c63389SBarry Smith SNESCheckObjectiveDomainError(snes, objA);
256*76c63389SBarry Smith } else objA = fAnorm;
2571aa26658SKarl Rupp
2589f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
2594936d809SStefano Zampini 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));
260*76c63389SBarry Smith SNESCheckFunctionDomainError(snes, fnorm);
261*76c63389SBarry Smith if (objective) {
262*76c63389SBarry Smith PetscCall(SNESComputeObjective(snes, X, &obj));
263*76c63389SBarry Smith SNESCheckObjectiveDomainError(snes, obj);
264*76c63389SBarry Smith } else obj = fnorm;
26519653cdaSPeter Brune selectRestart = PETSC_FALSE;
26623b3e82cSAsbjørn Nilsen Riseth
26713a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
2684936d809SStefano Zampini PetscCall(SNESNGMRESSelectRestart_Private(snes, l, obj, objM, objA, dnorm, objmin, dminnorm, &selectRestart));
26923b3e82cSAsbjørn Nilsen Riseth
27028ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */
2711aa26658SKarl Rupp if (selectRestart) restart_count++;
2721aa26658SKarl Rupp else restart_count = 0;
27313a62661SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
27413a62661SPeter Brune if (k_restart > ngmres->restart_periodic) {
27563a3b9bcSJacob Faibussowitsch if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart));
27613a62661SPeter Brune restart_count = ngmres->restart_it;
27713a62661SPeter Brune }
27813a62661SPeter Brune }
27923b3e82cSAsbjørn Nilsen Riseth
280b3c6a99cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
28123b3e82cSAsbjørn Nilsen Riseth
28228ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */
28328ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) {
28448a46eb9SPierre Jolivet if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart));
28528ed4a04SPeter Brune restart_count = 0;
28619653cdaSPeter Brune k_restart = 1;
28719653cdaSPeter Brune l = 1;
288b3c6a99cSPeter Brune ivec = 0;
28998b3e84cSPeter Brune /* q_{00} = nu */
2909566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, FM, fMnorm, XM));
291d2e16ddcSPeter Brune } else {
29298b3e84cSPeter Brune /* select the current size of the subspace */
2931e633543SBarry Smith if (l < ngmres->msize) l++;
29419653cdaSPeter Brune k_restart++;
29598b3e84cSPeter Brune /* place the current entry in the list of previous entries */
29638774f0aSPeter Brune if (ngmres->candidate) {
2974936d809SStefano Zampini objmin = PetscMin(objmin, objM);
2989566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fMnorm, XM));
299d2e16ddcSPeter Brune } else {
3004936d809SStefano Zampini objmin = PetscMin(objmin, obj);
3019566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, F, fnorm, X));
30219653cdaSPeter Brune }
303d2e16ddcSPeter Brune }
30419653cdaSPeter Brune
3059566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
306087dfb9eSxuemin snes->iter = k;
307f109b39eSPeter Brune snes->norm = fnorm;
30838606ef4SRené Chenard snes->ynorm = ynorm;
30938606ef4SRené Chenard snes->xnorm = xnorm;
3109566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
3119566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
3122d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, 0, 0, fnorm));
3139566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
3143ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
315a312c225SMatthew G Knepley }
316a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT;
3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
318a312c225SMatthew G Knepley }
319a312c225SMatthew G Knepley
32023b3e82cSAsbjørn Nilsen Riseth /*@
3210b4b7b1cSBarry Smith SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M inside a `SNESNGMRES` solve
32223b3e82cSAsbjørn Nilsen Riseth
32323b3e82cSAsbjørn Nilsen Riseth Input Parameters:
324f6dfbefdSBarry Smith + snes - the `SNES` context.
325f6dfbefdSBarry Smith - flg - boolean value deciding whether to use the option or not, default is `PETSC_FALSE`
32623b3e82cSAsbjørn Nilsen Riseth
327f6dfbefdSBarry Smith Options Database Key:
328f6dfbefdSBarry Smith . -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M
32923b3e82cSAsbjørn Nilsen Riseth
3300b4b7b1cSBarry Smith Level: advanced
33123b3e82cSAsbjørn Nilsen Riseth
33223b3e82cSAsbjørn Nilsen Riseth Notes:
33323b3e82cSAsbjø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.
3340b4b7b1cSBarry Smith To help the solver do that, remove the current stored solutions and residuals whenever F_M increases.
33523b3e82cSAsbjørn Nilsen Riseth
336f6dfbefdSBarry Smith This option must be used with the `SNESNGMRES` `SNESNGMRESRestartType` of `SNES_NGMRES_RESTART_DIFFERENCE`
33723b3e82cSAsbjørn Nilsen Riseth
338420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartType()`
33923b3e82cSAsbjørn Nilsen Riseth @*/
SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg)340d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes, PetscBool flg)
341d71ae5a4SJacob Faibussowitsch {
34223b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES, PetscBool);
34323b3e82cSAsbjørn Nilsen Riseth
34423b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin;
3459566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", &f));
3469566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, flg));
3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
34823b3e82cSAsbjørn Nilsen Riseth }
34923b3e82cSAsbjørn Nilsen Riseth
SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg)35066976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes, PetscBool flg)
351d71ae5a4SJacob Faibussowitsch {
35223b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
35323b3e82cSAsbjørn Nilsen Riseth
35423b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin;
35523b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = flg;
3563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
35723b3e82cSAsbjørn Nilsen Riseth }
35823b3e82cSAsbjørn Nilsen Riseth
SNESNGMRESGetRestartFmRise(SNES snes,PetscBool * flg)359d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes, PetscBool *flg)
360d71ae5a4SJacob Faibussowitsch {
36123b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES, PetscBool *);
36223b3e82cSAsbjørn Nilsen Riseth
36323b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin;
3649566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", &f));
3659566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, flg));
3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
36723b3e82cSAsbjørn Nilsen Riseth }
36823b3e82cSAsbjørn Nilsen Riseth
SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool * flg)36966976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes, PetscBool *flg)
370d71ae5a4SJacob Faibussowitsch {
37123b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
37223b3e82cSAsbjørn Nilsen Riseth
37323b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin;
37423b3e82cSAsbjørn Nilsen Riseth *flg = ngmres->restart_fm_rise;
3753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
37623b3e82cSAsbjørn Nilsen Riseth }
37723b3e82cSAsbjørn Nilsen Riseth
37813a62661SPeter Brune /*@
379f6dfbefdSBarry Smith SNESNGMRESSetRestartType - Sets the restart type for `SNESNGMRES`.
38013a62661SPeter Brune
381c3339decSBarry Smith Logically Collective
38213a62661SPeter Brune
38313a62661SPeter Brune Input Parameters:
38413a62661SPeter Brune + snes - the iterative context
385ceaaa498SBarry Smith - rtype - restart type, see `SNESNGMRESRestartType`
38613a62661SPeter Brune
387da81f932SPierre Jolivet Options Database Keys:
38813a62661SPeter Brune + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
389420bcc1bSBarry Smith - -snes_ngmres_restart <30> - sets the number of iterations before restart for periodic
39013a62661SPeter Brune
39113a62661SPeter Brune Level: intermediate
39213a62661SPeter Brune
393420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartFmRise()`,
394420bcc1bSBarry Smith `SNESNGMRESSetSelectType()`
39513a62661SPeter Brune @*/
SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype)396d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype)
397d71ae5a4SJacob Faibussowitsch {
39813a62661SPeter Brune PetscFunctionBegin;
39913a62661SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
400cac4c232SBarry Smith PetscTryMethod(snes, "SNESNGMRESSetRestartType_C", (SNES, SNESNGMRESRestartType), (snes, rtype));
4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
40213a62661SPeter Brune }
40313a62661SPeter Brune
40413a62661SPeter Brune /*@
405f6dfbefdSBarry Smith SNESNGMRESSetSelectType - Sets the selection type for `SNESNGMRES`. This determines how the candidate solution and
40613a62661SPeter Brune combined solution are used to create the next iterate.
40713a62661SPeter Brune
408c3339decSBarry Smith Logically Collective
40913a62661SPeter Brune
41013a62661SPeter Brune Input Parameters:
41113a62661SPeter Brune + snes - the iterative context
412ceaaa498SBarry Smith - stype - selection type, see `SNESNGMRESSelectType`
41313a62661SPeter Brune
414f6dfbefdSBarry Smith Options Database Key:
41567b8a455SSatish Balay . -snes_ngmres_select_type<difference,none,linesearch> - select type
41613a62661SPeter Brune
41713a62661SPeter Brune Level: intermediate
41813a62661SPeter Brune
419f6dfbefdSBarry Smith Note:
420a99ef635SJonas Heinzmann The default line search used is the `SNESLINESEARCHSECANT` line search and it requires two additional function evaluations.
42113a62661SPeter Brune
422420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNGMRES`, `SNESNGMRESSelectType`, `SNES_NGMRES_SELECT_NONE`, `SNES_NGMRES_SELECT_DIFFERENCE`, `SNES_NGMRES_SELECT_LINESEARCH`,
423420bcc1bSBarry Smith `SNESNGMRESSetRestartType()`
42413a62661SPeter Brune @*/
SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype)425d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetSelectType(SNES snes, SNESNGMRESSelectType stype)
426d71ae5a4SJacob Faibussowitsch {
42713a62661SPeter Brune PetscFunctionBegin;
42813a62661SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
429cac4c232SBarry Smith PetscTryMethod(snes, "SNESNGMRESSetSelectType_C", (SNES, SNESNGMRESSelectType), (snes, stype));
4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
43113a62661SPeter Brune }
43213a62661SPeter Brune
SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype)43366976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype)
434d71ae5a4SJacob Faibussowitsch {
43513a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
4365fd66863SKarl Rupp
43713a62661SPeter Brune PetscFunctionBegin;
43813a62661SPeter Brune ngmres->select_type = stype;
4393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
44013a62661SPeter Brune }
44113a62661SPeter Brune
SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype)44266976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype)
443d71ae5a4SJacob Faibussowitsch {
44413a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
4455fd66863SKarl Rupp
44613a62661SPeter Brune PetscFunctionBegin;
44713a62661SPeter Brune ngmres->restart_type = rtype;
4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
44913a62661SPeter Brune }
45013a62661SPeter Brune
451dfbf837cSBarry Smith /*MC
4520b4b7b1cSBarry Smith SNESNGMRES - An implementation of the Nonlinear Generalized Minimum Residual method, Nonlinear GMRES, or N-GMRES {cite}`ow1`, {cite}`bruneknepleysmithtu15` for solving
4530b4b7b1cSBarry Smith nonlinear systems with `SNES`.
454a312c225SMatthew G Knepley
455dfbf837cSBarry Smith Level: beginner
456dfbf837cSBarry Smith
457f6dfbefdSBarry Smith Options Database Keys:
45813a62661SPeter Brune + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
45938774f0aSPeter Brune . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
460f6dfbefdSBarry Smith . -snes_ngmres_candidate - Use `SNESNGMRES` variant which combines candidate solutions instead of actual solutions
46113a62661SPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals
46213a62661SPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart
46313a62661SPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination
46413a62661SPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart
46513a62661SPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart
46613a62661SPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart
4670b4b7b1cSBarry Smith . -snes_ngmres_restart_fm_rise - Restart on residual rise from $x_M$ step
4680b4b7b1cSBarry Smith . -snes_ngmres_monitor - Prints relevant information about the nonlinear GNMRES iterations
4695c3e6ab7SPeter Brune . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
4704936d809SStefano Zampini - -snes_ngmres_additive_snes_linesearch_type - line search type used to select between the candidate and combined solution with additive select type
4711867fe5bSPeter Brune
4721867fe5bSPeter Brune Notes:
4731867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
4741867fe5bSPeter Brune optimization problem at each iteration.
4751867fe5bSPeter Brune
476f6dfbefdSBarry Smith Very similar to the `SNESANDERSON` algorithm.
4774f02bc6aSBarry Smith
4780b4b7b1cSBarry Smith Unlike the linear GMRES algorithm this algorithm does not compute a Krylov subspace using the Arnoldi process. Instead it stores a
4790b4b7b1cSBarry Smith collection of previous solutions and the residuals $ F(x) - b $ at those solutions.
4800b4b7b1cSBarry Smith
4810b4b7b1cSBarry Smith This algorithm ignores any Jacobian provided with `SNESSetJacobian()`
4820b4b7b1cSBarry Smith
48362842358SBarry Smith Only supports left non-linear preconditioning.
48462842358SBarry Smith
485420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESANDERSON`, `SNESNGMRESSetSelectType()`, `SNESNGMRESSetRestartType()`,
486f8d70eaaSPierre Jolivet `SNESNGMRESSetRestartFmRise()`, `SNESNGMRESSelectType`, `SNESNGMRESRestartType`
487dfbf837cSBarry Smith M*/
488a312c225SMatthew G Knepley
SNESCreate_NGMRES(SNES snes)489d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
490d71ae5a4SJacob Faibussowitsch {
491a312c225SMatthew G Knepley SNES_NGMRES *ngmres;
492d8d34be6SBarry Smith SNESLineSearch linesearch;
493a312c225SMatthew G Knepley
494a312c225SMatthew G Knepley PetscFunctionBegin;
495a312c225SMatthew G Knepley snes->ops->destroy = SNESDestroy_NGMRES;
496a312c225SMatthew G Knepley snes->ops->setup = SNESSetUp_NGMRES;
497a312c225SMatthew G Knepley snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
498a312c225SMatthew G Knepley snes->ops->view = SNESView_NGMRES;
499a312c225SMatthew G Knepley snes->ops->solve = SNESSolve_NGMRES;
500a312c225SMatthew G Knepley snes->ops->reset = SNESReset_NGMRES;
501a312c225SMatthew G Knepley
502efd4aadfSBarry Smith snes->usesnpc = PETSC_TRUE;
5032c155ee1SBarry Smith snes->usesksp = PETSC_FALSE;
504efd4aadfSBarry Smith snes->npcside = PC_RIGHT;
5052c155ee1SBarry Smith
5064fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE;
5074fc747eaSLawrence Mitchell
5084dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ngmres));
509a312c225SMatthew G Knepley snes->data = (void *)ngmres;
510d2e16ddcSPeter Brune ngmres->msize = 30;
51119653cdaSPeter Brune
51277e5a1f9SBarry Smith PetscCall(SNESParametersInitialize(snes));
51377e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_funcs, 30000);
51477e5a1f9SBarry Smith PetscObjectParameterSetDefault(snes, max_its, 10000);
5150e444f03SPeter Brune
51638774f0aSPeter Brune ngmres->candidate = PETSC_FALSE;
517d2e16ddcSPeter Brune
5189566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch));
51948a46eb9SPierre Jolivet if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));
520d8d34be6SBarry Smith
5210298fd71SBarry Smith ngmres->additive_linesearch = NULL;
522077c4231SPeter Brune ngmres->approxfunc = PETSC_FALSE;
52328ed4a04SPeter Brune ngmres->restart_it = 2;
52413a62661SPeter Brune ngmres->restart_periodic = 30;
525f109b39eSPeter Brune ngmres->gammaA = 2.0;
526f109b39eSPeter Brune ngmres->gammaC = 2.0;
527cac108bcSPeter Brune ngmres->deltaB = 0.9;
528cac108bcSPeter Brune ngmres->epsilonB = 0.1;
52923b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = PETSC_FALSE;
530e7058c64SPeter Brune
53113a62661SPeter Brune ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
53213a62661SPeter Brune ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE;
53313a62661SPeter Brune
5349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", SNESNGMRESSetSelectType_NGMRES));
5359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", SNESNGMRESSetRestartType_NGMRES));
5369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", SNESNGMRESSetRestartFmRise_NGMRES));
5379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", SNESNGMRESGetRestartFmRise_NGMRES));
5383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
539a312c225SMatthew G Knepley }
540