Actual source code: fas.c
1: /* Defines the basic SNES object */
2: #include <../src/snes/impls/fas/fasimpls.h>
4: const char *const SNESFASTypes[] = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "SNESFASType", "SNES_FAS", NULL};
6: static PetscErrorCode SNESReset_FAS(SNES snes)
7: {
8: SNES_FAS *fas = (SNES_FAS *)snes->data;
10: PetscFunctionBegin;
11: PetscCall(SNESDestroy(&fas->smoothu));
12: PetscCall(SNESDestroy(&fas->smoothd));
13: PetscCall(MatDestroy(&fas->inject));
14: PetscCall(MatDestroy(&fas->interpolate));
15: PetscCall(MatDestroy(&fas->restrct));
16: PetscCall(VecDestroy(&fas->rscale));
17: PetscCall(VecDestroy(&fas->Xg));
18: PetscCall(VecDestroy(&fas->Fg));
19: if (fas->next) PetscCall(SNESReset(fas->next));
20: PetscFunctionReturn(PETSC_SUCCESS);
21: }
23: static PetscErrorCode SNESDestroy_FAS(SNES snes)
24: {
25: SNES_FAS *fas = (SNES_FAS *)snes->data;
27: PetscFunctionBegin;
28: /* recursively resets and then destroys */
29: PetscCall(SNESReset_FAS(snes));
30: PetscCall(SNESDestroy(&fas->next));
31: PetscCall(PetscFree(fas));
32: PetscFunctionReturn(PETSC_SUCCESS);
33: }
35: static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth)
36: {
37: SNESLineSearch linesearch;
38: SNESLineSearch slinesearch;
39: void *lsprectx, *lspostctx;
40: PetscErrorCode (*precheck)(SNESLineSearch, Vec, Vec, PetscBool *, void *);
41: PetscErrorCode (*postcheck)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
43: PetscFunctionBegin;
44: if (!snes->linesearch) PetscFunctionReturn(PETSC_SUCCESS);
45: PetscCall(SNESGetLineSearch(snes, &linesearch));
46: PetscCall(SNESGetLineSearch(smooth, &slinesearch));
47: PetscCall(SNESLineSearchGetPreCheck(linesearch, &precheck, &lsprectx));
48: PetscCall(SNESLineSearchGetPostCheck(linesearch, &postcheck, &lspostctx));
49: PetscCall(SNESLineSearchSetPreCheck(slinesearch, precheck, lsprectx));
50: PetscCall(SNESLineSearchSetPostCheck(slinesearch, postcheck, lspostctx));
51: PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch));
52: PetscFunctionReturn(PETSC_SUCCESS);
53: }
55: static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth)
56: {
57: SNES_FAS *fas = (SNES_FAS *)snes->data;
59: PetscFunctionBegin;
60: PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth));
61: PetscCall(SNESSetFromOptions(smooth));
62: PetscCall(SNESFASSetUpLineSearch_Private(snes, smooth));
64: PetscCall(PetscObjectReference((PetscObject)snes->vec_sol));
65: PetscCall(PetscObjectReference((PetscObject)snes->vec_sol_update));
66: PetscCall(PetscObjectReference((PetscObject)snes->vec_func));
67: smooth->vec_sol = snes->vec_sol;
68: smooth->vec_sol_update = snes->vec_sol_update;
69: smooth->vec_func = snes->vec_func;
71: if (fas->eventsmoothsetup) PetscCall(PetscLogEventBegin(fas->eventsmoothsetup, smooth, 0, 0, 0));
72: PetscCall(SNESSetUp(smooth));
73: if (fas->eventsmoothsetup) PetscCall(PetscLogEventEnd(fas->eventsmoothsetup, smooth, 0, 0, 0));
74: PetscFunctionReturn(PETSC_SUCCESS);
75: }
77: static PetscErrorCode SNESSetUp_FAS(SNES snes)
78: {
79: SNES_FAS *fas = (SNES_FAS *)snes->data;
80: PetscInt dm_levels;
81: SNES next;
82: PetscBool isFine, hasCreateRestriction, hasCreateInjection;
84: PetscFunctionBegin;
85: PetscCall(SNESFASCycleIsFine(snes, &isFine));
86: if (fas->usedmfornumberoflevels && isFine) {
87: PetscCall(DMGetRefineLevel(snes->dm, &dm_levels));
88: dm_levels++;
89: if (dm_levels > fas->levels) {
90: /* reset the number of levels */
91: PetscCall(SNESFASSetLevels(snes, dm_levels, NULL));
92: PetscCall(SNESSetFromOptions(snes));
93: }
94: }
95: PetscCall(SNESFASCycleGetCorrection(snes, &next));
96: if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
98: PetscCall(SNESSetWorkVecs(snes, 2)); /* work vectors used for intergrid transfers */
100: /* set up the smoothers if they haven't already been set up */
101: if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
103: if (snes->dm) {
104: /* set the smoother DMs properly */
105: if (fas->smoothu) PetscCall(SNESSetDM(fas->smoothu, snes->dm));
106: PetscCall(SNESSetDM(fas->smoothd, snes->dm));
107: /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
108: if (next) {
109: /* for now -- assume the DM and the evaluation functions have been set externally */
110: if (!next->dm) {
111: PetscCall(DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm));
112: PetscCall(SNESSetDM(next, next->dm));
113: }
114: /* set the interpolation and restriction from the DM */
115: if (!fas->interpolate) {
116: PetscCall(DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale));
117: if (!fas->restrct) {
118: PetscCall(DMHasCreateRestriction(next->dm, &hasCreateRestriction));
119: /* DM can create restrictions, use that */
120: if (hasCreateRestriction) {
121: PetscCall(DMCreateRestriction(next->dm, snes->dm, &fas->restrct));
122: } else {
123: PetscCall(PetscObjectReference((PetscObject)fas->interpolate));
124: fas->restrct = fas->interpolate;
125: }
126: }
127: }
128: /* set the injection from the DM */
129: if (!fas->inject) {
130: PetscCall(DMHasCreateInjection(next->dm, &hasCreateInjection));
131: if (hasCreateInjection) PetscCall(DMCreateInjection(next->dm, snes->dm, &fas->inject));
132: }
133: }
134: }
136: /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
137: if (fas->galerkin) {
138: if (next) PetscCall(SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next));
139: if (fas->smoothd && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes));
140: if (fas->smoothu && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes));
141: }
143: /* sets the down (pre) smoother's default norm and sets it from options */
144: if (fas->smoothd) {
145: if (fas->level == 0) {
146: PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_ALWAYS));
147: } else {
148: PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY));
149: }
150: PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd));
151: }
153: /* sets the up (post) smoother's default norm and sets it from options */
154: if (fas->smoothu) {
155: if (fas->level != fas->levels - 1) {
156: PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE));
157: } else {
158: PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY));
159: }
160: PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu));
161: }
163: if (next) {
164: /* gotta set up the solution vector for this to work */
165: if (!next->vec_sol) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_sol));
166: if (!next->vec_rhs) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_rhs));
167: PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next));
168: PetscCall(SNESFASSetUpLineSearch_Private(snes, next));
169: PetscCall(SNESSetUp(next));
170: }
172: /* setup FAS work vectors */
173: if (fas->galerkin) {
174: PetscCall(VecDuplicate(snes->vec_sol, &fas->Xg));
175: PetscCall(VecDuplicate(snes->vec_sol, &fas->Fg));
176: }
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: static PetscErrorCode SNESSetFromOptions_FAS(SNES snes, PetscOptionItems *PetscOptionsObject)
181: {
182: SNES_FAS *fas = (SNES_FAS *)snes->data;
183: PetscInt levels = 1;
184: PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE, continuationflg = PETSC_FALSE;
185: SNESFASType fastype;
186: const char *optionsprefix;
187: SNESLineSearch linesearch;
188: PetscInt m, n_up, n_down;
189: SNES next;
190: PetscBool isFine;
192: PetscFunctionBegin;
193: PetscCall(SNESFASCycleIsFine(snes, &isFine));
194: PetscOptionsHeadBegin(PetscOptionsObject, "SNESFAS Options-----------------------------------");
196: /* number of levels -- only process most options on the finest level */
197: if (isFine) {
198: PetscCall(PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg));
199: if (!flg && snes->dm) {
200: PetscCall(DMGetRefineLevel(snes->dm, &levels));
201: levels++;
202: fas->usedmfornumberoflevels = PETSC_TRUE;
203: }
204: PetscCall(SNESFASSetLevels(snes, levels, NULL));
205: fastype = fas->fastype;
206: PetscCall(PetscOptionsEnum("-snes_fas_type", "FAS correction type", "SNESFASSetType", SNESFASTypes, (PetscEnum)fastype, (PetscEnum *)&fastype, &flg));
207: if (flg) PetscCall(SNESFASSetType(snes, fastype));
209: PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
210: PetscCall(PetscOptionsInt("-snes_fas_cycles", "Number of cycles", "SNESFASSetCycles", fas->n_cycles, &m, &flg));
211: if (flg) PetscCall(SNESFASSetCycles(snes, m));
212: PetscCall(PetscOptionsBool("-snes_fas_continuation", "Corrected grid-sequence continuation", "SNESFASSetContinuation", fas->continuation, &continuationflg, &flg));
213: if (flg) PetscCall(SNESFASSetContinuation(snes, continuationflg));
215: PetscCall(PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin", "SNESFASSetGalerkin", fas->galerkin, &galerkinflg, &flg));
216: if (flg) PetscCall(SNESFASSetGalerkin(snes, galerkinflg));
218: if (fas->fastype == SNES_FAS_FULL) {
219: PetscCall(PetscOptionsBool("-snes_fas_full_downsweep", "Smooth on the initial down sweep for full FAS cycles", "SNESFASFullSetDownSweep", fas->full_downsweep, &fas->full_downsweep, &flg));
220: if (flg) PetscCall(SNESFASFullSetDownSweep(snes, fas->full_downsweep));
221: PetscCall(PetscOptionsBool("-snes_fas_full_total", "Use total restriction and interpolaton on the indial down and up sweeps for the full FAS cycle", "SNESFASFullSetUseTotal", fas->full_total, &fas->full_total, &flg));
222: if (flg) PetscCall(SNESFASFullSetTotal(snes, fas->full_total));
223: }
225: PetscCall(PetscOptionsInt("-snes_fas_smoothup", "Number of post-smoothing steps", "SNESFASSetNumberSmoothUp", fas->max_up_it, &n_up, &upflg));
227: PetscCall(PetscOptionsInt("-snes_fas_smoothdown", "Number of pre-smoothing steps", "SNESFASSetNumberSmoothDown", fas->max_down_it, &n_down, &downflg));
229: {
230: PetscViewer viewer;
231: PetscViewerFormat format;
232: PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_fas_monitor", &viewer, &format, &monflg));
233: if (monflg) {
234: PetscViewerAndFormat *vf;
235: PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
236: PetscCall(PetscOptionsRestoreViewer(&viewer));
237: PetscCall(SNESFASSetMonitor(snes, vf, PETSC_TRUE));
238: }
239: }
240: flg = PETSC_FALSE;
241: monflg = PETSC_TRUE;
242: PetscCall(PetscOptionsBool("-snes_fas_log", "Log times for each FAS level", "SNESFASSetLog", monflg, &monflg, &flg));
243: if (flg) PetscCall(SNESFASSetLog(snes, monflg));
244: }
246: PetscOptionsHeadEnd();
248: /* setup from the determined types if there is no pointwise procedure or smoother defined */
249: if (upflg) PetscCall(SNESFASSetNumberSmoothUp(snes, n_up));
250: if (downflg) PetscCall(SNESFASSetNumberSmoothDown(snes, n_down));
252: /* set up the default line search for coarse grid corrections */
253: if (fas->fastype == SNES_FAS_ADDITIVE) {
254: if (!snes->linesearch) {
255: PetscCall(SNESGetLineSearch(snes, &linesearch));
256: PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
257: }
258: }
260: /* recursive option setting for the smoothers */
261: PetscCall(SNESFASCycleGetCorrection(snes, &next));
262: if (next) PetscCall(SNESSetFromOptions(next));
263: PetscFunctionReturn(PETSC_SUCCESS);
264: }
266: #include <petscdraw.h>
267: static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
268: {
269: SNES_FAS *fas = (SNES_FAS *)snes->data;
270: PetscBool isFine, iascii, isdraw;
271: PetscInt i;
272: SNES smoothu, smoothd, levelsnes;
274: PetscFunctionBegin;
275: PetscCall(SNESFASCycleIsFine(snes, &isFine));
276: if (isFine) {
277: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
278: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
279: if (iascii) {
280: PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, levels=%" PetscInt_FMT ", cycles=%" PetscInt_FMT "\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles));
281: if (fas->galerkin) {
282: PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid function evaluation\n"));
283: } else {
284: PetscCall(PetscViewerASCIIPrintf(viewer, " Not using Galerkin computed coarse grid function evaluation\n"));
285: }
286: for (i = 0; i < fas->levels; i++) {
287: PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
288: PetscCall(SNESFASCycleGetSmootherUp(levelsnes, &smoothu));
289: PetscCall(SNESFASCycleGetSmootherDown(levelsnes, &smoothd));
290: if (!i) {
291: PetscCall(PetscViewerASCIIPrintf(viewer, " Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i));
292: } else {
293: PetscCall(PetscViewerASCIIPrintf(viewer, " Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
294: }
295: PetscCall(PetscViewerASCIIPushTab(viewer));
296: if (smoothd) {
297: PetscCall(SNESView(smoothd, viewer));
298: } else {
299: PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n"));
300: }
301: PetscCall(PetscViewerASCIIPopTab(viewer));
302: if (i && (smoothd == smoothu)) {
303: PetscCall(PetscViewerASCIIPrintf(viewer, " Up solver (post-smoother) same as down solver (pre-smoother)\n"));
304: } else if (i) {
305: PetscCall(PetscViewerASCIIPrintf(viewer, " Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
306: PetscCall(PetscViewerASCIIPushTab(viewer));
307: if (smoothu) {
308: PetscCall(SNESView(smoothu, viewer));
309: } else {
310: PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n"));
311: }
312: PetscCall(PetscViewerASCIIPopTab(viewer));
313: }
314: }
315: } else if (isdraw) {
316: PetscDraw draw;
317: PetscReal x, w, y, bottom, th, wth;
318: SNES_FAS *curfas = fas;
319: PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
320: PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
321: PetscCall(PetscDrawStringGetSize(draw, &wth, &th));
322: bottom = y - th;
323: while (curfas) {
324: if (!curfas->smoothu) {
325: PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
326: if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer));
327: PetscCall(PetscDrawPopCurrentPoint(draw));
328: } else {
329: w = 0.5 * PetscMin(1.0 - x, x);
330: PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom));
331: if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer));
332: PetscCall(PetscDrawPopCurrentPoint(draw));
333: PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom));
334: if (curfas->smoothu) PetscCall(SNESView(curfas->smoothu, viewer));
335: PetscCall(PetscDrawPopCurrentPoint(draw));
336: }
337: /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
338: bottom -= 5 * th;
339: if (curfas->next) curfas = (SNES_FAS *)curfas->next->data;
340: else curfas = NULL;
341: }
342: }
343: }
344: PetscFunctionReturn(PETSC_SUCCESS);
345: }
347: /*
348: Defines the action of the downsmoother
349: */
350: static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
351: {
352: SNESConvergedReason reason;
353: Vec FPC;
354: SNES smoothd;
355: PetscBool flg;
356: SNES_FAS *fas = (SNES_FAS *)snes->data;
358: PetscFunctionBegin;
359: PetscCall(SNESFASCycleGetSmootherDown(snes, &smoothd));
360: PetscCall(SNESSetInitialFunction(smoothd, F));
361: if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothd, B, X, 0));
362: PetscCall(SNESSolve(smoothd, B, X));
363: if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothd, B, X, 0));
364: /* check convergence reason for the smoother */
365: PetscCall(SNESGetConvergedReason(smoothd, &reason));
366: if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
367: snes->reason = SNES_DIVERGED_INNER;
368: PetscFunctionReturn(PETSC_SUCCESS);
369: }
371: PetscCall(SNESGetFunction(smoothd, &FPC, NULL, NULL));
372: PetscCall(SNESGetAlwaysComputesFinalResidual(smoothd, &flg));
373: if (!flg) PetscCall(SNESComputeFunction(smoothd, X, FPC));
374: PetscCall(VecCopy(FPC, F));
375: if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm));
376: PetscFunctionReturn(PETSC_SUCCESS);
377: }
379: /*
380: Defines the action of the upsmoother
381: */
382: static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
383: {
384: SNESConvergedReason reason;
385: Vec FPC;
386: SNES smoothu;
387: PetscBool flg;
388: SNES_FAS *fas = (SNES_FAS *)snes->data;
390: PetscFunctionBegin;
391: PetscCall(SNESFASCycleGetSmootherUp(snes, &smoothu));
392: if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothu, 0, 0, 0));
393: PetscCall(SNESSolve(smoothu, B, X));
394: if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothu, 0, 0, 0));
395: /* check convergence reason for the smoother */
396: PetscCall(SNESGetConvergedReason(smoothu, &reason));
397: if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
398: snes->reason = SNES_DIVERGED_INNER;
399: PetscFunctionReturn(PETSC_SUCCESS);
400: }
401: PetscCall(SNESGetFunction(smoothu, &FPC, NULL, NULL));
402: PetscCall(SNESGetAlwaysComputesFinalResidual(smoothu, &flg));
403: if (!flg) PetscCall(SNESComputeFunction(smoothu, X, FPC));
404: PetscCall(VecCopy(FPC, F));
405: if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm));
406: PetscFunctionReturn(PETSC_SUCCESS);
407: }
409: /*@
410: SNESFASCreateCoarseVec - create a `Vec` corresponding to a state vector on one level coarser than the current level
412: Collective
414: Input Parameter:
415: . snes - `SNESFAS` object
417: Output Parameter:
418: . Xcoarse - vector on level one coarser than the current level
420: Level: developer
422: .seealso: [](ch_snes), `SNESFASSetRestriction()`, `SNESFASRestrict()`, `SNESFAS`
423: @*/
424: PetscErrorCode SNESFASCreateCoarseVec(SNES snes, Vec *Xcoarse)
425: {
426: SNES_FAS *fas;
428: PetscFunctionBegin;
430: PetscAssertPointer(Xcoarse, 2);
431: fas = (SNES_FAS *)snes->data;
432: if (fas->rscale) {
433: PetscCall(VecDuplicate(fas->rscale, Xcoarse));
434: } else if (fas->interpolate) {
435: PetscCall(MatCreateVecs(fas->interpolate, Xcoarse, NULL));
436: } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set rscale or interpolation");
437: PetscFunctionReturn(PETSC_SUCCESS);
438: }
440: /*@
441: SNESFASRestrict - restrict a `Vec` to the next coarser level
443: Collective
445: Input Parameters:
446: + fine - `SNES` from which to restrict
447: - Xfine - vector to restrict
449: Output Parameter:
450: . Xcoarse - result of restriction
452: Level: developer
454: .seealso: [](ch_snes), `SNES`, `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASSetInjection()`, `SNESFASCreateCoarseVec()`
455: @*/
456: PetscErrorCode SNESFASRestrict(SNES fine, Vec Xfine, Vec Xcoarse)
457: {
458: SNES_FAS *fas;
460: PetscFunctionBegin;
464: fas = (SNES_FAS *)fine->data;
465: if (fas->inject) {
466: PetscCall(MatRestrict(fas->inject, Xfine, Xcoarse));
467: } else {
468: PetscCall(MatRestrict(fas->restrct, Xfine, Xcoarse));
469: PetscCall(VecPointwiseMult(Xcoarse, fas->rscale, Xcoarse));
470: }
471: PetscFunctionReturn(PETSC_SUCCESS);
472: }
474: /*
475: Performs a variant of FAS using the interpolated total coarse solution
477: fine problem: F(x) = b
478: coarse problem: F^c(x^c) = Rb, Initial guess Rx
479: interpolated solution: x^f = I x^c (total solution interpolation
480: */
481: static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new)
482: {
483: Vec X_c, B_c;
484: SNESConvergedReason reason;
485: SNES next;
486: Mat restrct, interpolate;
487: SNES_FAS *fasc;
489: PetscFunctionBegin;
490: PetscCall(SNESFASCycleGetCorrection(snes, &next));
491: if (next) {
492: fasc = (SNES_FAS *)next->data;
494: PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
495: PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
497: X_c = next->vec_sol;
499: if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
500: /* restrict the total solution: Rb */
501: PetscCall(SNESFASRestrict(snes, X, X_c));
502: B_c = next->vec_rhs;
503: if (snes->vec_rhs) {
504: /* restrict the total rhs defect: Rb */
505: PetscCall(MatRestrict(restrct, snes->vec_rhs, B_c));
506: } else {
507: PetscCall(VecSet(B_c, 0.));
508: }
509: if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
511: PetscCall(SNESSolve(next, B_c, X_c));
512: PetscCall(SNESGetConvergedReason(next, &reason));
513: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
514: snes->reason = SNES_DIVERGED_INNER;
515: PetscFunctionReturn(PETSC_SUCCESS);
516: }
517: /* x^f <- Ix^c*/
518: DM dmc, dmf;
520: PetscCall(SNESGetDM(next, &dmc));
521: PetscCall(SNESGetDM(snes, &dmf));
522: if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
523: PetscCall(DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new));
524: if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
525: PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse solution"));
526: PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view"));
527: PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution"));
528: PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view"));
529: }
530: PetscFunctionReturn(PETSC_SUCCESS);
531: }
533: /*
534: Performs the FAS coarse correction as:
536: fine problem: F(x) = b
537: coarse problem: F^c(x^c) = b^c
539: b^c = F^c(Rx) - R(F(x) - b)
540: */
541: static PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
542: {
543: Vec X_c, Xo_c, F_c, B_c;
544: SNESConvergedReason reason;
545: SNES next;
546: Mat restrct, interpolate;
547: SNES_FAS *fasc;
549: PetscFunctionBegin;
550: PetscCall(SNESFASCycleGetCorrection(snes, &next));
551: if (next) {
552: fasc = (SNES_FAS *)next->data;
554: PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
555: PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
557: X_c = next->vec_sol;
558: Xo_c = next->work[0];
559: F_c = next->vec_func;
560: B_c = next->vec_rhs;
562: if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
563: PetscCall(SNESFASRestrict(snes, X, Xo_c));
564: /* restrict the defect: R(F(x) - b) */
565: PetscCall(MatRestrict(restrct, F, B_c));
566: if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
568: if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0));
569: /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
570: PetscCall(SNESComputeFunction(next, Xo_c, F_c));
571: if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0));
573: /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
574: PetscCall(VecCopy(B_c, X_c));
575: PetscCall(VecCopy(F_c, B_c));
576: PetscCall(VecCopy(X_c, F_c));
577: /* set initial guess of the coarse problem to the projected fine solution */
578: PetscCall(VecCopy(Xo_c, X_c));
580: /* recurse to the next level */
581: PetscCall(SNESSetInitialFunction(next, F_c));
582: PetscCall(SNESSolve(next, B_c, X_c));
583: PetscCall(SNESGetConvergedReason(next, &reason));
584: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
585: snes->reason = SNES_DIVERGED_INNER;
586: PetscFunctionReturn(PETSC_SUCCESS);
587: }
588: /* correct as x <- x + I(x^c - Rx)*/
589: PetscCall(VecAXPY(X_c, -1.0, Xo_c));
591: if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
592: PetscCall(MatInterpolateAdd(interpolate, X_c, X, X_new));
593: if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
594: PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse correction"));
595: PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view"));
596: PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution"));
597: PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view"));
598: }
599: PetscFunctionReturn(PETSC_SUCCESS);
600: }
602: /*
603: The additive cycle is:
605: xhat = x
606: xhat = dS(x, b)
607: x = coarsecorrection(xhat, b_d)
608: x = x + nu*(xhat - x);
609: (optional) x = uS(x, b)
611: With the coarse RHS (defect correction) as below.
612: */
613: static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
614: {
615: Vec F, B, Xhat;
616: Vec X_c, Xo_c, F_c, B_c;
617: SNESConvergedReason reason;
618: PetscReal xnorm, fnorm, ynorm;
619: SNESLineSearchReason lsresult;
620: SNES next;
621: Mat restrct, interpolate;
622: SNES_FAS *fas = (SNES_FAS *)snes->data, *fasc;
624: PetscFunctionBegin;
625: PetscCall(SNESFASCycleGetCorrection(snes, &next));
626: F = snes->vec_func;
627: B = snes->vec_rhs;
628: Xhat = snes->work[1];
629: PetscCall(VecCopy(X, Xhat));
630: /* recurse first */
631: if (next) {
632: fasc = (SNES_FAS *)next->data;
633: PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
634: PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
635: if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0));
636: PetscCall(SNESComputeFunction(snes, Xhat, F));
637: if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0));
638: PetscCall(VecNorm(F, NORM_2, &fnorm));
639: X_c = next->vec_sol;
640: Xo_c = next->work[0];
641: F_c = next->vec_func;
642: B_c = next->vec_rhs;
644: PetscCall(SNESFASRestrict(snes, Xhat, Xo_c));
645: /* restrict the defect */
646: PetscCall(MatRestrict(restrct, F, B_c));
648: /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
649: if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0));
650: PetscCall(SNESComputeFunction(next, Xo_c, F_c));
651: if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0));
652: PetscCall(VecCopy(B_c, X_c));
653: PetscCall(VecCopy(F_c, B_c));
654: PetscCall(VecCopy(X_c, F_c));
655: /* set initial guess of the coarse problem to the projected fine solution */
656: PetscCall(VecCopy(Xo_c, X_c));
658: /* recurse */
659: PetscCall(SNESSetInitialFunction(next, F_c));
660: PetscCall(SNESSolve(next, B_c, X_c));
662: /* smooth on this level */
663: PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &fnorm));
665: PetscCall(SNESGetConvergedReason(next, &reason));
666: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
667: snes->reason = SNES_DIVERGED_INNER;
668: PetscFunctionReturn(PETSC_SUCCESS);
669: }
671: /* correct as x <- x + I(x^c - Rx)*/
672: PetscCall(VecAYPX(X_c, -1.0, Xo_c));
673: PetscCall(MatInterpolate(interpolate, X_c, Xhat));
675: /* additive correction of the coarse direction*/
676: PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat));
677: PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
678: PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm));
679: if (lsresult) {
680: if (++snes->numFailures >= snes->maxFailures) {
681: snes->reason = SNES_DIVERGED_LINE_SEARCH;
682: PetscFunctionReturn(PETSC_SUCCESS);
683: }
684: }
685: } else {
686: PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
687: }
688: PetscFunctionReturn(PETSC_SUCCESS);
689: }
691: /*
692: Defines the FAS cycle as:
694: fine problem: F(x) = b
695: coarse problem: F^c(x) = b^c
697: b^c = F^c(Rx) - R(F(x) - b)
699: correction:
701: x = x + I(x^c - Rx)
702: */
703: static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
704: {
705: Vec F, B;
706: SNES next;
708: PetscFunctionBegin;
709: F = snes->vec_func;
710: B = snes->vec_rhs;
711: /* pre-smooth -- just update using the pre-smoother */
712: PetscCall(SNESFASCycleGetCorrection(snes, &next));
713: PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
714: if (next) {
715: PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
716: PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
717: }
718: PetscFunctionReturn(PETSC_SUCCESS);
719: }
721: static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
722: {
723: SNES next;
724: SNES_FAS *fas = (SNES_FAS *)snes->data;
725: PetscBool isFine;
727: PetscFunctionBegin;
728: /* pre-smooth -- just update using the pre-smoother */
729: PetscCall(SNESFASCycleIsFine(snes, &isFine));
730: PetscCall(SNESFASCycleGetCorrection(snes, &next));
731: fas->full_stage = 0;
732: if (next) PetscCall(SNESFASCycleSetupPhase_Full(next));
733: PetscFunctionReturn(PETSC_SUCCESS);
734: }
736: static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
737: {
738: Vec F, B;
739: SNES_FAS *fas = (SNES_FAS *)snes->data;
740: PetscBool isFine;
741: SNES next;
743: PetscFunctionBegin;
744: F = snes->vec_func;
745: B = snes->vec_rhs;
746: PetscCall(SNESFASCycleIsFine(snes, &isFine));
747: PetscCall(SNESFASCycleGetCorrection(snes, &next));
749: if (isFine) PetscCall(SNESFASCycleSetupPhase_Full(snes));
751: if (fas->full_stage == 0) {
752: /* downsweep */
753: if (next) {
754: if (fas->level != 1) next->max_its += 1;
755: if (fas->full_downsweep) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
756: fas->full_downsweep = PETSC_TRUE;
757: if (fas->full_total) PetscCall(SNESFASInterpolatedCoarseSolution(snes, X, X));
758: else PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
759: fas->full_total = PETSC_FALSE;
760: PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
761: if (fas->level != 1) next->max_its -= 1;
762: } else {
763: /* The smoother on the coarse level is the coarse solver */
764: PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
765: }
766: fas->full_stage = 1;
767: } else if (fas->full_stage == 1) {
768: if (snes->iter == 0) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
769: if (next) {
770: PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
771: PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
772: }
773: }
774: /* final v-cycle */
775: if (isFine) {
776: if (next) {
777: PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
778: PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
779: }
780: }
781: PetscFunctionReturn(PETSC_SUCCESS);
782: }
784: static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
785: {
786: Vec F, B;
787: SNES next;
789: PetscFunctionBegin;
790: F = snes->vec_func;
791: B = snes->vec_rhs;
792: PetscCall(SNESFASCycleGetCorrection(snes, &next));
793: if (next) {
794: PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
795: PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
796: } else {
797: PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
798: }
799: PetscFunctionReturn(PETSC_SUCCESS);
800: }
802: PetscBool SNEScite = PETSC_FALSE;
803: const char SNESCitation[] = "@Article{bruneknepleysmithtu15,"
804: " title = {Composing Scalable Nonlinear Algebraic Solvers},"
805: " author = {Peter R. Brune and Matthew G. Knepley and Barry F. Smith and Xuemin Tu},"
806: " journal = {SIAM Review},"
807: " volume = {57},"
808: " number = {4},"
809: " pages = {535--565},"
810: " doi = {10.1137/130936725},"
811: " year = {2015}\n";
813: static PetscErrorCode SNESSolve_FAS(SNES snes)
814: {
815: PetscInt i;
816: Vec X, F;
817: PetscReal fnorm;
818: SNES_FAS *fas = (SNES_FAS *)snes->data, *ffas;
819: DM dm;
820: PetscBool isFine;
822: PetscFunctionBegin;
823: 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);
825: PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
826: snes->reason = SNES_CONVERGED_ITERATING;
827: X = snes->vec_sol;
828: F = snes->vec_func;
830: PetscCall(SNESFASCycleIsFine(snes, &isFine));
831: /* norm setup */
832: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
833: snes->iter = 0;
834: snes->norm = 0;
835: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
836: if (!snes->vec_func_init_set) {
837: if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0));
838: PetscCall(SNESComputeFunction(snes, X, F));
839: if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0));
840: } else snes->vec_func_init_set = PETSC_FALSE;
842: PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */
843: SNESCheckFunctionNorm(snes, fnorm);
844: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
845: snes->norm = fnorm;
846: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
847: PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
849: /* test convergence */
850: PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
851: PetscCall(SNESMonitor(snes, snes->iter, fnorm));
852: if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
854: if (isFine) {
855: /* propagate scale-dependent data up the hierarchy */
856: PetscCall(SNESGetDM(snes, &dm));
857: for (ffas = fas; ffas->next; ffas = (SNES_FAS *)ffas->next->data) {
858: DM dmcoarse;
859: PetscCall(SNESGetDM(ffas->next, &dmcoarse));
860: PetscCall(DMRestrict(dm, ffas->restrct, ffas->rscale, ffas->inject, dmcoarse));
861: dm = dmcoarse;
862: }
863: }
865: for (i = 0; i < snes->max_its; i++) {
866: /* Call general purpose update function */
867: PetscTryTypeMethod(snes, update, snes->iter);
869: if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
870: PetscCall(SNESFASCycle_Multiplicative(snes, X));
871: } else if (fas->fastype == SNES_FAS_ADDITIVE) {
872: PetscCall(SNESFASCycle_Additive(snes, X));
873: } else if (fas->fastype == SNES_FAS_FULL) {
874: PetscCall(SNESFASCycle_Full(snes, X));
875: } else if (fas->fastype == SNES_FAS_KASKADE) {
876: PetscCall(SNESFASCycle_Kaskade(snes, X));
877: } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported FAS type");
879: /* check for FAS cycle divergence */
880: if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
882: /* Monitor convergence */
883: PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
884: snes->iter = i + 1;
885: PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
886: PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
887: PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, snes->norm));
888: PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
889: if (snes->reason) break;
890: }
891: PetscFunctionReturn(PETSC_SUCCESS);
892: }
894: /*MC
895: SNESFAS - Full Approximation Scheme nonlinear multigrid solver {cite}`bruneknepleysmithtu15`.
897: The nonlinear problem is solved by correction using coarse versions
898: of the nonlinear problem. This problem is perturbed so that a projected
899: solution of the fine problem elicits no correction from the coarse problem.
901: Options Database Keys and Prefixes:
902: + -snes_fas_levels <l> - The number of levels
903: . -snes_fas_cycles <c> - The number of cycles -- 1 for V, 2 for W
904: . -snes_fas_type<additive,multiplicative,full,kaskade> - Additive or multiplicative cycle
905: . -snes_fas_galerkin <false,true> - Form coarse problems by projection back upon the fine problem
906: . -snes_fas_smoothup <u> - The number of iterations of the post-smoother
907: . -snes_fas_smoothdown <d> - The number of iterations of the pre-smoother
908: . -snes_fas_monitor - Monitor progress of all of the levels
909: . -snes_fas_full_downsweep <false,true> - call the downsmooth on the initial downsweep of full FAS
910: . -fas_levels_snes_ - prefix for `SNES` options for all smoothers
911: . -fas_levels_cycle_snes_ - prefix for `SNES` options for all cycles
912: . -fas_levels_i_snes_ - prefix `SNES` options for the smoothers on level i
913: . -fas_levels_i_cycle_snes_ - prefix for `SNES` options for the cycle on level i
914: - -fas_coarse_snes_ - prefix for `SNES` options for the coarsest smoother
916: Level: beginner
918: Note:
919: The organization of the FAS solver is slightly different from the organization of `PCMG`
920: As each level has smoother `SNES` instances(down and potentially up) and a cycle `SNES` instance.
921: The cycle `SNES` instance may be used for monitoring convergence on a particular level.
923: .seealso: [](ch_snes), `PCMG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESFASSetRestriction()`, `SNESFASSetInjection()`,
924: `SNESFASFullGetTotal()`, `SNESFASSetType()`, `SNESFASGetType()`, `SNESFASSetLevels()`, `SNESFASGetLevels()`, `SNESFASGetCycleSNES()`,
925: `SNESFASSetNumberSmoothUp()`, `SNESFASSetNumberSmoothDown()`, `SNESFASSetContinuation()`, `SNESFASSetCycles()`, `SNESFASSetMonitor()`,
926: `SNESFASSetLog()`, `SNESFASCycleSetCycles()`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()`,
927: `SNESFASCycleGetCorrection()`, `SNESFASCycleGetInterpolation()`, `SNESFASCycleGetRestriction()`, `SNESFASCycleGetInjection()`,
928: `SNESFASCycleGetRScale()`, `SNESFASCycleIsFine()`, `SNESFASSetInterpolation()`, `SNESFASGetInterpolation()`, `SNESFASSetRestriction()`,
929: `SNESFASGetRestriction()`, `SNESFASSetInjection()`, `SNESFASGetInjection()`, `SNESFASSetRScale()`,`SNESFASGetSmoother()`,
930: `SNESFASGetSmootherDown()`, `SNESFASGetSmootherUp()`, `SNESFASGetCoarseSolve()`, `SNESFASFullSetDownSweep()`, `SNESFASFullSetTotal()`, `SNESFASFullGetTotal()`
931: M*/
933: PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
934: {
935: SNES_FAS *fas;
937: PetscFunctionBegin;
938: snes->ops->destroy = SNESDestroy_FAS;
939: snes->ops->setup = SNESSetUp_FAS;
940: snes->ops->setfromoptions = SNESSetFromOptions_FAS;
941: snes->ops->view = SNESView_FAS;
942: snes->ops->solve = SNESSolve_FAS;
943: snes->ops->reset = SNESReset_FAS;
945: snes->usesksp = PETSC_FALSE;
946: snes->usesnpc = PETSC_FALSE;
948: if (!snes->tolerancesset) {
949: snes->max_funcs = 30000;
950: snes->max_its = 10000;
951: }
953: snes->alwayscomputesfinalresidual = PETSC_TRUE;
955: PetscCall(PetscNew(&fas));
957: snes->data = (void *)fas;
958: fas->level = 0;
959: fas->levels = 1;
960: fas->n_cycles = 1;
961: fas->max_up_it = 1;
962: fas->max_down_it = 1;
963: fas->smoothu = NULL;
964: fas->smoothd = NULL;
965: fas->next = NULL;
966: fas->previous = NULL;
967: fas->fine = snes;
968: fas->interpolate = NULL;
969: fas->restrct = NULL;
970: fas->inject = NULL;
971: fas->usedmfornumberoflevels = PETSC_FALSE;
972: fas->fastype = SNES_FAS_MULTIPLICATIVE;
973: fas->full_downsweep = PETSC_FALSE;
974: fas->full_total = PETSC_FALSE;
976: fas->eventsmoothsetup = 0;
977: fas->eventsmoothsolve = 0;
978: fas->eventresidual = 0;
979: fas->eventinterprestrict = 0;
980: PetscFunctionReturn(PETSC_SUCCESS);
981: }