Actual source code: snesngmres.c
petsc-3.11.4 2019-09-28
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_",0};
6: const char *const SNESNGMRESSelectTypes[] = {"NONE","DIFFERENCE","LINESEARCH","SNESNGMRESSelectType","SNES_NGMRES_SELECT_",0};
8: PetscErrorCode SNESReset_NGMRES(SNES snes)
9: {
10: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
14: VecDestroyVecs(ngmres->msize,&ngmres->Fdot);
15: VecDestroyVecs(ngmres->msize,&ngmres->Xdot);
16: SNESLineSearchDestroy(&ngmres->additive_linesearch);
17: return(0);
18: }
20: PetscErrorCode SNESDestroy_NGMRES(SNES snes)
21: {
23: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
26: SNESReset_NGMRES(snes);
27: PetscFree5(ngmres->h,ngmres->beta,ngmres->xi,ngmres->fnorms,ngmres->q);
28: PetscFree(ngmres->s);
29: PetscFree(ngmres->xnorms);
30: #if defined(PETSC_USE_COMPLEX)
31: PetscFree(ngmres->rwork);
32: #endif
33: PetscFree(ngmres->work);
34: PetscFree(snes->data);
35: return(0);
36: }
38: PetscErrorCode SNESSetUp_NGMRES(SNES snes)
39: {
40: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
41: const char *optionsprefix;
42: PetscInt msize,hsize;
44: DM dm;
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: SNESSetWorkVecs(snes,5);
53: if (!snes->vec_sol) {
54: SNESGetDM(snes,&dm);
55: DMCreateGlobalVector(dm,&snes->vec_sol);
56: }
58: if (!ngmres->Xdot) {VecDuplicateVecs(snes->vec_sol,ngmres->msize,&ngmres->Xdot);}
59: if (!ngmres->Fdot) {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: PetscMalloc5(hsize,&ngmres->h, msize,&ngmres->beta, msize,&ngmres->xi, msize,&ngmres->fnorms, hsize,&ngmres->q);
66: PetscMalloc1(msize,&ngmres->xnorms);
67: ngmres->nrhs = 1;
68: ngmres->lda = msize;
69: ngmres->ldb = msize;
70: PetscMalloc1(msize,&ngmres->s);
71: PetscMemzero(ngmres->h, hsize*sizeof(PetscScalar));
72: PetscMemzero(ngmres->q, hsize*sizeof(PetscScalar));
73: PetscMemzero(ngmres->xi, msize*sizeof(PetscScalar));
74: PetscMemzero(ngmres->beta,msize*sizeof(PetscScalar));
75: ngmres->lwork = 12*msize;
76: #if defined(PETSC_USE_COMPLEX)
77: PetscMalloc1(ngmres->lwork,&ngmres->rwork);
78: #endif
79: PetscMalloc1(ngmres->lwork,&ngmres->work);
80: }
82: /* linesearch setup */
83: SNESGetOptionsPrefix(snes,&optionsprefix);
85: if (ngmres->select_type == SNES_NGMRES_SELECT_LINESEARCH) {
86: SNESLineSearchCreate(PetscObjectComm((PetscObject)snes),&ngmres->additive_linesearch);
87: SNESLineSearchSetSNES(ngmres->additive_linesearch,snes);
88: SNESLineSearchSetType(ngmres->additive_linesearch,SNESLINESEARCHL2);
89: SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,"additive_");
90: SNESLineSearchAppendOptionsPrefix(ngmres->additive_linesearch,optionsprefix);
91: SNESLineSearchSetFromOptions(ngmres->additive_linesearch);
92: }
94: ngmres->setup_called = PETSC_TRUE;
95: return(0);
96: }
98: PetscErrorCode SNESSetFromOptions_NGMRES(PetscOptionItems *PetscOptionsObject,SNES snes)
99: {
100: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
102: PetscBool debug = PETSC_FALSE;
103: SNESLineSearch linesearch;
106: PetscOptionsHead(PetscOptionsObject,"SNES NGMRES options");
107: PetscOptionsEnum("-snes_ngmres_select_type","Select type","SNESNGMRESSetSelectType",SNESNGMRESSelectTypes,
108: (PetscEnum)ngmres->select_type,(PetscEnum*)&ngmres->select_type,NULL);
109: PetscOptionsEnum("-snes_ngmres_restart_type","Restart type","SNESNGMRESSetRestartType",SNESNGMRESRestartTypes,
110: (PetscEnum)ngmres->restart_type,(PetscEnum*)&ngmres->restart_type,NULL);
111: PetscOptionsBool("-snes_ngmres_candidate", "Use candidate storage", "SNES",ngmres->candidate,&ngmres->candidate,NULL);
112: PetscOptionsBool("-snes_ngmres_approxfunc","Linearly approximate the function", "SNES",ngmres->approxfunc,&ngmres->approxfunc,NULL);
113: PetscOptionsInt("-snes_ngmres_m", "Number of directions", "SNES",ngmres->msize,&ngmres->msize,NULL);
114: PetscOptionsInt("-snes_ngmres_restart", "Iterations before forced restart", "SNES",ngmres->restart_periodic,&ngmres->restart_periodic,NULL);
115: PetscOptionsInt("-snes_ngmres_restart_it", "Tolerance iterations before restart","SNES",ngmres->restart_it,&ngmres->restart_it,NULL);
116: PetscOptionsBool("-snes_ngmres_monitor", "Monitor actions of NGMRES", "SNES",ngmres->monitor ? PETSC_TRUE : PETSC_FALSE,&debug,NULL);
117: if (debug) {
118: ngmres->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
119: }
120: PetscOptionsReal("-snes_ngmres_gammaA", "Residual selection constant", "SNES",ngmres->gammaA,&ngmres->gammaA,NULL);
121: PetscOptionsReal("-snes_ngmres_gammaC", "Residual restart constant", "SNES",ngmres->gammaC,&ngmres->gammaC,NULL);
122: PetscOptionsReal("-snes_ngmres_epsilonB", "Difference selection constant", "SNES",ngmres->epsilonB,&ngmres->epsilonB,NULL);
123: PetscOptionsReal("-snes_ngmres_deltaB", "Difference residual selection constant", "SNES",ngmres->deltaB,&ngmres->deltaB,NULL);
124: PetscOptionsBool("-snes_ngmres_single_reduction", "Aggregate reductions", "SNES",ngmres->singlereduction,&ngmres->singlereduction,NULL);
125: PetscOptionsBool("-snes_ngmres_restart_fm_rise", "Restart on F_M residual rise", "SNESNGMRESSetRestartFmRise",ngmres->restart_fm_rise,&ngmres->restart_fm_rise,NULL);
126: PetscOptionsTail();
127: if ((ngmres->gammaA > ngmres->gammaC) && (ngmres->gammaC > 2.)) ngmres->gammaC = ngmres->gammaA;
129: /* set the default type of the line search if the user hasn't already. */
130: if (!snes->linesearch) {
131: SNESGetLineSearch(snes,&linesearch);
132: SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);
133: }
134: return(0);
135: }
137: PetscErrorCode SNESView_NGMRES(SNES snes,PetscViewer viewer)
138: {
139: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
140: PetscBool iascii;
144: PetscObjectTypeCompare((PetscObject) viewer,PETSCVIEWERASCII,&iascii);
145: if (iascii) {
146: PetscViewerASCIIPrintf(viewer," Number of stored past updates: %d\n", ngmres->msize);
147: PetscViewerASCIIPrintf(viewer," Residual selection: gammaA=%1.0e, gammaC=%1.0e\n",ngmres->gammaA,ngmres->gammaC);
148: PetscViewerASCIIPrintf(viewer," Difference restart: epsilonB=%1.0e, deltaB=%1.0e\n",ngmres->epsilonB,ngmres->deltaB);
149: PetscViewerASCIIPrintf(viewer," Restart on F_M residual increase: %s\n",ngmres->restart_fm_rise?"TRUE":"FALSE");
150: }
151: return(0);
152: }
154: PetscErrorCode SNESSolve_NGMRES(SNES snes)
155: {
156: SNES_NGMRES *ngmres = (SNES_NGMRES*) snes->data;
157: /* present solution, residual, and preconditioned residual */
158: Vec X,F,B,D,Y;
160: /* candidate linear combination answers */
161: Vec XA,FA,XM,FM;
163: /* coefficients and RHS to the minimization problem */
164: PetscReal fnorm,fMnorm,fAnorm;
165: PetscReal xnorm,xMnorm,xAnorm;
166: PetscReal ynorm,yMnorm,yAnorm;
167: PetscInt k,k_restart,l,ivec,restart_count = 0;
169: /* solution selection data */
170: PetscBool selectRestart;
171: /*
172: These two variables are initialized to prevent compilers/analyzers from producing false warnings about these variables being passed
173: to SNESNGMRESSelect_Private() without being set when SNES_NGMRES_RESTART_DIFFERENCE, the values are not used in the subroutines in that case
174: so the code is correct as written.
175: */
176: PetscReal dnorm = 0.0,dminnorm = 0.0;
177: PetscReal fminnorm;
179: SNESConvergedReason reason;
180: SNESLineSearchReason lssucceed;
181: PetscErrorCode ierr;
184: if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
186: PetscCitationsRegister(SNESCitation,&SNEScite);
187: /* variable initialization */
188: snes->reason = SNES_CONVERGED_ITERATING;
189: X = snes->vec_sol;
190: F = snes->vec_func;
191: B = snes->vec_rhs;
192: XA = snes->vec_sol_update;
193: FA = snes->work[0];
194: D = snes->work[1];
196: /* work for the line search */
197: Y = snes->work[2];
198: XM = snes->work[3];
199: FM = snes->work[4];
201: PetscObjectSAWsTakeAccess((PetscObject)snes);
202: snes->iter = 0;
203: snes->norm = 0.;
204: PetscObjectSAWsGrantAccess((PetscObject)snes);
206: /* initialization */
208: if (snes->npc && snes->npcside== PC_LEFT) {
209: SNESApplyNPC(snes,X,NULL,F);
210: SNESGetConvergedReason(snes->npc,&reason);
211: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
212: snes->reason = SNES_DIVERGED_INNER;
213: return(0);
214: }
215: VecNorm(F,NORM_2,&fnorm);
216: } else {
217: if (!snes->vec_func_init_set) {
218: SNESComputeFunction(snes,X,F);
219: } else snes->vec_func_init_set = PETSC_FALSE;
221: VecNorm(F,NORM_2,&fnorm);
222: SNESCheckFunctionNorm(snes,fnorm);
223: }
224: fminnorm = fnorm;
226: PetscObjectSAWsTakeAccess((PetscObject)snes);
227: snes->norm = fnorm;
228: PetscObjectSAWsGrantAccess((PetscObject)snes);
229: SNESLogConvergenceHistory(snes,fnorm,0);
230: SNESMonitor(snes,0,fnorm);
231: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
232: if (snes->reason) return(0);
233: SNESNGMRESUpdateSubspace_Private(snes,0,0,F,fnorm,X);
235: k_restart = 1;
236: l = 1;
237: ivec = 0;
238: for (k=1; k < snes->max_its+1; k++) {
239: /* Computation of x^M */
240: if (snes->npc && snes->npcside== PC_RIGHT) {
241: VecCopy(X,XM);
242: SNESSetInitialFunction(snes->npc,F);
244: PetscLogEventBegin(SNES_NPCSolve,snes->npc,XM,B,0);
245: SNESSolve(snes->npc,B,XM);
246: PetscLogEventEnd(SNES_NPCSolve,snes->npc,XM,B,0);
248: SNESGetConvergedReason(snes->npc,&reason);
249: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
250: snes->reason = SNES_DIVERGED_INNER;
251: return(0);
252: }
253: SNESGetNPCFunction(snes,FM,&fMnorm);
254: } else {
255: /* no preconditioner -- just take gradient descent with line search */
256: VecCopy(F,Y);
257: VecCopy(F,FM);
258: VecCopy(X,XM);
260: fMnorm = fnorm;
262: SNESLineSearchApply(snes->linesearch,XM,FM,&fMnorm,Y);
263: SNESLineSearchGetReason(snes->linesearch,&lssucceed);
264: if (lssucceed) {
265: if (++snes->numFailures >= snes->maxFailures) {
266: snes->reason = SNES_DIVERGED_LINE_SEARCH;
267: return(0);
268: }
269: }
270: }
272: SNESNGMRESFormCombinedSolution_Private(snes,ivec,l,XM,FM,fMnorm,X,XA,FA);
273: /* r = F(x) */
274: if (fminnorm > fMnorm) fminnorm = fMnorm; /* the minimum norm is now of F^M */
276: /* differences for selection and restart */
277: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE || ngmres->select_type == SNES_NGMRES_SELECT_DIFFERENCE) {
278: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,&dnorm,&dminnorm,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);
279: } else {
280: SNESNGMRESNorms_Private(snes,l,X,F,XM,FM,XA,FA,D,NULL,NULL,&xMnorm,NULL,&yMnorm,&xAnorm,&fAnorm,&yAnorm);
281: }
282: SNESCheckFunctionNorm(snes,fnorm);
284: /* combination (additive) or selection (multiplicative) of the N-GMRES solution */
285: SNESNGMRESSelect_Private(snes,k_restart,XM,FM,xMnorm,fMnorm,yMnorm,XA,FA,xAnorm,fAnorm,yAnorm,dnorm,fminnorm,dminnorm,X,F,Y,&xnorm,&fnorm,&ynorm);
286: selectRestart = PETSC_FALSE;
288: if (ngmres->restart_type == SNES_NGMRES_RESTART_DIFFERENCE) {
289: SNESNGMRESSelectRestart_Private(snes,l,fMnorm,fAnorm,dnorm,fminnorm,dminnorm,&selectRestart);
291: /* if the restart conditions persist for more than restart_it iterations, restart. */
292: if (selectRestart) restart_count++;
293: else restart_count = 0;
294: } else if (ngmres->restart_type == SNES_NGMRES_RESTART_PERIODIC) {
295: if (k_restart > ngmres->restart_periodic) {
296: if (ngmres->monitor) PetscViewerASCIIPrintf(ngmres->monitor,"periodic restart after %D iterations\n",k_restart);
297: restart_count = ngmres->restart_it;
298: }
299: }
301: ivec = k_restart % ngmres->msize; /* replace the last used part of the subspace */
303: /* restart after restart conditions have persisted for a fixed number of iterations */
304: if (restart_count >= ngmres->restart_it) {
305: if (ngmres->monitor) {
306: PetscViewerASCIIPrintf(ngmres->monitor,"Restarted at iteration %d\n",k_restart);
307: }
308: restart_count = 0;
309: k_restart = 1;
310: l = 1;
311: ivec = 0;
312: /* q_{00} = nu */
313: SNESNGMRESUpdateSubspace_Private(snes,0,0,FM,fMnorm,XM);
314: } else {
315: /* select the current size of the subspace */
316: if (l < ngmres->msize) l++;
317: k_restart++;
318: /* place the current entry in the list of previous entries */
319: if (ngmres->candidate) {
320: if (fminnorm > fMnorm) fminnorm = fMnorm;
321: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,FM,fMnorm,XM);
322: } else {
323: if (fminnorm > fnorm) fminnorm = fnorm;
324: SNESNGMRESUpdateSubspace_Private(snes,ivec,l,F,fnorm,X);
325: }
326: }
328: PetscObjectSAWsTakeAccess((PetscObject)snes);
329: snes->iter = k;
330: snes->norm = fnorm;
331: PetscObjectSAWsGrantAccess((PetscObject)snes);
332: SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
333: SNESMonitor(snes,snes->iter,snes->norm);
334: (*snes->ops->converged)(snes,snes->iter,0,0,fnorm,&snes->reason,snes->cnvP);
335: if (snes->reason) return(0);
336: }
337: snes->reason = SNES_DIVERGED_MAX_IT;
338: return(0);
339: }
341: /*@
342: SNESNGMRESSetRestartFmRise - Increase the restart count if the step x_M increases the residual F_M
344: Input Parameters:
345: + snes - the SNES context.
346: - flg - boolean value deciding whether to use the option or not
348: Options Database:
349: + -snes_ngmres_restart_fm_rise - Increase the restart count if the step x_M increases the residual F_M
351: Level: intermediate
353: Notes:
354: If the proposed step x_M increases the residual F_M, it might be trying to get out of a stagnation area.
355: To help the solver do that, reset the Krylov subspace whenever F_M increases.
357: This option must be used with SNES_NGMRES_RESTART_DIFFERENCE
359: The default is FALSE.
360: .seealso: SNES_NGMRES_RESTART_DIFFERENCE
361: @*/
362: PetscErrorCode SNESNGMRESSetRestartFmRise(SNES snes,PetscBool flg)
363: {
364: PetscErrorCode (*f)(SNES,PetscBool);
368: PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",&f);
369: if (f) {(f)(snes,flg);}
370: return(0);
371: }
373: PetscErrorCode SNESNGMRESSetRestartFmRise_NGMRES(SNES snes,PetscBool flg)
374: {
375: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
378: ngmres->restart_fm_rise = flg;
379: return(0);
380: }
382: PetscErrorCode SNESNGMRESGetRestartFmRise(SNES snes,PetscBool *flg)
383: {
384: PetscErrorCode (*f)(SNES,PetscBool*);
388: PetscObjectQueryFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",&f);
389: if (f) {(f)(snes,flg);}
390: return(0);
391: }
393: PetscErrorCode SNESNGMRESGetRestartFmRise_NGMRES(SNES snes,PetscBool *flg)
394: {
395: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
398: *flg = ngmres->restart_fm_rise;
399: return(0);
400: }
403: /*@
404: SNESNGMRESSetRestartType - Sets the restart type for SNESNGMRES.
406: Logically Collective on SNES
408: Input Parameters:
409: + snes - the iterative context
410: - rtype - restart type
412: Options Database:
413: + -snes_ngmres_restart_type<difference,periodic,none> - set the restart type
414: - -snes_ngmres_restart[30] - sets the number of iterations before restart for periodic
416: Level: intermediate
418: SNESNGMRESRestartTypes:
419: + SNES_NGMRES_RESTART_NONE - never restart
420: . SNES_NGMRES_RESTART_DIFFERENCE - restart based upon difference criteria
421: - SNES_NGMRES_RESTART_PERIODIC - restart after a fixed number of iterations
423: Notes:
424: The default line search used is the L2 line search and it requires two additional function evaluations.
426: .keywords: SNES, SNESNGMRES, restart, type, set SNESLineSearch
427: @*/
428: PetscErrorCode SNESNGMRESSetRestartType(SNES snes,SNESNGMRESRestartType rtype)
429: {
434: PetscTryMethod(snes,"SNESNGMRESSetRestartType_C",(SNES,SNESNGMRESRestartType),(snes,rtype));
435: return(0);
436: }
438: /*@
439: SNESNGMRESSetSelectType - Sets the selection type for SNESNGMRES. This determines how the candidate solution and
440: combined solution are used to create the next iterate.
442: Logically Collective on SNES
444: Input Parameters:
445: + snes - the iterative context
446: - stype - selection type
448: Options Database:
449: . -snes_ngmres_select_type<difference,none,linesearch>
451: Level: intermediate
453: SNESNGMRESSelectTypes:
454: + SNES_NGMRES_SELECT_NONE - choose the combined solution all the time
455: . SNES_NGMRES_SELECT_DIFFERENCE - choose based upon the selection criteria
456: - SNES_NGMRES_SELECT_LINESEARCH - choose based upon line search combination
458: Notes:
459: The default line search used is the L2 line search and it requires two additional function evaluations.
461: .keywords: SNES, SNESNGMRES, selection, type, set SNESLineSearch
462: @*/
463: PetscErrorCode SNESNGMRESSetSelectType(SNES snes,SNESNGMRESSelectType stype)
464: {
469: PetscTryMethod(snes,"SNESNGMRESSetSelectType_C",(SNES,SNESNGMRESSelectType),(snes,stype));
470: return(0);
471: }
473: PetscErrorCode SNESNGMRESSetSelectType_NGMRES(SNES snes,SNESNGMRESSelectType stype)
474: {
475: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
478: ngmres->select_type = stype;
479: return(0);
480: }
482: PetscErrorCode SNESNGMRESSetRestartType_NGMRES(SNES snes,SNESNGMRESRestartType rtype)
483: {
484: SNES_NGMRES *ngmres = (SNES_NGMRES*)snes->data;
487: ngmres->restart_type = rtype;
488: return(0);
489: }
491: /*MC
492: SNESNGMRES - The Nonlinear Generalized Minimum Residual method.
494: Level: beginner
496: Options Database:
497: + -snes_ngmres_select_type<difference,none,linesearch> - choose the select between candidate and combined solution
498: . -snes_ngmres_restart_type<difference,none,periodic> - choose the restart conditions
499: . -snes_ngmres_candidate - Use NGMRES variant which combines candidate solutions instead of actual solutions
500: . -snes_ngmres_m - Number of stored previous solutions and residuals
501: . -snes_ngmres_restart_it - Number of iterations the restart conditions hold before restart
502: . -snes_ngmres_gammaA - Residual tolerance for solution select between the candidate and combination
503: . -snes_ngmres_gammaC - Residual tolerance for restart
504: . -snes_ngmres_epsilonB - Difference tolerance between subsequent solutions triggering restart
505: . -snes_ngmres_deltaB - Difference tolerance between residuals triggering restart
506: . -snes_ngmres_restart_fm_rise - Restart on residual rise from x_M step
507: . -snes_ngmres_monitor - Prints relevant information about the ngmres iteration
508: . -snes_linesearch_type <basic,l2,cp> - Line search type used for the default smoother
509: - -additive_snes_linesearch_type - linesearch type used to select between the candidate and combined solution with additive select type
511: Notes:
513: The N-GMRES method combines m previous solutions into a minimum-residual solution by solving a small linearized
514: optimization problem at each iteration.
516: Very similar to the SNESANDERSON algorithm.
518: References:
519: + 1. - C. W. Oosterlee and T. Washio, "Krylov Subspace Acceleration of Nonlinear Multigrid with Application to Recirculating Flows",
520: SIAM Journal on Scientific Computing, 21(5), 2000.
521: - 2. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
522: SIAM Review, 57(4), 2015
525: .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
526: M*/
528: PETSC_EXTERN PetscErrorCode SNESCreate_NGMRES(SNES snes)
529: {
530: SNES_NGMRES *ngmres;
534: snes->ops->destroy = SNESDestroy_NGMRES;
535: snes->ops->setup = SNESSetUp_NGMRES;
536: snes->ops->setfromoptions = SNESSetFromOptions_NGMRES;
537: snes->ops->view = SNESView_NGMRES;
538: snes->ops->solve = SNESSolve_NGMRES;
539: snes->ops->reset = SNESReset_NGMRES;
541: snes->usesnpc = PETSC_TRUE;
542: snes->usesksp = PETSC_FALSE;
543: snes->npcside = PC_RIGHT;
545: snes->alwayscomputesfinalresidual = PETSC_TRUE;
547: PetscNewLog(snes,&ngmres);
548: snes->data = (void*) ngmres;
549: ngmres->msize = 30;
551: if (!snes->tolerancesset) {
552: snes->max_funcs = 30000;
553: snes->max_its = 10000;
554: }
556: ngmres->candidate = PETSC_FALSE;
558: ngmres->additive_linesearch = NULL;
559: ngmres->approxfunc = PETSC_FALSE;
560: ngmres->restart_it = 2;
561: ngmres->restart_periodic = 30;
562: ngmres->gammaA = 2.0;
563: ngmres->gammaC = 2.0;
564: ngmres->deltaB = 0.9;
565: ngmres->epsilonB = 0.1;
566: ngmres->restart_fm_rise = PETSC_FALSE;
568: ngmres->restart_type = SNES_NGMRES_RESTART_DIFFERENCE;
569: ngmres->select_type = SNES_NGMRES_SELECT_DIFFERENCE;
571: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetSelectType_C",SNESNGMRESSetSelectType_NGMRES);
572: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartType_C",SNESNGMRESSetRestartType_NGMRES);
573: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESSetRestartFmRise_C",SNESNGMRESSetRestartFmRise_NGMRES);
574: PetscObjectComposeFunction((PetscObject)snes,"SNESNGMRESGetRestartFmRise_C",SNESNGMRESGetRestartFmRise_NGMRES);
575: return(0);
576: }