Actual source code: snesngmres.c

  1: #include <../src/snes/impls/ngmres/snesngmres.h>
  2: #include <petscblaslapack.h>
  3: #include <petscdm.h>

  5: const char *const SNESNGMRESRestartTypes[] = {"NONE", "PERIODIC", "DIFFERENCE", "SNESNGMRESRestartType", "SNES_NGMRES_RESTART_", NULL};
  6: const char *const SNESNGMRESSelectTypes[]  = {"NONE", "DIFFERENCE", "LINESEARCH", "SNESNGMRESSelectType", "SNES_NGMRES_SELECT_", NULL};

  8: PetscErrorCode SNESReset_NGMRES(SNES snes)
  9: {
 10:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

 12:   PetscFunctionBegin;
 13:   PetscCall(VecDestroyVecs(ngmres->msize, &ngmres->Fdot));
 14:   PetscCall(VecDestroyVecs(ngmres->msize, &ngmres->Xdot));
 15:   PetscCall(SNESLineSearchDestroy(&ngmres->additive_linesearch));
 16:   PetscFunctionReturn(PETSC_SUCCESS);
 17: }

 19: PetscErrorCode SNESDestroy_NGMRES(SNES snes)
 20: {
 21:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

 23:   PetscFunctionBegin;
 24:   PetscCall(SNESReset_NGMRES(snes));
 25:   PetscCall(PetscFree4(ngmres->h, ngmres->beta, ngmres->xi, ngmres->q));
 26:   PetscCall(PetscFree3(ngmres->xnorms, ngmres->fnorms, ngmres->s));
 27: #if defined(PETSC_USE_COMPLEX)
 28:   PetscCall(PetscFree(ngmres->rwork));
 29: #endif
 30:   PetscCall(PetscFree(ngmres->work));
 31:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", NULL));
 32:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", NULL));
 33:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", NULL));
 34:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", NULL));
 35:   PetscCall(PetscFree(snes->data));
 36:   PetscFunctionReturn(PETSC_SUCCESS);
 37: }

 39: PetscErrorCode SNESSetUp_NGMRES(SNES snes)
 40: {
 41:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
 42:   PetscInt     msize, hsize;
 43:   DM           dm;

 45:   PetscFunctionBegin;
 46:   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
 47:     SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNGMRES does not support left preconditioning with unpreconditioned function");
 48:   }
 49:   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
 50:   PetscCall(SNESSetWorkVecs(snes, 5));

 52:   if (!snes->vec_sol) {
 53:     PetscCall(SNESGetDM(snes, &dm));
 54:     PetscCall(DMCreateGlobalVector(dm, &snes->vec_sol));
 55:   }

 57:   if (!ngmres->Xdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Xdot));
 58:   if (!ngmres->Fdot) PetscCall(VecDuplicateVecs(snes->vec_sol, ngmres->msize, &ngmres->Fdot));
 59:   if (!ngmres->setup_called) {
 60:     msize = ngmres->msize; /* restart size */
 61:     hsize = msize * msize;

 63:     /* explicit least squares minimization solve */
 64:     PetscCall(PetscCalloc4(hsize, &ngmres->h, msize, &ngmres->beta, msize, &ngmres->xi, hsize, &ngmres->q));
 65:     PetscCall(PetscMalloc3(msize, &ngmres->xnorms, msize, &ngmres->fnorms, msize, &ngmres->s));
 66:     ngmres->nrhs  = 1;
 67:     ngmres->lda   = msize;
 68:     ngmres->ldb   = msize;
 69:     ngmres->lwork = 12 * msize;
 70: #if defined(PETSC_USE_COMPLEX)
 71:     PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->rwork));
 72: #endif
 73:     PetscCall(PetscMalloc1(ngmres->lwork, &ngmres->work));
 74:   }

 76:   ngmres->setup_called = PETSC_TRUE;
 77:   PetscFunctionReturn(PETSC_SUCCESS);
 78: }

 80: static PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes, PetscOptionItems *PetscOptionsObject)
 81: {
 82:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
 83:   PetscBool    debug  = PETSC_FALSE;

 85:   PetscFunctionBegin;
 86:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NGMRES options");
 87:   PetscCall(PetscOptionsEnum("-snes_ngmres_select_type", "Select type", "SNESNGMRESSetSelectType", SNESNGMRESSelectTypes, (PetscEnum)ngmres->select_type, (PetscEnum *)&ngmres->select_type, NULL));
 88:   PetscCall(PetscOptionsEnum("-snes_ngmres_restart_type", "Restart type", "SNESNGMRESSetRestartType", SNESNGMRESRestartTypes, (PetscEnum)ngmres->restart_type, (PetscEnum *)&ngmres->restart_type, NULL));
 89:   PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES", ngmres->candidate, &ngmres->candidate, NULL));
 90:   PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc", "Linearly approximate the function", "SNES", ngmres->approxfunc, &ngmres->approxfunc, NULL));
 91:   PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, NULL));
 92:   PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES", ngmres->restart_periodic, &ngmres->restart_periodic, NULL));
 93:   PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart", "SNES", ngmres->restart_it, &ngmres->restart_it, NULL));
 94:   PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
 95:   if (debug) ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
 96:   PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, NULL));
 97:   PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, NULL));
 98:   PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, NULL));
 99:   PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, NULL));
