xref: /petsc/src/snes/impls/fas/fas.c (revision 2d1571503b53ddc44d5587ea811639ee03673e35)
1421d9b32SPeter Brune /* Defines the basic SNES object */
2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I  "petscsnes.h"  I*/
3421d9b32SPeter Brune 
49e5d0892SLisandro Dalcin const char *const SNESFASTypes[] = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "SNESFASType", "SNES_FAS", NULL};
507144faaSPeter Brune 
6d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_FAS(SNES snes)
7d71ae5a4SJacob Faibussowitsch {
8421d9b32SPeter Brune   SNES_FAS *fas = (SNES_FAS *)snes->data;
9421d9b32SPeter Brune 
10421d9b32SPeter Brune   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&fas->smoothu));
129566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&fas->smoothd));
139566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->inject));
149566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->interpolate));
159566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->restrct));
169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fas->rscale));
179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fas->Xg));
189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fas->Fg));
199566063dSJacob Faibussowitsch   if (fas->next) PetscCall(SNESReset(fas->next));
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21421d9b32SPeter Brune }
22421d9b32SPeter Brune 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_FAS(SNES snes)
24d71ae5a4SJacob Faibussowitsch {
25421d9b32SPeter Brune   SNES_FAS *fas = (SNES_FAS *)snes->data;
26421d9b32SPeter Brune 
27421d9b32SPeter Brune   PetscFunctionBegin;
28421d9b32SPeter Brune   /* recursively resets and then destroys */
299566063dSJacob Faibussowitsch   PetscCall(SNESReset_FAS(snes));
309566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&fas->next));
319566063dSJacob Faibussowitsch   PetscCall(PetscFree(fas));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33421d9b32SPeter Brune }
34421d9b32SPeter Brune 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth)
36d71ae5a4SJacob Faibussowitsch {
37f833ba53SLisandro Dalcin   SNESLineSearch linesearch;
38f833ba53SLisandro Dalcin   SNESLineSearch slinesearch;
39f833ba53SLisandro Dalcin   void          *lsprectx, *lspostctx;
40f833ba53SLisandro Dalcin   PetscErrorCode (*precheck)(SNESLineSearch, Vec, Vec, PetscBool *, void *);
41f833ba53SLisandro Dalcin   PetscErrorCode (*postcheck)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *);
42f833ba53SLisandro Dalcin 
43f833ba53SLisandro Dalcin   PetscFunctionBegin;
443ba16761SJacob Faibussowitsch   if (!snes->linesearch) PetscFunctionReturn(PETSC_SUCCESS);
459566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
469566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(smooth, &slinesearch));
479566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetPreCheck(linesearch, &precheck, &lsprectx));
489566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetPostCheck(linesearch, &postcheck, &lspostctx));
499566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetPreCheck(slinesearch, precheck, lsprectx));
509566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetPostCheck(slinesearch, postcheck, lspostctx));
519566063dSJacob Faibussowitsch   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch));
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53f833ba53SLisandro Dalcin }
54f833ba53SLisandro Dalcin 
55d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth)
56d71ae5a4SJacob Faibussowitsch {
57f833ba53SLisandro Dalcin   SNES_FAS *fas = (SNES_FAS *)snes->data;
58f833ba53SLisandro Dalcin 
59f833ba53SLisandro Dalcin   PetscFunctionBegin;
609566063dSJacob Faibussowitsch   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth));
619566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(smooth));
629566063dSJacob Faibussowitsch   PetscCall(SNESFASSetUpLineSearch_Private(snes, smooth));
63f833ba53SLisandro Dalcin 
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)snes->vec_sol));
659566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)snes->vec_sol_update));
669566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)snes->vec_func));
67f833ba53SLisandro Dalcin   smooth->vec_sol        = snes->vec_sol;
68f833ba53SLisandro Dalcin   smooth->vec_sol_update = snes->vec_sol_update;
69f833ba53SLisandro Dalcin   smooth->vec_func       = snes->vec_func;
70f833ba53SLisandro Dalcin 
719566063dSJacob Faibussowitsch   if (fas->eventsmoothsetup) PetscCall(PetscLogEventBegin(fas->eventsmoothsetup, smooth, 0, 0, 0));
729566063dSJacob Faibussowitsch   PetscCall(SNESSetUp(smooth));
739566063dSJacob Faibussowitsch   if (fas->eventsmoothsetup) PetscCall(PetscLogEventEnd(fas->eventsmoothsetup, smooth, 0, 0, 0));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75f833ba53SLisandro Dalcin }
76f833ba53SLisandro Dalcin 
77d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_FAS(SNES snes)
78d71ae5a4SJacob Faibussowitsch {
7948bfdf8aSPeter Brune   SNES_FAS *fas = (SNES_FAS *)snes->data;
80d1adcc6fSPeter Brune   PetscInt  dm_levels;
81ab8d36c9SPeter Brune   SNES      next;
82f833ba53SLisandro Dalcin   PetscBool isFine, hasCreateRestriction, hasCreateInjection;
83eff52c0eSPeter Brune 
846b2b7091SBarry Smith   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
86ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
879566063dSJacob Faibussowitsch     PetscCall(DMGetRefineLevel(snes->dm, &dm_levels));
88d1adcc6fSPeter Brune     dm_levels++;
89cc05f883SPeter Brune     if (dm_levels > fas->levels) {
903dccd265SPeter Brune       /* reset the number of levels */
919566063dSJacob Faibussowitsch       PetscCall(SNESFASSetLevels(snes, dm_levels, NULL));
929566063dSJacob Faibussowitsch       PetscCall(SNESSetFromOptions(snes));
93d1adcc6fSPeter Brune     }
94d1adcc6fSPeter Brune   }
959566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
96ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
973dccd265SPeter Brune 
989566063dSJacob Faibussowitsch   PetscCall(SNESSetWorkVecs(snes, 2)); /* work vectors used for intergrid transfers */
99cc05f883SPeter Brune 
100ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
10148a46eb9SPierre Jolivet   if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
102ab8d36c9SPeter Brune 
10379d9a41aSPeter Brune   if (snes->dm) {
104ab8d36c9SPeter Brune     /* set the smoother DMs properly */
1059566063dSJacob Faibussowitsch     if (fas->smoothu) PetscCall(SNESSetDM(fas->smoothu, snes->dm));
1069566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(fas->smoothd, snes->dm));
10779d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
108ab8d36c9SPeter Brune     if (next) {
10979d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
110ab8d36c9SPeter Brune       if (!next->dm) {
1119566063dSJacob Faibussowitsch         PetscCall(DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm));
1129566063dSJacob Faibussowitsch         PetscCall(SNESSetDM(next, next->dm));
11379d9a41aSPeter Brune       }
11479d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
11579d9a41aSPeter Brune       if (!fas->interpolate) {
1169566063dSJacob Faibussowitsch         PetscCall(DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale));
117bccf9bb3SJed Brown         if (!fas->restrct) {
1189566063dSJacob Faibussowitsch           PetscCall(DMHasCreateRestriction(next->dm, &hasCreateRestriction));
1190a7266b2SPatrick Farrell           /* DM can create restrictions, use that */
1200a7266b2SPatrick Farrell           if (hasCreateRestriction) {
1219566063dSJacob Faibussowitsch             PetscCall(DMCreateRestriction(next->dm, snes->dm, &fas->restrct));
1220a7266b2SPatrick Farrell           } else {
1239566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)fas->interpolate));
12479d9a41aSPeter Brune             fas->restrct = fas->interpolate;
12579d9a41aSPeter Brune           }
126bccf9bb3SJed Brown         }
1270a7266b2SPatrick Farrell       }
12879d9a41aSPeter Brune       /* set the injection from the DM */
12979d9a41aSPeter Brune       if (!fas->inject) {
1309566063dSJacob Faibussowitsch         PetscCall(DMHasCreateInjection(next->dm, &hasCreateInjection));
13148a46eb9SPierre Jolivet         if (hasCreateInjection) PetscCall(DMCreateInjection(next->dm, snes->dm, &fas->inject));
13279d9a41aSPeter Brune       }
13379d9a41aSPeter Brune     }
13423e68893SLawrence Mitchell   }
135f833ba53SLisandro Dalcin 
13679d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
13779d9a41aSPeter Brune   if (fas->galerkin) {
1381baa6e33SBarry Smith     if (next) PetscCall(SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next));
13948a46eb9SPierre Jolivet     if (fas->smoothd && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes));
14048a46eb9SPierre Jolivet     if (fas->smoothu && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes));
14179d9a41aSPeter Brune   }
14279d9a41aSPeter Brune 
143534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
144534ebe21SPeter Brune   if (fas->smoothd) {
145*2d157150SStefano Zampini     if (fas->level == 0) {
146*2d157150SStefano Zampini       PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_ALWAYS));
147534ebe21SPeter Brune     } else {
1489566063dSJacob Faibussowitsch       PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY));
149534ebe21SPeter Brune     }
1509566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd));
151534ebe21SPeter Brune   }
152534ebe21SPeter Brune 
153534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
154534ebe21SPeter Brune   if (fas->smoothu) {
155534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
1569566063dSJacob Faibussowitsch       PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE));
157534ebe21SPeter Brune     } else {
1589566063dSJacob Faibussowitsch       PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY));
159534ebe21SPeter Brune     }
1609566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu));
161534ebe21SPeter Brune   }
162d06165b7SPeter Brune 
163ab8d36c9SPeter Brune   if (next) {
16479d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
1659566063dSJacob Faibussowitsch     if (!next->vec_sol) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_sol));
1669566063dSJacob Faibussowitsch     if (!next->vec_rhs) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_rhs));
1679566063dSJacob Faibussowitsch     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next));
1689566063dSJacob Faibussowitsch     PetscCall(SNESFASSetUpLineSearch_Private(snes, next));
1699566063dSJacob Faibussowitsch     PetscCall(SNESSetUp(next));
17079d9a41aSPeter Brune   }
171f833ba53SLisandro Dalcin 
1726273346dSPeter Brune   /* setup FAS work vectors */
1736273346dSPeter Brune   if (fas->galerkin) {
1749566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(snes->vec_sol, &fas->Xg));
1759566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(snes->vec_sol, &fas->Fg));
1766273346dSPeter Brune   }
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178421d9b32SPeter Brune }
179421d9b32SPeter Brune 
180d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_FAS(SNES snes, PetscOptionItems *PetscOptionsObject)
181d71ae5a4SJacob Faibussowitsch {
182ee78dd50SPeter Brune   SNES_FAS      *fas    = (SNES_FAS *)snes->data;
183ee78dd50SPeter Brune   PetscInt       levels = 1;
18487f44e3fSPeter Brune   PetscBool      flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE, continuationflg = PETSC_FALSE;
18507144faaSPeter Brune   SNESFASType    fastype;
186fde0ff24SPeter Brune   const char    *optionsprefix;
187f1c6b773SPeter Brune   SNESLineSearch linesearch;
18866585501SPeter Brune   PetscInt       m, n_up, n_down;
189ab8d36c9SPeter Brune   SNES           next;
190ab8d36c9SPeter Brune   PetscBool      isFine;
191421d9b32SPeter Brune 
192421d9b32SPeter Brune   PetscFunctionBegin;
1939566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
194d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNESFAS Options-----------------------------------");
195ee78dd50SPeter Brune 
196ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
197ab8d36c9SPeter Brune   if (isFine) {
1989566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg));
199c732cbdbSBarry Smith     if (!flg && snes->dm) {
2009566063dSJacob Faibussowitsch       PetscCall(DMGetRefineLevel(snes->dm, &levels));
201c732cbdbSBarry Smith       levels++;
202d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
203c732cbdbSBarry Smith     }
2049566063dSJacob Faibussowitsch     PetscCall(SNESFASSetLevels(snes, levels, NULL));
20507144faaSPeter Brune     fastype = fas->fastype;
2069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-snes_fas_type", "FAS correction type", "SNESFASSetType", SNESFASTypes, (PetscEnum)fastype, (PetscEnum *)&fastype, &flg));
2071baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetType(snes, fastype));
208ee78dd50SPeter Brune 
2099566063dSJacob Faibussowitsch     PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
2109566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_cycles", "Number of cycles", "SNESFASSetCycles", fas->n_cycles, &m, &flg));
2111baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetCycles(snes, m));
2129566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-snes_fas_continuation", "Corrected grid-sequence continuation", "SNESFASSetContinuation", fas->continuation, &continuationflg, &flg));
2131baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetContinuation(snes, continuationflg));
214fde0ff24SPeter Brune 
2159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin", "SNESFASSetGalerkin", fas->galerkin, &galerkinflg, &flg));
2161baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetGalerkin(snes, galerkinflg));
217ee78dd50SPeter Brune 
218928e959bSPeter Brune     if (fas->fastype == SNES_FAS_FULL) {
2199566063dSJacob Faibussowitsch       PetscCall(PetscOptionsBool("-snes_fas_full_downsweep", "Smooth on the initial down sweep for full FAS cycles", "SNESFASFullSetDownSweep", fas->full_downsweep, &fas->full_downsweep, &flg));
2209566063dSJacob Faibussowitsch       if (flg) PetscCall(SNESFASFullSetDownSweep(snes, fas->full_downsweep));
2219566063dSJacob Faibussowitsch       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));
2229566063dSJacob Faibussowitsch       if (flg) PetscCall(SNESFASFullSetTotal(snes, fas->full_total));
223928e959bSPeter Brune     }
224928e959bSPeter Brune 
2259566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_smoothup", "Number of post-smoothing steps", "SNESFASSetNumberSmoothUp", fas->max_up_it, &n_up, &upflg));
226162d76ddSPeter Brune 
2279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_smoothdown", "Number of pre-smoothing steps", "SNESFASSetNumberSmoothDown", fas->max_down_it, &n_down, &downflg));
228162d76ddSPeter Brune 
229d142ab34SLawrence Mitchell     {
230d142ab34SLawrence Mitchell       PetscViewer       viewer;
231d142ab34SLawrence Mitchell       PetscViewerFormat format;
2329566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_fas_monitor", &viewer, &format, &monflg));
233d142ab34SLawrence Mitchell       if (monflg) {
234d142ab34SLawrence Mitchell         PetscViewerAndFormat *vf;
2359566063dSJacob Faibussowitsch         PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf));
2369566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)viewer));
2379566063dSJacob Faibussowitsch         PetscCall(SNESFASSetMonitor(snes, vf, PETSC_TRUE));
238d142ab34SLawrence Mitchell       }
239d142ab34SLawrence Mitchell     }
2400dd27c6cSPeter Brune     flg    = PETSC_FALSE;
2410dd27c6cSPeter Brune     monflg = PETSC_TRUE;
2429566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-snes_fas_log", "Log times for each FAS level", "SNESFASSetLog", monflg, &monflg, &flg));
2439566063dSJacob Faibussowitsch     if (flg) PetscCall(SNESFASSetLog(snes, monflg));
244ab8d36c9SPeter Brune   }
245ee78dd50SPeter Brune 
246d0609cedSBarry Smith   PetscOptionsHeadEnd();
247f833ba53SLisandro Dalcin 
2488cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
2491baa6e33SBarry Smith   if (upflg) PetscCall(SNESFASSetNumberSmoothUp(snes, n_up));
2501baa6e33SBarry Smith   if (downflg) PetscCall(SNESFASSetNumberSmoothDown(snes, n_down));
251eff52c0eSPeter Brune 
2529e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
2539e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
2549e764e56SPeter Brune     if (!snes->linesearch) {
2559566063dSJacob Faibussowitsch       PetscCall(SNESGetLineSearch(snes, &linesearch));
2569566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
2579e764e56SPeter Brune     }
2589e764e56SPeter Brune   }
2599e764e56SPeter Brune 
260ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
2619566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
2629566063dSJacob Faibussowitsch   if (next) PetscCall(SNESSetFromOptions(next));
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
264421d9b32SPeter Brune }
265421d9b32SPeter Brune 
2669804daf3SBarry Smith #include <petscdraw.h>
267d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
268d71ae5a4SJacob Faibussowitsch {
269421d9b32SPeter Brune   SNES_FAS *fas = (SNES_FAS *)snes->data;
270656ede7eSPeter Brune   PetscBool isFine, iascii, isdraw;
271ab8d36c9SPeter Brune   PetscInt  i;
272ab8d36c9SPeter Brune   SNES      smoothu, smoothd, levelsnes;
273421d9b32SPeter Brune 
274421d9b32SPeter Brune   PetscFunctionBegin;
2759566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
276ab8d36c9SPeter Brune   if (isFine) {
2779566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2789566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
279421d9b32SPeter Brune     if (iascii) {
28063a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%" PetscInt_FMT ", cycles=%" PetscInt_FMT "\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles));
281ab8d36c9SPeter Brune       if (fas->galerkin) {
2829566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Using Galerkin computed coarse grid function evaluation\n"));
283421d9b32SPeter Brune       } else {
2849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Not using Galerkin computed coarse grid function evaluation\n"));
285421d9b32SPeter Brune       }
286ab8d36c9SPeter Brune       for (i = 0; i < fas->levels; i++) {
2879566063dSJacob Faibussowitsch         PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
2889566063dSJacob Faibussowitsch         PetscCall(SNESFASCycleGetSmootherUp(levelsnes, &smoothu));
2899566063dSJacob Faibussowitsch         PetscCall(SNESFASCycleGetSmootherDown(levelsnes, &smoothd));
290ab8d36c9SPeter Brune         if (!i) {
29163a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i));
292421d9b32SPeter Brune         } else {
29363a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
294421d9b32SPeter Brune         }
2959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
296166b3ea4SJed Brown         if (smoothd) {
2979566063dSJacob Faibussowitsch           PetscCall(SNESView(smoothd, viewer));
298166b3ea4SJed Brown         } else {
2999566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n"));
300166b3ea4SJed Brown         }
3019566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
302ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
3039566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  Up solver (post-smoother) same as down solver (pre-smoother)\n"));
304ab8d36c9SPeter Brune         } else if (i) {
30563a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "  Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
3069566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
307166b3ea4SJed Brown           if (smoothu) {
3089566063dSJacob Faibussowitsch             PetscCall(SNESView(smoothu, viewer));
309166b3ea4SJed Brown           } else {
3109566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n"));
311166b3ea4SJed Brown           }
3129566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
313ab8d36c9SPeter Brune         }
314ab8d36c9SPeter Brune       }
315656ede7eSPeter Brune     } else if (isdraw) {
316656ede7eSPeter Brune       PetscDraw draw;
317b4375e8dSPeter Brune       PetscReal x, w, y, bottom, th, wth;
318656ede7eSPeter Brune       SNES_FAS *curfas = fas;
3199566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
3209566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
3219566063dSJacob Faibussowitsch       PetscCall(PetscDrawStringGetSize(draw, &wth, &th));
322656ede7eSPeter Brune       bottom = y - th;
323656ede7eSPeter Brune       while (curfas) {
324b4375e8dSPeter Brune         if (!curfas->smoothu) {
3259566063dSJacob Faibussowitsch           PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
3269566063dSJacob Faibussowitsch           if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer));
3279566063dSJacob Faibussowitsch           PetscCall(PetscDrawPopCurrentPoint(draw));
328b4375e8dSPeter Brune         } else {
329b4375e8dSPeter Brune           w = 0.5 * PetscMin(1.0 - x, x);
3309566063dSJacob Faibussowitsch           PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom));
3319566063dSJacob Faibussowitsch           if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer));
3329566063dSJacob Faibussowitsch           PetscCall(PetscDrawPopCurrentPoint(draw));
3339566063dSJacob Faibussowitsch           PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom));
3349566063dSJacob Faibussowitsch           if (curfas->smoothu) PetscCall(SNESView(curfas->smoothu, viewer));
3359566063dSJacob Faibussowitsch           PetscCall(PetscDrawPopCurrentPoint(draw));
336b4375e8dSPeter Brune         }
337656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
338656ede7eSPeter Brune         bottom -= 5 * th;
3391aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS *)curfas->next->data;
3400298fd71SBarry Smith         else curfas = NULL;
341656ede7eSPeter Brune       }
342421d9b32SPeter Brune     }
343ab8d36c9SPeter Brune   }
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
345421d9b32SPeter Brune }
346421d9b32SPeter Brune 
34739bd7f45SPeter Brune /*
34839bd7f45SPeter Brune Defines the action of the downsmoother
34939bd7f45SPeter Brune  */
350d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
351d71ae5a4SJacob Faibussowitsch {
352742fe5e2SPeter Brune   SNESConvergedReason reason;
353ab8d36c9SPeter Brune   Vec                 FPC;
354ab8d36c9SPeter Brune   SNES                smoothd;
3556cbb2f26SLawrence Mitchell   PetscBool           flg;
3560dd27c6cSPeter Brune   SNES_FAS           *fas = (SNES_FAS *)snes->data;
3576e111a19SKarl Rupp 
358421d9b32SPeter Brune   PetscFunctionBegin;
3599566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetSmootherDown(snes, &smoothd));
3609566063dSJacob Faibussowitsch   PetscCall(SNESSetInitialFunction(smoothd, F));
3619566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothd, B, X, 0));
3629566063dSJacob Faibussowitsch   PetscCall(SNESSolve(smoothd, B, X));
3639566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothd, B, X, 0));
364742fe5e2SPeter Brune   /* check convergence reason for the smoother */
3659566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(smoothd, &reason));
3661ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
367742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
3683ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
369742fe5e2SPeter Brune   }
3706cbb2f26SLawrence Mitchell 
3719566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(smoothd, &FPC, NULL, NULL));
3729566063dSJacob Faibussowitsch   PetscCall(SNESGetAlwaysComputesFinalResidual(smoothd, &flg));
37348a46eb9SPierre Jolivet   if (!flg) PetscCall(SNESComputeFunction(smoothd, X, FPC));
3749566063dSJacob Faibussowitsch   PetscCall(VecCopy(FPC, F));
3759566063dSJacob Faibussowitsch   if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm));
3763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37739bd7f45SPeter Brune }
37839bd7f45SPeter Brune 
37939bd7f45SPeter Brune /*
38007144faaSPeter Brune Defines the action of the upsmoother
38139bd7f45SPeter Brune  */
382d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
383d71ae5a4SJacob Faibussowitsch {
38439bd7f45SPeter Brune   SNESConvergedReason reason;
385ab8d36c9SPeter Brune   Vec                 FPC;
386ab8d36c9SPeter Brune   SNES                smoothu;
3876cbb2f26SLawrence Mitchell   PetscBool           flg;
3880dd27c6cSPeter Brune   SNES_FAS           *fas = (SNES_FAS *)snes->data;
389ab8d36c9SPeter Brune 
3906e111a19SKarl Rupp   PetscFunctionBegin;
3919566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetSmootherUp(snes, &smoothu));
3929566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothu, 0, 0, 0));
3939566063dSJacob Faibussowitsch   PetscCall(SNESSolve(smoothu, B, X));
3949566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothu, 0, 0, 0));
39539bd7f45SPeter Brune   /* check convergence reason for the smoother */
3969566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(smoothu, &reason));
3971ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
39839bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
3993ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
40039bd7f45SPeter Brune   }
4019566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(smoothu, &FPC, NULL, NULL));
4029566063dSJacob Faibussowitsch   PetscCall(SNESGetAlwaysComputesFinalResidual(smoothu, &flg));
40348a46eb9SPierre Jolivet   if (!flg) PetscCall(SNESComputeFunction(smoothu, X, FPC));
4049566063dSJacob Faibussowitsch   PetscCall(VecCopy(FPC, F));
4059566063dSJacob Faibussowitsch   if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm));
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40739bd7f45SPeter Brune }
40839bd7f45SPeter Brune 
409938e4a01SJed Brown /*@
410f6dfbefdSBarry Smith    SNESFASCreateCoarseVec - create `Vec` corresponding to a state vector on one level coarser than current level
411938e4a01SJed Brown 
412938e4a01SJed Brown    Collective
413938e4a01SJed Brown 
4144165533cSJose E. Roman    Input Parameter:
415f6dfbefdSBarry Smith .  snes - `SNESFAS` object
416938e4a01SJed Brown 
4174165533cSJose E. Roman    Output Parameter:
418938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
419938e4a01SJed Brown 
420938e4a01SJed Brown    Level: developer
421938e4a01SJed Brown 
422db781477SPatrick Sanan .seealso: `SNESFASSetRestriction()`, `SNESFASRestrict()`
423938e4a01SJed Brown @*/
424d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCreateCoarseVec(SNES snes, Vec *Xcoarse)
425d71ae5a4SJacob Faibussowitsch {
426f833ba53SLisandro Dalcin   SNES_FAS *fas;
427938e4a01SJed Brown 
428938e4a01SJed Brown   PetscFunctionBegin;
429f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
430f833ba53SLisandro Dalcin   PetscValidPointer(Xcoarse, 2);
431f833ba53SLisandro Dalcin   fas = (SNES_FAS *)snes->data;
4321aa26658SKarl Rupp   if (fas->rscale) {
4339566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(fas->rscale, Xcoarse));
43490b9e4c9SPablo Brubeck   } else if (fas->interpolate) {
43590b9e4c9SPablo Brubeck     PetscCall(MatCreateVecs(fas->interpolate, Xcoarse, NULL));
43690b9e4c9SPablo Brubeck   } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set rscale or interpolation");
4373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
438938e4a01SJed Brown }
439938e4a01SJed Brown 
440e9923e8dSJed Brown /*@
441f6dfbefdSBarry Smith    SNESFASRestrict - restrict a `Vec` to the next coarser level
442e9923e8dSJed Brown 
443e9923e8dSJed Brown    Collective
444e9923e8dSJed Brown 
4454165533cSJose E. Roman    Input Parameters:
446f6dfbefdSBarry Smith +  fine - `SNES` from which to restrict
447e9923e8dSJed Brown -  Xfine - vector to restrict
448e9923e8dSJed Brown 
4494165533cSJose E. Roman    Output Parameter:
450e9923e8dSJed Brown .  Xcoarse - result of restriction
451e9923e8dSJed Brown 
452e9923e8dSJed Brown    Level: developer
453e9923e8dSJed Brown 
454f6dfbefdSBarry Smith .seealso: `SNES`, `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASSetInjection()`
455e9923e8dSJed Brown @*/
456d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASRestrict(SNES fine, Vec Xfine, Vec Xcoarse)
457d71ae5a4SJacob Faibussowitsch {
458f833ba53SLisandro Dalcin   SNES_FAS *fas;
459e9923e8dSJed Brown 
460e9923e8dSJed Brown   PetscFunctionBegin;
461f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(fine, SNES_CLASSID, 1, SNESFAS);
462e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine, VEC_CLASSID, 2);
463e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse, VEC_CLASSID, 3);
464f833ba53SLisandro Dalcin   fas = (SNES_FAS *)fine->data;
465e9923e8dSJed Brown   if (fas->inject) {
4669566063dSJacob Faibussowitsch     PetscCall(MatRestrict(fas->inject, Xfine, Xcoarse));
467e9923e8dSJed Brown   } else {
4689566063dSJacob Faibussowitsch     PetscCall(MatRestrict(fas->restrct, Xfine, Xcoarse));
4699566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Xcoarse, fas->rscale, Xcoarse));
470e9923e8dSJed Brown   }
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
472e9923e8dSJed Brown }
473e9923e8dSJed Brown 
47439bd7f45SPeter Brune /*
47539bd7f45SPeter Brune 
4766dfbafefSToby Isaac Performs a variant of FAS using the interpolated total coarse solution
4776dfbafefSToby Isaac 
4786dfbafefSToby Isaac fine problem:   F(x) = b
4796dfbafefSToby Isaac coarse problem: F^c(x^c) = Rb, Initial guess Rx
4806dfbafefSToby Isaac interpolated solution: x^f = I x^c (total solution interpolation
4816dfbafefSToby Isaac 
4826dfbafefSToby Isaac  */
483d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new)
484d71ae5a4SJacob Faibussowitsch {
4856dfbafefSToby Isaac   Vec                 X_c, B_c;
4866dfbafefSToby Isaac   SNESConvergedReason reason;
4876dfbafefSToby Isaac   SNES                next;
4886dfbafefSToby Isaac   Mat                 restrct, interpolate;
4896dfbafefSToby Isaac   SNES_FAS           *fasc;
4906dfbafefSToby Isaac 
4916dfbafefSToby Isaac   PetscFunctionBegin;
4929566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
4936dfbafefSToby Isaac   if (next) {
4946dfbafefSToby Isaac     fasc = (SNES_FAS *)next->data;
4956dfbafefSToby Isaac 
4969566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
4979566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
4986dfbafefSToby Isaac 
4996dfbafefSToby Isaac     X_c = next->vec_sol;
5006dfbafefSToby Isaac 
5019566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
5026dfbafefSToby Isaac     /* restrict the total solution: Rb */
5039566063dSJacob Faibussowitsch     PetscCall(SNESFASRestrict(snes, X, X_c));
5046dfbafefSToby Isaac     B_c = next->vec_rhs;
5056dfbafefSToby Isaac     if (snes->vec_rhs) {
5066dfbafefSToby Isaac       /* restrict the total rhs defect: Rb */
5079566063dSJacob Faibussowitsch       PetscCall(MatRestrict(restrct, snes->vec_rhs, B_c));
5086dfbafefSToby Isaac     } else {
5099566063dSJacob Faibussowitsch       PetscCall(VecSet(B_c, 0.));
5106dfbafefSToby Isaac     }
5119566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
5126dfbafefSToby Isaac 
5139566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next, B_c, X_c));
5149566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next, &reason));
5156dfbafefSToby Isaac     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
5166dfbafefSToby Isaac       snes->reason = SNES_DIVERGED_INNER;
5173ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
5186dfbafefSToby Isaac     }
5196dfbafefSToby Isaac     /* x^f <- Ix^c*/
5206dfbafefSToby Isaac     DM dmc, dmf;
5216dfbafefSToby Isaac 
5229566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(next, &dmc));
5239566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dmf));
5249566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
5259566063dSJacob Faibussowitsch     PetscCall(DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new));
5269566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
5279566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse solution"));
5289566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view"));
5299566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution"));
5309566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view"));
5316dfbafefSToby Isaac   }
5323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5336dfbafefSToby Isaac }
5346dfbafefSToby Isaac 
5356dfbafefSToby Isaac /*
5366dfbafefSToby Isaac 
53739bd7f45SPeter Brune Performs the FAS coarse correction as:
53839bd7f45SPeter Brune 
539b5527d98SMatthew G. Knepley fine problem:   F(x) = b
540b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c
54139bd7f45SPeter Brune 
542b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
54339bd7f45SPeter Brune 
54439bd7f45SPeter Brune  */
545d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
546d71ae5a4SJacob Faibussowitsch {
54739bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
54839bd7f45SPeter Brune   SNESConvergedReason reason;
549ab8d36c9SPeter Brune   SNES                next;
550ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
5510dd27c6cSPeter Brune   SNES_FAS           *fasc;
5525fd66863SKarl Rupp 
55339bd7f45SPeter Brune   PetscFunctionBegin;
5549566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
555ab8d36c9SPeter Brune   if (next) {
5560dd27c6cSPeter Brune     fasc = (SNES_FAS *)next->data;
5570dd27c6cSPeter Brune 
5589566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
5599566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
560ab8d36c9SPeter Brune 
561ab8d36c9SPeter Brune     X_c  = next->vec_sol;
562ab8d36c9SPeter Brune     Xo_c = next->work[0];
563ab8d36c9SPeter Brune     F_c  = next->vec_func;
564ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
565efe1f98aSPeter Brune 
5669566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
5679566063dSJacob Faibussowitsch     PetscCall(SNESFASRestrict(snes, X, Xo_c));
5685609cbf2SMatthew G. Knepley     /* restrict the defect: R(F(x) - b) */
5699566063dSJacob Faibussowitsch     PetscCall(MatRestrict(restrct, F, B_c));
5709566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
5710dd27c6cSPeter Brune 
5729566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0));
5735609cbf2SMatthew G. Knepley     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
5749566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(next, Xo_c, F_c));
5759566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0));
5760dd27c6cSPeter Brune 
5770dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
5789566063dSJacob Faibussowitsch     PetscCall(VecCopy(B_c, X_c));
5799566063dSJacob Faibussowitsch     PetscCall(VecCopy(F_c, B_c));
5809566063dSJacob Faibussowitsch     PetscCall(VecCopy(X_c, F_c));
581ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
5829566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xo_c, X_c));
583c90fad12SPeter Brune 
584c90fad12SPeter Brune     /* recurse to the next level */
5859566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next, F_c));
5869566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next, B_c, X_c));
5879566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next, &reason));
588742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
589742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
5903ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
591742fe5e2SPeter Brune     }
592fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
5939566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X_c, -1.0, Xo_c));
5940dd27c6cSPeter Brune 
5959566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0));
5969566063dSJacob Faibussowitsch     PetscCall(MatInterpolateAdd(interpolate, X_c, X, X_new));
5979566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0));
5989566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse correction"));
5999566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view"));
6009566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution"));
6019566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view"));
602293a7e31SPeter Brune   }
6033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60439bd7f45SPeter Brune }
60539bd7f45SPeter Brune 
60639bd7f45SPeter Brune /*
60739bd7f45SPeter Brune 
60839bd7f45SPeter Brune The additive cycle looks like:
60939bd7f45SPeter Brune 
61007144faaSPeter Brune xhat = x
61107144faaSPeter Brune xhat = dS(x, b)
61207144faaSPeter Brune x = coarsecorrection(xhat, b_d)
61307144faaSPeter Brune x = x + nu*(xhat - x);
61439bd7f45SPeter Brune (optional) x = uS(x, b)
61539bd7f45SPeter Brune 
61639bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
61739bd7f45SPeter Brune 
61839bd7f45SPeter Brune  */
619d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
620d71ae5a4SJacob Faibussowitsch {
62107144faaSPeter Brune   Vec                  F, B, Xhat;
62222c1e704SPeter Brune   Vec                  X_c, Xo_c, F_c, B_c;
62307144faaSPeter Brune   SNESConvergedReason  reason;
62422c1e704SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
625422a814eSBarry Smith   SNESLineSearchReason lsresult;
626ab8d36c9SPeter Brune   SNES                 next;
627ab8d36c9SPeter Brune   Mat                  restrct, interpolate;
6280dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data, *fasc;
6290dd27c6cSPeter Brune 
63039bd7f45SPeter Brune   PetscFunctionBegin;
6319566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
63239bd7f45SPeter Brune   F    = snes->vec_func;
63339bd7f45SPeter Brune   B    = snes->vec_rhs;
634e7f468e7SPeter Brune   Xhat = snes->work[1];
6359566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Xhat));
63607144faaSPeter Brune   /* recurse first */
637ab8d36c9SPeter Brune   if (next) {
6380dd27c6cSPeter Brune     fasc = (SNES_FAS *)next->data;
6399566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
6409566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
6419566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0));
6429566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, Xhat, F));
6439566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0));
6449566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
645ab8d36c9SPeter Brune     X_c  = next->vec_sol;
646ab8d36c9SPeter Brune     Xo_c = next->work[0];
647ab8d36c9SPeter Brune     F_c  = next->vec_func;
648ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
64939bd7f45SPeter Brune 
6509566063dSJacob Faibussowitsch     PetscCall(SNESFASRestrict(snes, Xhat, Xo_c));
65107144faaSPeter Brune     /* restrict the defect */
6529566063dSJacob Faibussowitsch     PetscCall(MatRestrict(restrct, F, B_c));
65307144faaSPeter Brune 
65407144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
6559566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0));
6569566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(next, Xo_c, F_c));
6579566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0));
6589566063dSJacob Faibussowitsch     PetscCall(VecCopy(B_c, X_c));
6599566063dSJacob Faibussowitsch     PetscCall(VecCopy(F_c, B_c));
6609566063dSJacob Faibussowitsch     PetscCall(VecCopy(X_c, F_c));
66107144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
6629566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xo_c, X_c));
66307144faaSPeter Brune 
66407144faaSPeter Brune     /* recurse */
6659566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next, F_c));
6669566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next, B_c, X_c));
66707144faaSPeter Brune 
66807144faaSPeter Brune     /* smooth on this level */
6699566063dSJacob Faibussowitsch     PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &fnorm));
67007144faaSPeter Brune 
6719566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next, &reason));
67207144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
67307144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
6743ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
67507144faaSPeter Brune     }
67607144faaSPeter Brune 
67707144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
6789566063dSJacob Faibussowitsch     PetscCall(VecAYPX(X_c, -1.0, Xo_c));
6799566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(interpolate, X_c, Xhat));
68007144faaSPeter Brune 
681ddebd997SPeter Brune     /* additive correction of the coarse direction*/
6829566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat));
6839566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
6849566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm));
685422a814eSBarry Smith     if (lsresult) {
6869e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6879e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6883ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
6899e764e56SPeter Brune       }
6909e764e56SPeter Brune     }
69107144faaSPeter Brune   } else {
6929566063dSJacob Faibussowitsch     PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
69307144faaSPeter Brune   }
6943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69539bd7f45SPeter Brune }
69639bd7f45SPeter Brune 
69739bd7f45SPeter Brune /*
69839bd7f45SPeter Brune 
69939bd7f45SPeter Brune Defines the FAS cycle as:
70039bd7f45SPeter Brune 
701b5527d98SMatthew G. Knepley fine problem: F(x) = b
70239bd7f45SPeter Brune coarse problem: F^c(x) = b^c
70339bd7f45SPeter Brune 
704b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
70539bd7f45SPeter Brune 
70639bd7f45SPeter Brune correction:
70739bd7f45SPeter Brune 
70839bd7f45SPeter Brune x = x + I(x^c - Rx)
70939bd7f45SPeter Brune 
71039bd7f45SPeter Brune  */
711d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
712d71ae5a4SJacob Faibussowitsch {
71339bd7f45SPeter Brune   Vec  F, B;
71434d65b3cSPeter Brune   SNES next;
71539bd7f45SPeter Brune 
71639bd7f45SPeter Brune   PetscFunctionBegin;
71739bd7f45SPeter Brune   F = snes->vec_func;
71839bd7f45SPeter Brune   B = snes->vec_rhs;
71939bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7209566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
7219566063dSJacob Faibussowitsch   PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
72234d65b3cSPeter Brune   if (next) {
7239566063dSJacob Faibussowitsch     PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
7249566063dSJacob Faibussowitsch     PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
725fe6f9142SPeter Brune   }
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
727421d9b32SPeter Brune }
728421d9b32SPeter Brune 
729d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
730d71ae5a4SJacob Faibussowitsch {
7318c94862eSPeter Brune   SNES      next;
7328c94862eSPeter Brune   SNES_FAS *fas = (SNES_FAS *)snes->data;
7338c94862eSPeter Brune   PetscBool isFine;
7348c94862eSPeter Brune 
7358c94862eSPeter Brune   PetscFunctionBegin;
7368c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7379566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
7389566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
7398c94862eSPeter Brune   fas->full_stage = 0;
7409566063dSJacob Faibussowitsch   if (next) PetscCall(SNESFASCycleSetupPhase_Full(next));
7413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7428c94862eSPeter Brune }
7438c94862eSPeter Brune 
744d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
745d71ae5a4SJacob Faibussowitsch {
7468c94862eSPeter Brune   Vec       F, B;
7478c94862eSPeter Brune   SNES_FAS *fas = (SNES_FAS *)snes->data;
7488c94862eSPeter Brune   PetscBool isFine;
7498c94862eSPeter Brune   SNES      next;
7508c94862eSPeter Brune 
7518c94862eSPeter Brune   PetscFunctionBegin;
7528c94862eSPeter Brune   F = snes->vec_func;
7538c94862eSPeter Brune   B = snes->vec_rhs;
7549566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
7559566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
7568c94862eSPeter Brune 
7571baa6e33SBarry Smith   if (isFine) PetscCall(SNESFASCycleSetupPhase_Full(snes));
7588c94862eSPeter Brune 
7598c94862eSPeter Brune   if (fas->full_stage == 0) {
760928e959bSPeter Brune     /* downsweep */
7618c94862eSPeter Brune     if (next) {
7628c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
7639566063dSJacob Faibussowitsch       if (fas->full_downsweep) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
764e140b65aSPatrick Farrell       fas->full_downsweep = PETSC_TRUE;
7659566063dSJacob Faibussowitsch       if (fas->full_total) PetscCall(SNESFASInterpolatedCoarseSolution(snes, X, X));
7669566063dSJacob Faibussowitsch       else PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
7676dfbafefSToby Isaac       fas->full_total = PETSC_FALSE;
7689566063dSJacob Faibussowitsch       PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
7698c94862eSPeter Brune       if (fas->level != 1) next->max_its -= 1;
7708c94862eSPeter Brune     } else {
771a3a80b83SMatthew G. Knepley       /* The smoother on the coarse level is the coarse solver */
7729566063dSJacob Faibussowitsch       PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
7738c94862eSPeter Brune     }
7748c94862eSPeter Brune     fas->full_stage = 1;
7758c94862eSPeter Brune   } else if (fas->full_stage == 1) {
7769566063dSJacob Faibussowitsch     if (snes->iter == 0) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
7778c94862eSPeter Brune     if (next) {
7789566063dSJacob Faibussowitsch       PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
7799566063dSJacob Faibussowitsch       PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
7808c94862eSPeter Brune     }
7818c94862eSPeter Brune   }
7828c94862eSPeter Brune   /* final v-cycle */
7838c94862eSPeter Brune   if (isFine) {
7848c94862eSPeter Brune     if (next) {
7859566063dSJacob Faibussowitsch       PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
7869566063dSJacob Faibussowitsch       PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
7878c94862eSPeter Brune     }
7888c94862eSPeter Brune   }
7893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7908c94862eSPeter Brune }
7918c94862eSPeter Brune 
792d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
793d71ae5a4SJacob Faibussowitsch {
79434d65b3cSPeter Brune   Vec  F, B;
79534d65b3cSPeter Brune   SNES next;
79634d65b3cSPeter Brune 
79734d65b3cSPeter Brune   PetscFunctionBegin;
79834d65b3cSPeter Brune   F = snes->vec_func;
79934d65b3cSPeter Brune   B = snes->vec_rhs;
8009566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
80134d65b3cSPeter Brune   if (next) {
8029566063dSJacob Faibussowitsch     PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
8039566063dSJacob Faibussowitsch     PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
80434d65b3cSPeter Brune   } else {
8059566063dSJacob Faibussowitsch     PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
80634d65b3cSPeter Brune   }
8073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
80834d65b3cSPeter Brune }
80934d65b3cSPeter Brune 
810fffbeea8SBarry Smith PetscBool  SNEScite       = PETSC_FALSE;
811fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
812fffbeea8SBarry Smith                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
813fffbeea8SBarry Smith                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
814fffbeea8SBarry Smith                             "  year = 2013,\n"
815fffbeea8SBarry Smith                             "  type = Preprint,\n"
816fffbeea8SBarry Smith                             "  number = {ANL/MCS-P2010-0112},\n"
817fffbeea8SBarry Smith                             "  institution = {Argonne National Laboratory}\n}\n";
818fffbeea8SBarry Smith 
819d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_FAS(SNES snes)
820d71ae5a4SJacob Faibussowitsch {
821f833ba53SLisandro Dalcin   PetscInt  i;
822ddb5aff1SPeter Brune   Vec       X, F;
823fe6f9142SPeter Brune   PetscReal fnorm;
824b17ce1afSJed Brown   SNES_FAS *fas = (SNES_FAS *)snes->data, *ffas;
825b17ce1afSJed Brown   DM        dm;
826e70c42e5SPeter Brune   PetscBool isFine;
827b17ce1afSJed Brown 
828421d9b32SPeter Brune   PetscFunctionBegin;
8290b121fc5SBarry Smith   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
830c579b300SPatrick Farrell 
8319566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
832fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
833fa9694d7SPeter Brune   X            = snes->vec_sol;
834f5a6d4f9SBarry Smith   F            = snes->vec_func;
835293a7e31SPeter Brune 
8369566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
837293a7e31SPeter Brune   /* norm setup */
8389566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
839fe6f9142SPeter Brune   snes->iter = 0;
840f833ba53SLisandro Dalcin   snes->norm = 0;
8419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
842e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
8439566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0));
8449566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, X, F));
8459566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0));
8461aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
847e4ed7901SPeter Brune 
8489566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
849422a814eSBarry Smith   SNESCheckFunctionNorm(snes, fnorm);
8509566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
851fe6f9142SPeter Brune   snes->norm = fnorm;
8529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
8539566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
854fe6f9142SPeter Brune 
855fe6f9142SPeter Brune   /* test convergence */
856*2d157150SStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
857*2d157150SStefano Zampini   PetscCall(SNESMonitor(snes, snes->iter, fnorm));
8583ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
859e4ed7901SPeter Brune 
860b9c2fdf1SPeter Brune   if (isFine) {
861b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
8629566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dm));
863b17ce1afSJed Brown     for (ffas = fas; ffas->next; ffas = (SNES_FAS *)ffas->next->data) {
864b17ce1afSJed Brown       DM dmcoarse;
8659566063dSJacob Faibussowitsch       PetscCall(SNESGetDM(ffas->next, &dmcoarse));
8669566063dSJacob Faibussowitsch       PetscCall(DMRestrict(dm, ffas->restrct, ffas->rscale, ffas->inject, dmcoarse));
867b17ce1afSJed Brown       dm = dmcoarse;
868b17ce1afSJed Brown     }
869b9c2fdf1SPeter Brune   }
870b17ce1afSJed Brown 
871f833ba53SLisandro Dalcin   for (i = 0; i < snes->max_its; i++) {
872fe6f9142SPeter Brune     /* Call general purpose update function */
873dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
874f833ba53SLisandro Dalcin 
87507144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
8769566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Multiplicative(snes, X));
8778c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
8789566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Additive(snes, X));
8798c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
8809566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Full(snes, X));
88134d65b3cSPeter Brune     } else if (fas->fastype == SNES_FAS_KASKADE) {
8829566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Kaskade(snes, X));
8836c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported FAS type");
884742fe5e2SPeter Brune 
885742fe5e2SPeter Brune     /* check for FAS cycle divergence */
8863ba16761SJacob Faibussowitsch     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
887b9c2fdf1SPeter Brune 
888c90fad12SPeter Brune     /* Monitor convergence */
8899566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
890c90fad12SPeter Brune     snes->iter = i + 1;
8919566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
8929566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
893*2d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, snes->norm));
8949566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
895c90fad12SPeter Brune     if (snes->reason) break;
896fe6f9142SPeter Brune   }
8973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
898421d9b32SPeter Brune }
89940244768SBarry Smith 
90040244768SBarry Smith /*MC
90140244768SBarry Smith 
90240244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
90340244768SBarry Smith 
90440244768SBarry Smith    The nonlinear problem is solved by correction using coarse versions
90540244768SBarry Smith    of the nonlinear problem.  This problem is perturbed so that a projected
90640244768SBarry Smith    solution of the fine problem elicits no correction from the coarse problem.
90740244768SBarry Smith 
908f6dfbefdSBarry Smith    Options Database Keys and Prefixes:
90940244768SBarry Smith +   -snes_fas_levels -  The number of levels
91040244768SBarry Smith .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
91140244768SBarry Smith .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  Additive or multiplicative cycle
912f6dfbefdSBarry Smith .   -snes_fas_galerkin<`PETSC_FALSE`> -  Form coarse problems by projection back upon the fine problem
91340244768SBarry Smith .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
91440244768SBarry Smith .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
91540244768SBarry Smith .   -snes_fas_monitor -  Monitor progress of all of the levels
916f6dfbefdSBarry Smith .   -snes_fas_full_downsweep<`PETSC_FALSE`> - call the downsmooth on the initial downsweep of full FAS
917f6dfbefdSBarry Smith .   -fas_levels_snes_ -  `SNES` options for all smoothers
918f6dfbefdSBarry Smith .   -fas_levels_cycle_snes_ -  `SNES` options for all cycles
919f6dfbefdSBarry Smith .   -fas_levels_i_snes_ -  `SNES` options for the smoothers on level i
920f6dfbefdSBarry Smith .   -fas_levels_i_cycle_snes_ - `SNES` options for the cycle on level i
921f6dfbefdSBarry Smith -   -fas_coarse_snes_ -  `SNES` options for the coarsest smoother
92240244768SBarry Smith 
923f6dfbefdSBarry Smith    Note:
924f6dfbefdSBarry Smith    The organization of the FAS solver is slightly different from the organization of `PCMG`
925f6dfbefdSBarry Smith    As each level has smoother `SNES` instances(down and potentially up) and a cycle `SNES` instance.
926f6dfbefdSBarry Smith    The cycle `SNES` instance may be used for monitoring convergence on a particular level.
92740244768SBarry Smith 
92840244768SBarry Smith    Level: beginner
92940244768SBarry Smith 
93040244768SBarry Smith    References:
931606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
93240244768SBarry Smith    SIAM Review, 57(4), 2015
93340244768SBarry Smith 
934f6dfbefdSBarry Smith .seealso: `PCMG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESFASSetRestriction()`, `SNESFASSetInjection()`,
935f6dfbefdSBarry Smith           `SNESFASFullGetTotal()`
93640244768SBarry Smith M*/
93740244768SBarry Smith 
938d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
939d71ae5a4SJacob Faibussowitsch {
94040244768SBarry Smith   SNES_FAS *fas;
94140244768SBarry Smith 
94240244768SBarry Smith   PetscFunctionBegin;
94340244768SBarry Smith   snes->ops->destroy        = SNESDestroy_FAS;
94440244768SBarry Smith   snes->ops->setup          = SNESSetUp_FAS;
94540244768SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
94640244768SBarry Smith   snes->ops->view           = SNESView_FAS;
94740244768SBarry Smith   snes->ops->solve          = SNESSolve_FAS;
94840244768SBarry Smith   snes->ops->reset          = SNESReset_FAS;
94940244768SBarry Smith 
95040244768SBarry Smith   snes->usesksp = PETSC_FALSE;
951efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
95240244768SBarry Smith 
95340244768SBarry Smith   if (!snes->tolerancesset) {
95440244768SBarry Smith     snes->max_funcs = 30000;
95540244768SBarry Smith     snes->max_its   = 10000;
95640244768SBarry Smith   }
95740244768SBarry Smith 
9584fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
9594fc747eaSLawrence Mitchell 
9604dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&fas));
96140244768SBarry Smith 
96240244768SBarry Smith   snes->data                  = (void *)fas;
96340244768SBarry Smith   fas->level                  = 0;
96440244768SBarry Smith   fas->levels                 = 1;
96540244768SBarry Smith   fas->n_cycles               = 1;
96640244768SBarry Smith   fas->max_up_it              = 1;
96740244768SBarry Smith   fas->max_down_it            = 1;
96840244768SBarry Smith   fas->smoothu                = NULL;
96940244768SBarry Smith   fas->smoothd                = NULL;
97040244768SBarry Smith   fas->next                   = NULL;
97140244768SBarry Smith   fas->previous               = NULL;
97240244768SBarry Smith   fas->fine                   = snes;
97340244768SBarry Smith   fas->interpolate            = NULL;
97440244768SBarry Smith   fas->restrct                = NULL;
97540244768SBarry Smith   fas->inject                 = NULL;
97640244768SBarry Smith   fas->usedmfornumberoflevels = PETSC_FALSE;
97740244768SBarry Smith   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
97840244768SBarry Smith   fas->full_downsweep         = PETSC_FALSE;
9796dfbafefSToby Isaac   fas->full_total             = PETSC_FALSE;
98040244768SBarry Smith 
98140244768SBarry Smith   fas->eventsmoothsetup    = 0;
98240244768SBarry Smith   fas->eventsmoothsolve    = 0;
98340244768SBarry Smith   fas->eventresidual       = 0;
98440244768SBarry Smith   fas->eventinterprestrict = 0;
9853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98640244768SBarry Smith }
987