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:   const char  *optionsprefix;
 43:   PetscInt     msize, hsize;
 44:   DM           dm;

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

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

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

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

 77:   /* linesearch setup */
 78:   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));

 80:   if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
 81:     PetscCall(SNESLineSearchCreate(PetscObjectComm((PetscObject)snes), &ngmres->additive_linesearch));
 82:     PetscCall(SNESLineSearchSetSNES(ngmres->additive_linesearch, snes));
 83:     if (!((PetscObject)ngmres->additive_linesearch)->type_name) PetscCall(SNESLineSearchSetType(ngmres->additive_linesearch, SNESLINESEARCHL2));
 84:     PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, "additive_"));
 85:     PetscCall(SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch, optionsprefix));
 86:     PetscCall(SNESLineSearchSetFromOptions(ngmres->additive_linesearch));
 87:   }

 89:   ngmres->setup_called = PETSC_TRUE;
 90:   PetscFunctionReturn(PETSC_SUCCESS);
 91: }

 93: static PetscErrorCode SNESSetFromOptions_NGMRES(SNES snes, PetscOptionItems *PetscOptionsObject)
 94: {
 95:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
 96:   PetscBool    debug  = PETSC_FALSE;

 98:   PetscFunctionBegin;
 99:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NGMRES options");
100:   PetscCall(PetscOptionsEnum("-snes_ngmres_select_type", "Select type", "SNESNGMRESSetSelectType", SNESNGMRESSelectTypes, (PetscEnum)ngmres->select_type, (PetscEnum *)&ngmres->select_type, NULL));
101:   PetscCall(PetscOptionsEnum("-snes_ngmres_restart_type", "Restart type", "SNESNGMRESSetRestartType", SNESNGMRESRestartTypes, (PetscEnum)ngmres->restart_type, (PetscEnum *)&ngmres->restart_type, NULL));
102:   PetscCall(PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES", ngmres->candidate, &ngmres->candidate, NULL));
103:   PetscCall(PetscOptionsBool("-snes_ngmres_approxfunc", "Linearly approximate the function", "SNES", ngmres->approxfunc, &ngmres->approxfunc, NULL));
104:   PetscCall(PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES", ngmres->msize, &ngmres->msize, NULL));
105:   PetscCall(PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES", ngmres->restart_periodic, &ngmres->restart_periodic, NULL));
106:   PetscCall(PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart", "SNES", ngmres->restart_it, &ngmres->restart_it, NULL));
107:   PetscCall(PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES", ngmres->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
108:   if (debug) ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
109:   PetscCall(PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES", ngmres->gammaA, &ngmres->gammaA, NULL));
110:   PetscCall(PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES", ngmres->gammaC, &ngmres->gammaC, NULL));
111:   PetscCall(PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES", ngmres->epsilonB, &ngmres->epsilonB, NULL));
112:   PetscCall(PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES", ngmres->deltaB, &ngmres->deltaB, NULL));
113:   PetscCall(PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES", ngmres->singlereduction, &ngmres->singlereduction, NULL));
114:   PetscCall(PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise", ngmres->restart_fm_rise, &ngmres->restart_fm_rise, NULL));
115:   PetscOptionsHeadEnd();
116:   if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: PetscErrorCode SNESView_NGMRES(SNES snes, PetscViewer viewer)
121: {
122:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
123:   PetscBool    iascii;

125:   PetscFunctionBegin;
126:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
127:   if (iascii) {
128:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of stored past updates: %" PetscInt_FMT "\n", ngmres->msize));
129:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual selection: gammaA=%1.0e, gammaC=%1.0e\n", (double)ngmres->gammaA, (double)ngmres->gammaC));
130:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n", (double)ngmres->epsilonB, (double)ngmres->deltaB));
131:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Restart on F_M residual increase: %s\n", PetscBools[ngmres->restart_fm_rise]));
132:   }
133:   PetscFunctionReturn(PETSC_SUCCESS);
134: }

136: static PetscErrorCode SNESSolve_NGMRES(SNES snes)
137: {
138:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;
139:   /* present solution, residual, and preconditioned residual */
140:   Vec X, F, B, D, Y;

142:   /* candidate linear combination answers */
143:   Vec XA, FA, XM, FM;

145:   /* coefficients and RHS to the minimization problem */
146:   PetscReal fnorm, fMnorm, fAnorm;
147:   PetscReal xnorm, xMnorm, xAnorm;
148:   PetscReal ynorm, yMnorm, yAnorm;
149:   PetscInt  k, k_restart, l, ivec, restart_count = 0;

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

161:   SNESConvergedReason  reason;
162:   SNESLineSearchReason lssucceed;

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->vec_sol_update;
174:   FA           = snes->work[0];
175:   D            = snes->work[1];

177:   /* work for the line search */
178:   Y  = snes->work[2];
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:   fminnorm = fnorm;

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

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

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

228:       PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, XM, B, 0));
229:       PetscCall(SNESSolve(snes->npc, B, XM));
230:       PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, XM, B, 0));

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

244:       fMnorm = fnorm;

246:       PetscCall(SNESLineSearchApply(snes->linesearch, XM, FM, &fMnorm, Y));
247:       PetscCall(SNESLineSearchGetReason(snes->linesearch, &lssucceed));
248:       if (lssucceed) {
249:         if (++snes->numFailures >= snes->maxFailures) {
250:           snes->reason = SNES_DIVERGED_LINE_SEARCH;
251:           PetscFunctionReturn(PETSC_SUCCESS);
252:         }
253:       }
254:     }

256:     PetscCall(SNESNGMRESFormCombinedSolution_Private(snes, ivec, l, XM, FM, fMnorm, X, XA, FA));
257:     /* r = F(x) */
258:     if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */

260:     /* differences for selection and restart */
261:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
262:       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, &dnorm, &dminnorm, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
263:     } else {
264:       PetscCall(SNESNGMRESNorms_Private(snes, l, X, F, XM, FM, XA, FA, D, NULL, NULL, &xMnorm, NULL, &yMnorm, &xAnorm, &fAnorm, &yAnorm));
265:     }
266:     SNESCheckFunctionNorm(snes, fnorm);

268:     /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
269:     PetscCall(SNESNGMRESSelect_Private(snes, k_restart, XM, FM, xMnorm, fMnorm, yMnorm, XA, FA, xAnorm, fAnorm, yAnorm, dnorm, fminnorm, dminnorm, X, F, Y, &xnorm, &fnorm, &ynorm));
270:     selectRestart = PETSC_FALSE;

272:     if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
273:       PetscCall(SNESNGMRESSelectRestart_Private(snes, l, fMnorm, fAnorm, dnorm, fminnorm, dminnorm, &selectRestart));

275:       /* if the restart conditions persist for more than restart_it iterations, restart. */
276:       if (selectRestart) restart_count++;
277:       else restart_count = 0;
278:     } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
279:       if (k_restart > ngmres->restart_periodic) {
280:         if (ngmres->monitor) PetscCall(PetscViewerASCIIPrintf(ngmres->monitor, "periodic restart after %" PetscInt_FMT " iterations\n", k_restart));
281:         restart_count = ngmres->restart_it;
282:       }
283:     }

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

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

310:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
311:     snes->iter = k;
312:     snes->norm = fnorm;
313:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
314:     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, snes->iter));
315:     PetscCall(SNESConverged(snes, snes->iter, 0, 0, fnorm));
316:     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
317:     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
318:   }
319:   snes->reason = SNES_DIVERGED_MAX_IT;
320:   PetscFunctionReturn(PETSC_SUCCESS);
321: }

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

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

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

333:   Level: intermediate

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

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

341: .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartType()`
342:   @*/
343: PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes, PetscBool flg)
344: {
345:   PetscErrorCode (*f)(SNES, PetscBool);

347:   PetscFunctionBegin;
348:   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", &f));
349:   if (f) PetscCall((f)(snes, flg));
350:   PetscFunctionReturn(PETSC_SUCCESS);
351: }

353: static PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes, PetscBool flg)
354: {
355:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

357:   PetscFunctionBegin;
358:   ngmres->restart_fm_rise = flg;
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }

362: PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes, PetscBool *flg)
363: {
364:   PetscErrorCode (*f)(SNES, PetscBool *);

366:   PetscFunctionBegin;
367:   PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", &f));
368:   if (f) PetscCall((f)(snes, flg));
369:   PetscFunctionReturn(PETSC_SUCCESS);
370: }

372: static PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes, PetscBool *flg)
373: {
374:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

376:   PetscFunctionBegin;
377:   *flg = ngmres->restart_fm_rise;
378:   PetscFunctionReturn(PETSC_SUCCESS);
379: }

381: /*@
382:   SNESNGMRESSetRestartType - Sets the restart type for `SNESNGMRES`.

384:   Logically Collective

386:   Input Parameters:
387: + snes  - the iterative context
388: - rtype - restart type, see `SNESNGMRESRestartType`

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

394:   Level: intermediate

396: .seealso: [](ch_snes), `SNES`, `SNES_NGMRES_RESTART_DIFFERENCE`, `SNESNGMRES`, `SNESNGMRESRestartType`, `SNESNGMRESSetRestartFmRise()`,
397:           `SNESNGMRESSetSelectType()`
398: @*/
399: PetscErrorCode SNESNGMRESSetRestartType(SNES snes, SNESNGMRESRestartType rtype)
400: {
401:   PetscFunctionBegin;
403:   PetscTryMethod(snes, "SNESNGMRESSetRestartType_C", (SNES, SNESNGMRESRestartType), (snes, rtype));
404:   PetscFunctionReturn(PETSC_SUCCESS);
405: }

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

411:   Logically Collective

413:   Input Parameters:
414: + snes  - the iterative context
415: - stype - selection type, see `SNESNGMRESSelectType`

417:   Options Database Key:
418: . -snes_ngmres_select_type<difference,none,linesearch> - select type

420:   Level: intermediate

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

425: .seealso: [](ch_snes), `SNES`, `SNESNGMRES`, `SNESNGMRESSelectType`, `SNES_NGMRES_SELECT_NONE`, `SNES_NGMRES_SELECT_DIFFERENCE`, `SNES_NGMRES_SELECT_LINESEARCH`,
426:           `SNESNGMRESSetRestartType()`
427: @*/
428: PetscErrorCode SNESNGMRESSetSelectType(SNES snes, SNESNGMRESSelectType stype)
429: {
430:   PetscFunctionBegin;
432:   PetscTryMethod(snes, "SNESNGMRESSetSelectType_C", (SNES, SNESNGMRESSelectType), (snes, stype));
433:   PetscFunctionReturn(PETSC_SUCCESS);
434: }

436: static PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes, SNESNGMRESSelectType stype)
437: {
438:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

440:   PetscFunctionBegin;
441:   ngmres->select_type = stype;
442:   PetscFunctionReturn(PETSC_SUCCESS);
443: }

445: static PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes, SNESNGMRESRestartType rtype)
446: {
447:   SNES_NGMRES *ngmres = (SNES_NGMRES *)snes->data;

449:   PetscFunctionBegin;
450:   ngmres->restart_type = rtype;
451:   PetscFunctionReturn(PETSC_SUCCESS);
452: }

454: /*MC
455:   SNESNGMRES - The Nonlinear Generalized Minimum Residual method.

457:    Level: beginner

459:    Options Database Keys:
460: +  -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
461: .  -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
462: .  -snes_ngmres_candidate        - Use `SNESNGMRES` variant which combines candidate solutions instead of actual solutions
463: .  -snes_ngmres_m                - Number of stored previous solutions and residuals
464: .  -snes_ngmres_restart_it       - Number of iterations the restart conditions hold before restart
465: .  -snes_ngmres_gammaA           - Residual tolerance for solution select between the candidate and combination
466: .  -snes_ngmres_gammaC           - Residual tolerance for restart
467: .  -snes_ngmres_epsilonB         - Difference tolerance between subsequent solutions triggering restart
468: .  -snes_ngmres_deltaB           - Difference tolerance between residuals triggering restart
469: .  -snes_ngmres_restart_fm_rise  - Restart on residual rise from x_M step
470: .  -snes_ngmres_monitor          - Prints relevant information about the ngmres iteration
471: .  -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
472: -  -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type

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

478:    Very similar to the `SNESANDERSON` algorithm.

480:    References:
481: +  * - C. W. Oosterlee and T. Washio, "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows",
482:    SIAM Journal on Scientific Computing, 21(5), 2000.
483: -  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
484:    SIAM Review, 57(4), 2015

486: .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESANDERSON`, `SNESNGMRESSetSelectType()`, `SNESNGMRESSetRestartType()`,
487:           `SNESNGMRESSetRestartFmRise()`, `SNESNGMRESSelectType`, ``SNESNGMRESRestartType`
488: M*/

490: PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
491: {
492:   SNES_NGMRES   *ngmres;
493:   SNESLineSearch linesearch;

495:   PetscFunctionBegin;
496:   snes->ops->destroy        = SNESDestroy_NGMRES;
497:   snes->ops->setup          = SNESSetUp_NGMRES;
498:   snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
499:   snes->ops->view           = SNESView_NGMRES;
500:   snes->ops->solve          = SNESSolve_NGMRES;
501:   snes->ops->reset          = SNESReset_NGMRES;

503:   snes->usesnpc = PETSC_TRUE;
504:   snes->usesksp = PETSC_FALSE;
505:   snes->npcside = PC_RIGHT;

507:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

509:   PetscCall(PetscNew(&ngmres));
510:   snes->data    = (void *)ngmres;
511:   ngmres->msize = 30;

513:   if (!snes->tolerancesset) {
514:     snes->max_funcs = 30000;
515:     snes->max_its   = 10000;
516:   }

518:   ngmres->candidate = PETSC_FALSE;

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

523:   ngmres->additive_linesearch = NULL;
524:   ngmres->approxfunc          = PETSC_FALSE;
525:   ngmres->restart_it          = 2;
526:   ngmres->restart_periodic    = 30;
527:   ngmres->gammaA              = 2.0;
528:   ngmres->gammaC              = 2.0;
529:   ngmres->deltaB              = 0.9;
530:   ngmres->epsilonB            = 0.1;
531:   ngmres->restart_fm_rise     = PETSC_FALSE;

533:   ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
534:   ngmres->select_type  = SNES_NGMRES_SELECT_DIFFERENCE;

536:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetSelectType_C", SNESNGMRESSetSelectType_NGMRES));
537:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartType_C", SNESNGMRESSetRestartType_NGMRES));
538:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESSetRestartFmRise_C", SNESNGMRESSetRestartFmRise_NGMRES));
539:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNGMRESGetRestartFmRise_C", SNESNGMRESGetRestartFmRise_NGMRES));
540:   PetscFunctionReturn(PETSC_SUCCESS);
541: }