100:   PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise", ngmres->restart_fm_rise, &ngmres->restart_fm_rise, NULL));
101:   PetscOptionsHeadEnd();
102:   if (ngmres->gammaA > ngmres->gammaC && ngmres->gammaC > 2.) ngmres->gammaC = ngmres->gammaA;
103:   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
104:     PetscCall(SNESNGMRESGetAdditiveLineSearch_Private(snes, &ngmres->additive_linesearch));
105:     PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch));
106:   }
107:   PetscFunctionReturn(PETSC_SUCCESS);
108: }

110: PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
111: {
112:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
113:   PetscBool    iascii;

115:   PetscFunctionBegin;
116:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
117:   if (iascii) {
118:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of stored past updates: %" PetscInt_FMT "\n", ngmres->msize));
119:     if (ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
120:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", (double)ngmres->gammaA, (double)ngmres->gammaC));
121:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", (double)ngmres->epsilonB, (double)ngmres->deltaB));
122:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Restart on F_M residual increase: %s\n", PetscBools[ngmres->restart_fm_rise]));
123:     }
124:     if (ngmres->additive_linesearch) {
125:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Additive line-search details:\n"));
126:       PetscCall(SNESLineSearchView(ngmres->additive_linesearch, viewer));
127:     }
128:   }
129:   PetscFunctionReturn(PETSC_SUCCESS);
130: }

132: static PetscErrorCode SNESSolve_NGMRES(SNES snes)
133: {
134:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
135:   /* present solution, residual, and preconditioned residual */
136:   Vec X, F, B, D, Y;

138:   /* candidate linear combination answers */
139:   Vec XA, FA, XM, FM;

141:   /* coefficients and RHS to the minimization problem */
142:   PetscReal fnorm, fMnorm, fAnorm;
143:   PetscReal xnorm, xMnorm, xAnorm;
144:   PetscReal ynorm, yMnorm, yAnorm;
145:   PetscInt  k, k_restart, l, ivec, restart_count = 0;

147:   /* support for objective functions minimization */
148:   PetscReal objmin, objM, objA, obj;

150:   /* solution selection data */
151:   PetscBool selectRestart;
152:   /*
153:       These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed
154:       to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case
155:       so the code is correct as written.
156:   */
157:   PetscReal dnorm = 0.0, dminnorm = 0.0;

159:   SNESConvergedReason  reason;
160:   SNESLineSearchReason lssucceed;

162:   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);

164:   PetscFunctionBegin;
165:   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);

167:   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
168:   /* variable initialization */
169:   snes->reason = SNES_CONVERGED_ITERATING;
170:   X            = snes->vec_sol;
171:   F            = snes->vec_func;
172:   B            = snes->vec_rhs;
173:   XA           = snes->work[2];
174:   FA           = snes->work[0];
175:   D            = snes->work[1];

177:   /* work for the line search */
178:   Y  = snes->vec_sol_update;
179:   XM = snes->work[3];
180:   FM = snes->work[4];

182:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
183:   snes->iter = 0;
184:   snes->norm = 0.;
185:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));

187:   /* initialization */

