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; 4219653cdaSPeter Brune PetscInt msize, hsize; 43fced5a79SAsbjørn Nilsen Riseth DM dm; 44a312c225SMatthew G Knepley 45a312c225SMatthew G Knepley PetscFunctionBegin; 46efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) { 4746159c86SPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNGMRES does not support left preconditioning with unpreconditioned function"); 4846159c86SPeter Brune } 49efd4aadfSBarry Smith if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED; 509566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 5)); 51fced5a79SAsbjørn Nilsen Riseth 52fced5a79SAsbjørn Nilsen Riseth if (!snes->vec_sol) { 539566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 549566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol)); 55fced5a79SAsbjørn Nilsen Riseth } 56fced5a79SAsbjørn Nilsen Riseth 579566063dSJacob Faibussowitsch if (!ngmres->Xdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Xdot)); 589566063dSJacob Faibussowitsch if (!ngmres->Fdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Fdot)); 5978440776SJed Brown if (!ngmres->setup_called) { 60087dfb9eSxuemin msize = ngmres->msize; /* restart size */ 6119653cdaSPeter Brune hsize = msize * msize; 62087dfb9eSxuemin 6398b3e84cSPeter Brune /* explicit least squares minimization solve */ 649566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(hsize, &ngmres->h, msize, &ngmres->beta, msize, &ngmres->xi, hsize, &ngmres->q)); 659566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(msize, &ngmres->xnorms, msize, &ngmres->fnorms, msize, &ngmres->s)); 6619653cdaSPeter Brune ngmres->nrhs = 1; 676497c311SBarry Smith PetscCall(PetscBLASIntCast(msize, &ngmres->lda)); 686497c311SBarry Smith PetscCall(PetscBLASIntCast(msize, &ngmres->ldb)); 696497c311SBarry Smith PetscCall(PetscBLASIntCast(12 * msize, &ngmres->lwork)); 70dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->rwork)); 7219653cdaSPeter Brune #endif 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->work)); 7478440776SJed Brown } 75e7058c64SPeter Brune 7678440776SJed Brown ngmres->setup_called = PETSC_TRUE; 773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78a312c225SMatthew G Knepley } 79a312c225SMatthew G Knepley 8066976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes, PetscOptionItems *PetscOptionsObject) 81d71ae5a4SJacob Faibussowitsch { 82a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 8394ae4db5SBarry Smith PetscBool debug = PETSC_FALSE; 840adebc6cSBarry Smith 85a312c225SMatthew G Knepley PetscFunctionBegin; 86d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNES NGMRES options"); 87d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-snes_ngmres_select_type", "Select type", "SNESNGMRESSetSelectType", SNESNGMRESSelectTypes, (PetscEnum)ngmres->select_type, (PetscEnum *)&ngmres->select_type, NULL)); 88d0609cedSBarry Smith PetscCall(PetscOptionsEnum("-snes_ngmres_restart_type", "Restart type", "SNESNGMRESSetRestartType", SNESNGMRESRestartTypes, (PetscEnum)ngmres->restart_type, (PetscEnum *)&ngmres->restart_type, NULL)); 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES", ngmres->candidate, &ngmres->candidate, NULL)); 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc", "Linearly approximate the function", "SNES", ngmres->approxfunc, &ngmres->approxfunc, NULL)); 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, NULL)); 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES", ngmres->restart_periodic, &ngmres->restart_periodic, NULL)); 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart", "SNES", ngmres->restart_it, &ngmres->restart_it, NULL)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL)); 95ad540459SPierre Jolivet if (debug) ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes)); 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, NULL)); 979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, NULL)); 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, NULL)); 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, NULL)); 1009566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise", ngmres->restart_fm_rise, &ngmres->restart_fm_rise, NULL)); 101d0609cedSBarry Smith PetscOptionsHeadEnd(); 1024936d809SStefano Zampini if (ngmres->gammaA > ngmres->gammaC && ngmres->gammaC > 2.) ngmres->gammaC = ngmres->gammaA; 1034936d809SStefano Zampini if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) { 1044936d809SStefano Zampini PetscCall(SNESNGMRESGetAdditiveLineSearch_Private(snes, &ngmres->additive_linesearch)); 1054936d809SStefano Zampini PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch)); 1064936d809SStefano Zampini } 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 108a312c225SMatthew G Knepley } 109a312c225SMatthew G Knepley 110d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer) 111d71ae5a4SJacob Faibussowitsch { 112a312c225SMatthew G Knepley SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 113a312c225SMatthew G Knepley PetscBool iascii; 114a312c225SMatthew G Knepley 115a312c225SMatthew G Knepley PetscFunctionBegin; 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 117a312c225SMatthew G Knepley if (iascii) { 11863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of stored past updates: %" PetscInt_FMT "\n", ngmres->msize)); 1194936d809SStefano Zampini if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 12063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", (double)ngmres->gammaA, (double)ngmres->gammaC)); 12163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", (double)ngmres->epsilonB, (double)ngmres->deltaB)); 12263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Restart on F_M residual increase: %s\n", PetscBools[ngmres->restart_fm_rise])); 123a312c225SMatthew G Knepley } 1244936d809SStefano Zampini if (ngmres->additive_linesearch) { 1254936d809SStefano Zampini PetscCall(PetscViewerASCIIPrintf(viewer, " Additive line-search details:\n")); 1264936d809SStefano Zampini PetscCall(SNESLineSearchView(ngmres->additive_linesearch, viewer)); 1274936d809SStefano Zampini } 1284936d809SStefano Zampini } 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 130a312c225SMatthew G Knepley } 131a312c225SMatthew G Knepley 13266976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NGMRES(SNES snes) 133d71ae5a4SJacob Faibussowitsch { 134087dfb9eSxuemin SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 1356b72add0SBarry Smith Vec X, F, B, D, Y; /* present solution, residual, and preconditioned residual */ 1366b72add0SBarry Smith Vec XA, FA, XM, FM; /* candidate linear combination answers */ 1376b72add0SBarry Smith PetscReal fnorm, fMnorm, fAnorm; /* coefficients and RHS to the minimization problem */ 138b3c6a99cSPeter Brune PetscReal xnorm, xMnorm, xAnorm; 139b3c6a99cSPeter Brune PetscReal ynorm, yMnorm, yAnorm; 14038774f0aSPeter Brune PetscInt k, k_restart, l, ivec, restart_count = 0; 1416b72add0SBarry Smith PetscReal objmin, objM, objA, obj; /* support for objective functions minimization */ 1426b72add0SBarry Smith PetscBool selectRestart; /* solution selection data */ 1436b72add0SBarry Smith SNESConvergedReason reason; 1446b72add0SBarry Smith SNESLineSearchReason lssucceed; 1456b72add0SBarry Smith SNESObjectiveFn *objective; 14661ba4676SBarry Smith /* 14761ba4676SBarry Smith These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed 14861ba4676SBarry Smith to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case 14961ba4676SBarry Smith so the code is correct as written. 15061ba4676SBarry Smith */ 15161ba4676SBarry Smith PetscReal dnorm = 0.0, dminnorm = 0.0; 15219653cdaSPeter Brune 153a312c225SMatthew G Knepley PetscFunctionBegin; 1540b121fc5SBarry 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); 155c579b300SPatrick Farrell 1569566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 15798b3e84cSPeter Brune /* variable initialization */ 158a312c225SMatthew G Knepley snes->reason = SNES_CONVERGED_ITERATING; 159f109b39eSPeter Brune X = snes->vec_sol; 160f109b39eSPeter Brune F = snes->vec_func; 161f109b39eSPeter Brune B = snes->vec_rhs; 1625ef867c2SRené Chenard XA = snes->work[2]; 163f109b39eSPeter Brune FA = snes->work[0]; 164f109b39eSPeter Brune D = snes->work[1]; 165f109b39eSPeter Brune 166f109b39eSPeter Brune /* work for the line search */ 1675ef867c2SRené Chenard Y = snes->vec_sol_update; 1689f425c49SPeter Brune XM = snes->work[3]; 1699f425c49SPeter Brune FM = snes->work[4]; 170a312c225SMatthew G Knepley 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 172a312c225SMatthew G Knepley snes->iter = 0; 173a312c225SMatthew G Knepley snes->norm = 0.; 1749566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 17519653cdaSPeter Brune 17698b3e84cSPeter Brune /* initialization */ 17719653cdaSPeter Brune 178efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_LEFT) { 1799566063dSJacob Faibussowitsch PetscCall(SNESApplyNPC(snes, X, NULL, F)); 1809566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 1813a2ae377SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1823a2ae377SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1843a2ae377SPeter Brune } 1859566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 1863a2ae377SPeter Brune } else { 187e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 1889566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 1891aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 190c1c75074SPeter Brune 1919566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 192422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 1933a2ae377SPeter Brune } 1944936d809SStefano Zampini PetscCall(SNESGetObjective(snes, &objective, NULL)); 1954936d809SStefano Zampini objmin = fnorm; 1964936d809SStefano Zampini if (objective) PetscCall(SNESComputeObjective(snes, X, &objmin)); 1974936d809SStefano Zampini obj = objmin; 19819653cdaSPeter Brune 1999566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 200f109b39eSPeter Brune snes->norm = fnorm; 2019566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 2029566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 2032d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 2049566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, fnorm)); 2053ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 2063ba16761SJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, F, fnorm, X)); 207a312c225SMatthew G Knepley 20819653cdaSPeter Brune k_restart = 1; 20919653cdaSPeter Brune l = 1; 210b3c6a99cSPeter Brune ivec = 0; 21109c08436SPeter Brune for (k = 1; k < snes->max_its + 1; k++) { 212d52c4f7dSRené Chenard /* Call general purpose update function */ 213d52c4f7dSRené Chenard PetscTryTypeMethod(snes, update, snes->iter); 214d52c4f7dSRené Chenard 21598b3e84cSPeter Brune /* Computation of x^M */ 216efd4aadfSBarry Smith if (snes->npc && snes->npcside == PC_RIGHT) { 2179566063dSJacob Faibussowitsch PetscCall(VecCopy(X, XM)); 2189566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(snes->npc, F)); 21963e7833aSPeter Brune 2209566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0)); 2219566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes->npc, B, XM)); 2229566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0)); 22363e7833aSPeter Brune 2249566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(snes->npc, &reason)); 2258cc86e31SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2268cc86e31SPeter Brune snes->reason = SNES_DIVERGED_INNER; 2273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2288cc86e31SPeter Brune } 2299566063dSJacob Faibussowitsch PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm)); 2308cc86e31SPeter Brune } else { 231f109b39eSPeter Brune /* no preconditioner -- just take gradient descent with line search */ 2329566063dSJacob Faibussowitsch PetscCall(VecCopy(F, Y)); 2339566063dSJacob Faibussowitsch PetscCall(VecCopy(F, FM)); 2349566063dSJacob Faibussowitsch PetscCall(VecCopy(X, XM)); 2351aa26658SKarl Rupp 236e7058c64SPeter Brune fMnorm = fnorm; 2371aa26658SKarl Rupp 2389566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, XM, FM, &fMnorm, Y)); 2399566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed)); 240422a814eSBarry Smith if (lssucceed) { 241f109b39eSPeter Brune if (++snes->numFailures >= snes->maxFailures) { 242f109b39eSPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 244f109b39eSPeter Brune } 245f109b39eSPeter Brune } 2466634f59bSPeter Brune } 2474936d809SStefano Zampini if (objective) PetscCall(SNESComputeObjective(snes, XM, &objM)); 2484936d809SStefano Zampini else objM = fMnorm; 2494936d809SStefano Zampini objmin = PetscMin(objmin, objM); 25023b3e82cSAsbjørn Nilsen Riseth 2519566063dSJacob Faibussowitsch PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA)); 25219653cdaSPeter Brune 2539f425c49SPeter Brune /* differences for selection and restart */ 25413a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) { 2559566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm)); 25613a62661SPeter Brune } else { 2579566063dSJacob Faibussowitsch PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm)); 25813a62661SPeter Brune } 2594936d809SStefano Zampini if (objective) PetscCall(SNESComputeObjective(snes, XA, &objA)); 2604936d809SStefano Zampini else objA = fAnorm; 261422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 2621aa26658SKarl Rupp 2639f425c49SPeter Brune /* combination (additive) or selection (multiplicative) of the N-GMRES solution */ 2644936d809SStefano 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)); 2654936d809SStefano Zampini if (objective) PetscCall(SNESComputeObjective(snes, X, &obj)); 2664936d809SStefano Zampini else obj = fnorm; 26719653cdaSPeter Brune selectRestart = PETSC_FALSE; 26823b3e82cSAsbjørn Nilsen Riseth 26913a62661SPeter Brune if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) { 2704936d809SStefano Zampini PetscCall(SNESNGMRESSelectRestart_Private(snes, l, obj, objM, objA, dnorm, objmin, dminnorm, &selectRestart)); 27123b3e82cSAsbjørn Nilsen Riseth 27228ed4a04SPeter Brune /* if the restart conditions persist for more than restart_it iterations, restart. */ 2731aa26658SKarl Rupp if (selectRestart) restart_count++; 2741aa26658SKarl Rupp else restart_count = 0; 27513a62661SPeter Brune } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) { 27613a62661SPeter Brune if (k_restart > ngmres->restart_periodic) { 27763a3b9bcSJacob Faibussowitsch if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart)); 27813a62661SPeter Brune restart_count = ngmres->restart_it; 27913a62661SPeter Brune } 28013a62661SPeter Brune } 28123b3e82cSAsbjørn Nilsen Riseth 282b3c6a99cSPeter Brune ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */ 28323b3e82cSAsbjørn Nilsen Riseth 28428ed4a04SPeter Brune /* restart after restart conditions have persisted for a fixed number of iterations */ 28528ed4a04SPeter Brune if (restart_count >= ngmres->restart_it) { 28648a46eb9SPierre Jolivet if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart)); 28728ed4a04SPeter Brune restart_count = 0; 28819653cdaSPeter Brune k_restart = 1; 28919653cdaSPeter Brune l = 1; 290b3c6a99cSPeter Brune ivec = 0; 29198b3e84cSPeter Brune /* q_{00} = nu */ 2929566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, FM, fMnorm, XM)); 293d2e16ddcSPeter Brune } else { 29498b3e84cSPeter Brune /* select the current size of the subspace */ 2951e633543SBarry Smith if (l < ngmres->msize) l++; 29619653cdaSPeter Brune k_restart++; 29798b3e84cSPeter Brune /* place the current entry in the list of previous entries */ 29838774f0aSPeter Brune if (ngmres->candidate) { 2994936d809SStefano Zampini objmin = PetscMin(objmin, objM); 3009566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fMnorm, XM)); 301d2e16ddcSPeter Brune } else { 3024936d809SStefano Zampini objmin = PetscMin(objmin, obj); 3039566063dSJacob Faibussowitsch PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, F, fnorm, X)); 30419653cdaSPeter Brune } 305d2e16ddcSPeter Brune } 30619653cdaSPeter Brune 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 308087dfb9eSxuemin snes->iter = k; 309f109b39eSPeter Brune snes->norm = fnorm; 31038606ef4SRené Chenard snes->ynorm = ynorm; 31138606ef4SRené Chenard snes->xnorm = xnorm; 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 3139566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter)); 3142d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, 0, 0, fnorm)); 3159566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 3163ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 317a312c225SMatthew G Knepley } 318a312c225SMatthew G Knepley snes->reason = SNES_DIVERGED_MAX_IT; 3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 320a312c225SMatthew G Knepley } 321a312c225SMatthew G Knepley 32223b3e82cSAsbjørn Nilsen Riseth /*@ 323*0b4b7b1cSBarry Smith SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M inside a `SNESNGMRES` solve 32423b3e82cSAsbjørn Nilsen Riseth 32523b3e82cSAsbjørn Nilsen Riseth Input Parameters: 326f6dfbefdSBarry Smith + snes - the `SNES` context. 327f6dfbefdSBarry Smith - flg - boolean value deciding whether to use the option or not, default is `PETSC_FALSE` 32823b3e82cSAsbjørn Nilsen Riseth 329f6dfbefdSBarry Smith Options Database Key: 330f6dfbefdSBarry Smith . -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M 33123b3e82cSAsbjørn Nilsen Riseth 332*0b4b7b1cSBarry Smith Level: advanced 33323b3e82cSAsbjørn Nilsen Riseth 33423b3e82cSAsbjørn Nilsen Riseth Notes: 33523b3e82cSAsbjø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. 336*0b4b7b1cSBarry Smith To help the solver do that, remove the current stored solutions and residuals whenever F_M increases. 33723b3e82cSAsbjørn Nilsen Riseth 338f6dfbefdSBarry Smith This option must be used with the `SNESNGMRES` `SNESNGMRESRestartType` of `SNES_NGMRES_RESTART_DIFFERENCE` 33923b3e82cSAsbjørn Nilsen Riseth 340420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartType()` 34123b3e82cSAsbjørn Nilsen Riseth @*/ 342d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes, PetscBool flg) 343d71ae5a4SJacob Faibussowitsch { 34423b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES, PetscBool); 34523b3e82cSAsbjørn Nilsen Riseth 34623b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 3479566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", &f)); 3489566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, flg)); 3493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35023b3e82cSAsbjørn Nilsen Riseth } 35123b3e82cSAsbjørn Nilsen Riseth 35266976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes, PetscBool flg) 353d71ae5a4SJacob Faibussowitsch { 35423b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 35523b3e82cSAsbjørn Nilsen Riseth 35623b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 35723b3e82cSAsbjørn Nilsen Riseth ngmres->restart_fm_rise = flg; 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35923b3e82cSAsbjørn Nilsen Riseth } 36023b3e82cSAsbjørn Nilsen Riseth 361d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes, PetscBool *flg) 362d71ae5a4SJacob Faibussowitsch { 36323b3e82cSAsbjørn Nilsen Riseth PetscErrorCode (*f)(SNES, PetscBool *); 36423b3e82cSAsbjørn Nilsen Riseth 36523b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 3669566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", &f)); 3679566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, flg)); 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36923b3e82cSAsbjørn Nilsen Riseth } 37023b3e82cSAsbjørn Nilsen Riseth 37166976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes, PetscBool *flg) 372d71ae5a4SJacob Faibussowitsch { 37323b3e82cSAsbjørn Nilsen Riseth SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 37423b3e82cSAsbjørn Nilsen Riseth 37523b3e82cSAsbjørn Nilsen Riseth PetscFunctionBegin; 37623b3e82cSAsbjørn Nilsen Riseth *flg = ngmres->restart_fm_rise; 3773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37823b3e82cSAsbjørn Nilsen Riseth } 37923b3e82cSAsbjørn Nilsen Riseth 38013a62661SPeter Brune /*@ 381f6dfbefdSBarry Smith SNESNGMRESSetRestartType - Sets the restart type for `SNESNGMRES`. 38213a62661SPeter Brune 383c3339decSBarry Smith Logically Collective 38413a62661SPeter Brune 38513a62661SPeter Brune Input Parameters: 38613a62661SPeter Brune + snes - the iterative context 387ceaaa498SBarry Smith - rtype - restart type, see `SNESNGMRESRestartType` 38813a62661SPeter Brune 389da81f932SPierre Jolivet Options Database Keys: 39013a62661SPeter Brune + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type 391420bcc1bSBarry Smith - -snes_ngmres_restart <30> - sets the number of iterations before restart for periodic 39213a62661SPeter Brune 39313a62661SPeter Brune Level: intermediate 39413a62661SPeter Brune 395420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartFmRise()`, 396420bcc1bSBarry Smith `SNESNGMRESSetSelectType()` 39713a62661SPeter Brune @*/ 398d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype) 399d71ae5a4SJacob Faibussowitsch { 40013a62661SPeter Brune PetscFunctionBegin; 40113a62661SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 402cac4c232SBarry Smith PetscTryMethod(snes, "SNESNGMRESSetRestartType_C", (SNES, SNESNGMRESRestartType), (snes, rtype)); 4033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40413a62661SPeter Brune } 40513a62661SPeter Brune 40613a62661SPeter Brune /*@ 407f6dfbefdSBarry Smith SNESNGMRESSetSelectType - Sets the selection type for `SNESNGMRES`. This determines how the candidate solution and 40813a62661SPeter Brune combined solution are used to create the next iterate. 40913a62661SPeter Brune 410c3339decSBarry Smith Logically Collective 41113a62661SPeter Brune 41213a62661SPeter Brune Input Parameters: 41313a62661SPeter Brune + snes - the iterative context 414ceaaa498SBarry Smith - stype - selection type, see `SNESNGMRESSelectType` 41513a62661SPeter Brune 416f6dfbefdSBarry Smith Options Database Key: 41767b8a455SSatish Balay . -snes_ngmres_select_type<difference,none,linesearch> - select type 41813a62661SPeter Brune 41913a62661SPeter Brune Level: intermediate 42013a62661SPeter Brune 421f6dfbefdSBarry Smith Note: 422f6dfbefdSBarry Smith The default line search used is the `SNESLINESEARCHL2` line search and it requires two additional function evaluations. 42313a62661SPeter Brune 424420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNGMRES`, `SNESNGMRESSelectType`, `SNES_NGMRES_SELECT_NONE`, `SNES_NGMRES_SELECT_DIFFERENCE`, `SNES_NGMRES_SELECT_LINESEARCH`, 425420bcc1bSBarry Smith `SNESNGMRESSetRestartType()` 42613a62661SPeter Brune @*/ 427d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNGMRESSetSelectType(SNES snes, SNESNGMRESSelectType stype) 428d71ae5a4SJacob Faibussowitsch { 42913a62661SPeter Brune PetscFunctionBegin; 43013a62661SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 431cac4c232SBarry Smith PetscTryMethod(snes, "SNESNGMRESSetSelectType_C", (SNES, SNESNGMRESSelectType), (snes, stype)); 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43313a62661SPeter Brune } 43413a62661SPeter Brune 43566976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype) 436d71ae5a4SJacob Faibussowitsch { 43713a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 4385fd66863SKarl Rupp 43913a62661SPeter Brune PetscFunctionBegin; 44013a62661SPeter Brune ngmres->select_type = stype; 4413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44213a62661SPeter Brune } 44313a62661SPeter Brune 44466976f2fSJacob Faibussowitsch static PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype) 445d71ae5a4SJacob Faibussowitsch { 44613a62661SPeter Brune SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data; 4475fd66863SKarl Rupp 44813a62661SPeter Brune PetscFunctionBegin; 44913a62661SPeter Brune ngmres->restart_type = rtype; 4503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45113a62661SPeter Brune } 45213a62661SPeter Brune 453dfbf837cSBarry Smith /*MC 454*0b4b7b1cSBarry Smith SNESNGMRES - An implementation of the Nonlinear Generalized Minimum Residual method, Nonlinear GMRES, or N-GMRES {cite}`ow1`, {cite}`bruneknepleysmithtu15` for solving 455*0b4b7b1cSBarry Smith nonlinear systems with `SNES`. 456a312c225SMatthew G Knepley 457dfbf837cSBarry Smith Level: beginner 458dfbf837cSBarry Smith 459f6dfbefdSBarry Smith Options Database Keys: 46013a62661SPeter Brune + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution 46138774f0aSPeter Brune . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions 462f6dfbefdSBarry Smith . -snes_ngmres_candidate - Use `SNESNGMRES` variant which combines candidate solutions instead of actual solutions 46313a62661SPeter Brune . -snes_ngmres_m - Number of stored previous solutions and residuals 46413a62661SPeter Brune . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart 46513a62661SPeter Brune . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination 46613a62661SPeter Brune . -snes_ngmres_gammaC - Residual tolerance for restart 46713a62661SPeter Brune . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart 46813a62661SPeter Brune . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart 469*0b4b7b1cSBarry Smith . -snes_ngmres_restart_fm_rise - Restart on residual rise from $x_M$ step 470*0b4b7b1cSBarry Smith . -snes_ngmres_monitor - Prints relevant information about the nonlinear GNMRES iterations 4715c3e6ab7SPeter Brune . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother 4724936d809SStefano Zampini - -snes_ngmres_additive_snes_linesearch_type - line search type used to select between the candidate and combined solution with additive select type 4731867fe5bSPeter Brune 4741867fe5bSPeter Brune Notes: 4751867fe5bSPeter Brune The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized 4761867fe5bSPeter Brune optimization problem at each iteration. 4771867fe5bSPeter Brune 478f6dfbefdSBarry Smith Very similar to the `SNESANDERSON` algorithm. 4794f02bc6aSBarry Smith 480*0b4b7b1cSBarry Smith Unlike the linear GMRES algorithm this algorithm does not compute a Krylov subspace using the Arnoldi process. Instead it stores a 481*0b4b7b1cSBarry Smith collection of previous solutions and the residuals $ F(x) - b $ at those solutions. 482*0b4b7b1cSBarry Smith 483*0b4b7b1cSBarry Smith This algorithm ignores any Jacobian provided with `SNESSetJacobian()` 484*0b4b7b1cSBarry Smith 485420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESANDERSON`, `SNESNGMRESSetSelectType()`, `SNESNGMRESSetRestartType()`, 486f8d70eaaSPierre Jolivet `SNESNGMRESSetRestartFmRise()`, `SNESNGMRESSelectType`, `SNESNGMRESRestartType` 487dfbf837cSBarry Smith M*/ 488a312c225SMatthew G Knepley 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