189:   if (snes->npc && snes->npcside == PC_LEFT) {
190:     PetscCall(SNESApplyNPC(snes, X, NULL, F));
191:     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
192:     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
193:       snes->reason = SNES_DIVERGED_INNER;
194:       PetscFunctionReturn(PETSC_SUCCESS);
195:     }
196:     PetscCall(VecNorm(F, NORM_2, &fnorm));
197:   } else {
198:     if (!snes->vec_func_init_set) {
199:       PetscCall(SNESComputeFunction(snes, X, F));
200:     } else snes->vec_func_init_set = PETSC_FALSE;

202:     PetscCall(VecNorm(F, NORM_2, &fnorm));
203:     SNESCheckFunctionNorm(snes, fnorm);
204:   }
205:   PetscCall(SNESGetObjective(snes, &objective, NULL));
206:   objmin = fnorm;
207:   if (objective) PetscCall(SNESComputeObjective(snes, X, &objmin));
208:   obj = objmin;

210:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
211:   snes->norm = fnorm;
212:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
213:   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
214:   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
215:   PetscCall(SNESMonitor(snes, 0, fnorm));
216:   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
217:   PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, F, fnorm, X));

219:   k_restart = 1;
220:   l         = 1;
221:   ivec      = 0;
222:   for (k = 1; k < snes->max_its + 1; k++) {
223:     /* Call general purpose update function */
224:     PetscTryTypeMethod(snes, update, snes->iter);

226:     /* Computation of x^M */
227:     if (snes->npc && snes->npcside == PC_RIGHT) {
228:       PetscCall(VecCopy(X, XM));
229:       PetscCall(SNESSetInitialFunction(snes->npc, F));

231:       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0));
232:       PetscCall(SNESSolve(snes->npc, B, XM));
233:       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0));

235:       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
236:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
237:         snes->reason = SNES_DIVERGED_INNER;
238:         PetscFunctionReturn(PETSC_SUCCESS);
239:       }
240:       PetscCall(SNESGetNPCFunction(snes, FM, &fMnorm));
241:     } else {
242:       /* no preconditioner -- just take gradient descent with line search */
243:       PetscCall(VecCopy(F, Y));
244:       PetscCall(VecCopy(F, FM));
245:       PetscCall(VecCopy(X, XM));

247:       fMnorm = fnorm;

249:       PetscCall(SNESLineSearchApply(snes->linesearch, XM, FM, &fMnorm, Y));
250:       PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
251:       if (lssucceed) {
252:         if (++snes->numFailures >= snes->maxFailures) {
253:           snes->reason = SNES_DIVERGED_LINE_SEARCH;
254:           PetscFunctionReturn(PETSC_SUCCESS);
255:         }
256:       }
257:     }
258:     if (objective) PetscCall(SNESComputeObjective(snes, XM, &objM));
259:     else objM = fMnorm;
260:     objmin = PetscMin(objmin, objM);

262:     PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA));

264:     /* differences for selection and restart */
265:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
266:       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
267:     } else {
268:       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
269:     }
270:     if (objective) PetscCall(SNESComputeObjective(snes, XA, &objA));
271:     else objA = fAnorm;
272:     SNESCheckFunctionNorm(snes, fnorm);

274:     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
275:     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));
276:     if (objective) PetscCall(SNESComputeObjective(snes, X, &obj));
277:     else obj = fnorm;
278:     selectRestart = PETSC_FALSE;

280:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
281:       PetscCall(SNESNGMRESSelectRestart_Private(snes, l, obj, objM, objA, dnorm, objmin, dminnorm, &selectRestart));

283:       /* if the restart conditions persist for more than restart_it iterations, restart. */
284:       if (selectRestart) restart_count++;
285:       else restart_count = 0;
286:     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
287:       if (k_restart > ngmres->restart_periodic) {
288:         if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart));
289:         restart_count = ngmres->restart_it;
290:       }
291:     }

293:     ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */

295:     /* restart after restart conditions have persisted for a fixed number of iterations */
296:     if (restart_count >= ngmres->restart_it) {
297:       if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "Restarted at iteration %" PetscInt_FMT "\n", k_restart));
298:       restart_count = 0;
299:       k_restart     = 1;
300:       l             = 1;
301:       ivec          = 0;
302:       /* q_{00} = nu */
303:       PetscCall(SNESNGMRESUpdateSubspace_Private(snes, 0, 0, FM, fMnorm, XM));
304:     } else {
305:       /* select the current size of the subspace */
306:       if (l < ngmres->msize) l++;
307:       k_restart++;
308:       /* place the current entry in the list of previous entries */
309:       if (ngmres->candidate) {
310:         objmin = PetscMin(objmin, objM);
311:         PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, FM, fMnorm, XM));
312:       } else {
313:         objmin = PetscMin(objmin, obj);
314:         PetscCall(SNESNGMRESUpdateSubspace_Private(snes, ivec, l, F, fnorm, X));
315:       }
316:     }

318:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
319:     snes->iter  = k;
320:     snes->norm  = fnorm;
321:     snes->ynorm = ynorm;
322:     snes->xnorm = xnorm;
323:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
324:     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
325:     PetscCall(SNESConverged(snes, snes->iter, 0, 0, fnorm));
326:     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
327:     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
328:   }
329:   snes->reason = SNES_DIVERGED_MAX_IT;
330:   PetscFunctionReturn(PETSC_SUCCESS);
331: }

333: /*@
334:   SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M

336:   Input Parameters:
337: + snes - the `SNES` context.
338: - flg  - boolean value deciding whether to use the option or not, default is `PETSC_FALSE`

340:   Options Database Key:
341: . -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M

343:   Level: intermediate

345:   Notes:
346:   If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area.
347:   To help the solver do that, reset the Krylov subspace whenever F_M increases.

349:   This option must be used with the `SNESNGMRES` `SNESNGMRESRestartType` of `SNES_NGMRES_RESTART_DIFFERENCE`

351: .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartType()`
352:   @*/
353: PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes, PetscBool flg)
354: {
355:   PetscErrorCode (*f)(SNES, PetscBool);

357:   PetscFunctionBegin;
358:   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", &f));
359:   if (f) PetscCall((f)(snes, flg));
360:   PetscFunctionReturn(PETSC_SUCCESS);
361: }

363: static PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes, PetscBool flg)
364: {
365:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

367:   PetscFunctionBegin;
368:   ngmres->restart_fm_rise = flg;
369:   PetscFunctionReturn(PETSC_SUCCESS);
370: }

372: PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes, PetscBool *flg)
373: {
374:   PetscErrorCode (*f)(SNES, PetscBool *);

376:   PetscFunctionBegin;
377:   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", &f));
378:   if (f) PetscCall((f)(snes, flg));
379:   PetscFunctionReturn(PETSC_SUCCESS);
380: }

382: static PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes, PetscBool *flg)
383: {
384:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

386:   PetscFunctionBegin;
387:   *flg = ngmres->restart_fm_rise;
388:   PetscFunctionReturn(PETSC_SUCCESS);
389: }

391: /*@
392:   SNESNGMRESSetRestartType - Sets the restart type for `SNESNGMRES`.

394:   Logically Collective

396:   Input Parameters:
397: + snes  - the iterative context
398: - rtype - restart type, see `SNESNGMRESRestartType`

400:   Options Database Keys:
401: + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
402: - -snes_ngmres_restart <30>                           - sets the number of iterations before restart for periodic

404:   Level: intermediate

406: .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartFmRise()`,
407:           `SNESNGMRESSetSelectType()`
408: @*/
409: PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype)
410: {
411:   PetscFunctionBegin;
413:   PetscTryMethod(snes, "SNESNGMRESSetRestartType_C", (SNES, SNESNGMRESRestartType), (snes, rtype));
414:   PetscFunctionReturn(PETSC_SUCCESS);
415: }

417: /*@
418:   SNESNGMRESSetSelectType - Sets the selection type for `SNESNGMRES`.  This determines how the candidate solution and
419:   combined solution are used to create the next iterate.

421:   Logically Collective

423:   Input Parameters:
424: + snes  - the iterative context
425: - stype - selection type, see `SNESNGMRESSelectType`

427:   Options Database Key:
428: . -snes_ngmres_select_type<difference,none,linesearch> - select type

430:   Level: intermediate

432:   Note:
433:   The default line search used is the `SNESLINESEARCHL2` line search and it requires two additional function evaluations.

435: .seealso: [](ch_snes), `SNES`, `SNESNGMRES`, `SNESNGMRESSelectType`, `SNES_NGMRES_SELECT_NONE`, `SNES_NGMRES_SELECT_DIFFERENCE`, `SNES_NGMRES_SELECT_LINESEARCH`,
436:           `SNESNGMRESSetRestartType()`
437: @*/
438: PetscErrorCode SNESNGMRESSetSelectType(SNES snes, SNESNGMRESSelectType stype)
439: {
440:   PetscFunctionBegin;
442:   PetscTryMethod(snes, "SNESNGMRESSetSelectType_C", (SNES, SNESNGMRESSelectType), (snes, stype));
443:   PetscFunctionReturn(PETSC_SUCCESS);
444: }

446: static PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype)
447: {
448:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

450:   PetscFunctionBegin;
451:   ngmres->select_type = stype;
452:   PetscFunctionReturn(PETSC_SUCCESS);
453: }

455: static PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype)
456: {
457:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

459:   PetscFunctionBegin;
460:   ngmres->restart_type = rtype;
461:   PetscFunctionReturn(PETSC_SUCCESS);
462: }

464: /*MC
465:   SNESNGMRES - The Nonlinear Generalized Minimum Residual method {cite}`ow1`, {cite}`bruneknepleysmithtu15`

467:    Level: beginner

469:    Options Database Keys:
470: +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
471: .  -snes_ngmres_restart_type<difference,none,periodic>  - choose the restart conditions
472: .  -snes_ngmres_candidate                               - Use `SNESNGMRES` variant which combines candidate solutions instead of actual solutions
473: .  -snes_ngmres_m                                       - Number of stored previous solutions and residuals
474: .  -snes_ngmres_restart_it                              - Number of iterations the restart conditions hold before restart
475: .  -snes_ngmres_gammaA                                  - Residual tolerance for solution select between the candidate and combination
476: .  -snes_ngmres_gammaC                                  - Residual tolerance for restart
477: .  -snes_ngmres_epsilonB                                - Difference tolerance between subsequent solutions triggering restart
478: .  -snes_ngmres_deltaB                                  - Difference tolerance between residuals triggering restart
479: .  -snes_ngmres_restart_fm_rise                         - Restart on residual rise from x_M step
480: .  -snes_ngmres_monitor                                 - Prints relevant information about the ngmres iteration
481: .  -snes_linesearch_type <basic,l2,cp>                  - Line search type used for the default smoother
482: -  -snes_ngmres_additive_snes_linesearch_type           - line search type used to select between the candidate and combined solution with additive select type

484:    Notes:
485:    The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
486:    optimization problem at each iteration.

488:    Very similar to the `SNESANDERSON` algorithm.

490: .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESANDERSON`, `SNESNGMRESSetSelectType()`, `SNESNGMRESSetRestartType()`,
491:           `SNESNGMRESSetRestartFmRise()`, `SNESNGMRESSelectType`, ``SNESNGMRESRestartType`
492: M*/

494: PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
495: {
496:   SNES_NGMRES   *ngmres;
497:   SNESLineSearch linesearch;

499:   PetscFunctionBegin;
500:   snes->ops->destroy        = SNESDestroy_NGMRES;
501:   snes->ops->setup          = SNESSetUp_NGMRES;
502:   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
503:   snes->ops->view           = SNESView_NGMRES;
504:   snes->ops->solve          = SNESSolve_NGMRES;
505:   snes->ops->reset          = SNESReset_NGMRES;

507:   snes->usesnpc = PETSC_TRUE;
508:   snes->usesksp = PETSC_FALSE;
509:   snes->npcside = PC_RIGHT;

511:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

513:   PetscCall(PetscNew(&ngmres));
514:   snes->data    = (void *)ngmres;
515:   ngmres->msize = 30;

517:   if (!snes->tolerancesset) {
518:     snes->max_funcs = 30000;
519:     snes->max_its   = 10000;
520:   }

522:   ngmres->candidate = PETSC_FALSE;

524:   PetscCall(SNESGetLineSearch(snes, &linesearch));
525:   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC));

527:   ngmres->additive_linesearch = NULL;
528:   ngmres->approxfunc          = PETSC_FALSE;
529:   ngmres->restart_it          = 2;
530:   ngmres->restart_periodic    = 30;
531:   ngmres->gammaA              = 2.0;
532:   ngmres->gammaC              = 2.0;
533:   ngmres->deltaB              = 0.9;
534:   ngmres->epsilonB            = 0.1;
535:   ngmres->restart_fm_rise     = PETSC_FALSE;

537:   ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
538:   ngmres->select_type  = SNES_NGMRES_SELECT_DIFFERENCE;

540:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", SNESNGMRESSetSelectType_NGMRES));
541:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", SNESNGMRESSetRestartType_NGMRES));
542:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", SNESNGMRESSetRestartFmRise_NGMRES));
543:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", SNESNGMRESGetRestartFmRise_NGMRES));
544:   PetscFunctionReturn(PETSC_SUCCESS);
545: }