xref: /petsc/src/snes/impls/fas/fas.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
640244768SBarry Smith static PetscErrorCode SNESReset_FAS(SNES snes)
7421d9b32SPeter Brune {
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));
20421d9b32SPeter Brune   PetscFunctionReturn(0);
21421d9b32SPeter Brune }
22421d9b32SPeter Brune 
2340244768SBarry Smith static PetscErrorCode SNESDestroy_FAS(SNES snes)
24421d9b32SPeter Brune {
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));
32421d9b32SPeter Brune   PetscFunctionReturn(0);
33421d9b32SPeter Brune }
34421d9b32SPeter Brune 
35f833ba53SLisandro Dalcin static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth)
36f833ba53SLisandro Dalcin {
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;
44f833ba53SLisandro Dalcin   if (!snes->linesearch) PetscFunctionReturn(0);
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));
52f833ba53SLisandro Dalcin   PetscFunctionReturn(0);
53f833ba53SLisandro Dalcin }
54f833ba53SLisandro Dalcin 
55f833ba53SLisandro Dalcin static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth)
56f833ba53SLisandro Dalcin {
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));
74f833ba53SLisandro Dalcin   PetscFunctionReturn(0);
75f833ba53SLisandro Dalcin }
76f833ba53SLisandro Dalcin 
7740244768SBarry Smith static PetscErrorCode SNESSetUp_FAS(SNES snes)
78421d9b32SPeter Brune {
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 */
101ab8d36c9SPeter Brune   if (!fas->smoothd) {
1029566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
103ab8d36c9SPeter Brune   }
104ab8d36c9SPeter Brune 
10579d9a41aSPeter Brune   if (snes->dm) {
106ab8d36c9SPeter Brune     /* set the smoother DMs properly */
1079566063dSJacob Faibussowitsch     if (fas->smoothu) PetscCall(SNESSetDM(fas->smoothu, snes->dm));
1089566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(fas->smoothd, snes->dm));
10979d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
110ab8d36c9SPeter Brune     if (next) {
11179d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
112ab8d36c9SPeter Brune       if (!next->dm) {
1139566063dSJacob Faibussowitsch         PetscCall(DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm));
1149566063dSJacob Faibussowitsch         PetscCall(SNESSetDM(next, next->dm));
11579d9a41aSPeter Brune       }
11679d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
11779d9a41aSPeter Brune       if (!fas->interpolate) {
1189566063dSJacob Faibussowitsch         PetscCall(DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale));
119bccf9bb3SJed Brown         if (!fas->restrct) {
1209566063dSJacob Faibussowitsch           PetscCall(DMHasCreateRestriction(next->dm, &hasCreateRestriction));
1210a7266b2SPatrick Farrell           /* DM can create restrictions, use that */
1220a7266b2SPatrick Farrell           if (hasCreateRestriction) {
1239566063dSJacob Faibussowitsch             PetscCall(DMCreateRestriction(next->dm, snes->dm, &fas->restrct));
1240a7266b2SPatrick Farrell           } else {
1259566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)fas->interpolate));
12679d9a41aSPeter Brune             fas->restrct = fas->interpolate;
12779d9a41aSPeter Brune           }
128bccf9bb3SJed Brown         }
1290a7266b2SPatrick Farrell       }
13079d9a41aSPeter Brune       /* set the injection from the DM */
13179d9a41aSPeter Brune       if (!fas->inject) {
1329566063dSJacob Faibussowitsch         PetscCall(DMHasCreateInjection(next->dm, &hasCreateInjection));
133f833ba53SLisandro Dalcin         if (hasCreateInjection) {
1349566063dSJacob Faibussowitsch           PetscCall(DMCreateInjection(next->dm, snes->dm, &fas->inject));
13579d9a41aSPeter Brune         }
13679d9a41aSPeter Brune       }
13779d9a41aSPeter Brune     }
13823e68893SLawrence Mitchell   }
139f833ba53SLisandro Dalcin 
14079d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
14179d9a41aSPeter Brune   if (fas->galerkin) {
1421baa6e33SBarry Smith     if (next) PetscCall(SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next));
1431aa26658SKarl Rupp     if (fas->smoothd && fas->level != fas->levels - 1) {
1449566063dSJacob Faibussowitsch       PetscCall(SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes));
1451aa26658SKarl Rupp     }
1461aa26658SKarl Rupp     if (fas->smoothu && fas->level != fas->levels - 1) {
1479566063dSJacob Faibussowitsch       PetscCall(SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes));
1481aa26658SKarl Rupp     }
14979d9a41aSPeter Brune   }
15079d9a41aSPeter Brune 
151534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
152534ebe21SPeter Brune   if (fas->smoothd) {
153bc3f2f05SPeter Brune     if (fas->level == 0 && fas->levels != 1) {
1549566063dSJacob Faibussowitsch       PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE));
155534ebe21SPeter Brune     } else {
1569566063dSJacob Faibussowitsch       PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY));
157534ebe21SPeter Brune     }
1589566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd));
159534ebe21SPeter Brune   }
160534ebe21SPeter Brune 
161534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
162534ebe21SPeter Brune   if (fas->smoothu) {
163534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
1649566063dSJacob Faibussowitsch       PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE));
165534ebe21SPeter Brune     } else {
1669566063dSJacob Faibussowitsch       PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY));
167534ebe21SPeter Brune     }
1689566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu));
169534ebe21SPeter Brune   }
170d06165b7SPeter Brune 
171ab8d36c9SPeter Brune   if (next) {
17279d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
1739566063dSJacob Faibussowitsch     if (!next->vec_sol) PetscCall(SNESFASCreateCoarseVec(snes,&next->vec_sol));
1749566063dSJacob Faibussowitsch     if (!next->vec_rhs) PetscCall(SNESFASCreateCoarseVec(snes,&next->vec_rhs));
1759566063dSJacob Faibussowitsch     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next));
1769566063dSJacob Faibussowitsch     PetscCall(SNESFASSetUpLineSearch_Private(snes, next));
1779566063dSJacob Faibussowitsch     PetscCall(SNESSetUp(next));
17879d9a41aSPeter Brune   }
179f833ba53SLisandro Dalcin 
1806273346dSPeter Brune   /* setup FAS work vectors */
1816273346dSPeter Brune   if (fas->galerkin) {
1829566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(snes->vec_sol, &fas->Xg));
1839566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(snes->vec_sol, &fas->Fg));
1846273346dSPeter Brune   }
185421d9b32SPeter Brune   PetscFunctionReturn(0);
186421d9b32SPeter Brune }
187421d9b32SPeter Brune 
188*dbbe0bcdSBarry Smith static PetscErrorCode SNESSetFromOptions_FAS(SNES snes,PetscOptionItems *PetscOptionsObject)
189421d9b32SPeter Brune {
190ee78dd50SPeter Brune   SNES_FAS       *fas   = (SNES_FAS*) snes->data;
191ee78dd50SPeter Brune   PetscInt       levels = 1;
19287f44e3fSPeter Brune   PetscBool      flg    = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE;
19307144faaSPeter Brune   SNESFASType    fastype;
194fde0ff24SPeter Brune   const char     *optionsprefix;
195f1c6b773SPeter Brune   SNESLineSearch linesearch;
19666585501SPeter Brune   PetscInt       m, n_up, n_down;
197ab8d36c9SPeter Brune   SNES           next;
198ab8d36c9SPeter Brune   PetscBool      isFine;
199421d9b32SPeter Brune 
200421d9b32SPeter Brune   PetscFunctionBegin;
2019566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
202d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SNESFAS Options-----------------------------------");
203ee78dd50SPeter Brune 
204ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
205ab8d36c9SPeter Brune   if (isFine) {
2069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg));
207c732cbdbSBarry Smith     if (!flg && snes->dm) {
2089566063dSJacob Faibussowitsch       PetscCall(DMGetRefineLevel(snes->dm,&levels));
209c732cbdbSBarry Smith       levels++;
210d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
211c732cbdbSBarry Smith     }
2129566063dSJacob Faibussowitsch     PetscCall(SNESFASSetLevels(snes, levels, NULL));
21307144faaSPeter Brune     fastype = fas->fastype;
2149566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg));
2151baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetType(snes, fastype));
216ee78dd50SPeter Brune 
2179566063dSJacob Faibussowitsch     PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
2189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg));
2191baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetCycles(snes, m));
2209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg));
2211baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetContinuation(snes,continuationflg));
222fde0ff24SPeter Brune 
2239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg));
2241baa6e33SBarry Smith     if (flg) PetscCall(SNESFASSetGalerkin(snes, galerkinflg));
225ee78dd50SPeter Brune 
226928e959bSPeter Brune     if (fas->fastype == SNES_FAS_FULL) {
2279566063dSJacob 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));
2289566063dSJacob Faibussowitsch       if (flg) PetscCall(SNESFASFullSetDownSweep(snes,fas->full_downsweep));
2299566063dSJacob 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));
2309566063dSJacob Faibussowitsch       if (flg) PetscCall(SNESFASFullSetTotal(snes,fas->full_total));
231928e959bSPeter Brune     }
232928e959bSPeter Brune 
2339566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg));
234162d76ddSPeter Brune 
2359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg));
236162d76ddSPeter Brune 
237d142ab34SLawrence Mitchell     {
238d142ab34SLawrence Mitchell       PetscViewer viewer;
239d142ab34SLawrence Mitchell       PetscViewerFormat format;
2409566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_fas_monitor",&viewer,&format,&monflg));
241d142ab34SLawrence Mitchell       if (monflg) {
242d142ab34SLawrence Mitchell         PetscViewerAndFormat *vf;
2439566063dSJacob Faibussowitsch         PetscCall(PetscViewerAndFormatCreate(viewer,format,&vf));
2449566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)viewer));
2459566063dSJacob Faibussowitsch         PetscCall(SNESFASSetMonitor(snes,vf,PETSC_TRUE));
246d142ab34SLawrence Mitchell       }
247d142ab34SLawrence Mitchell     }
2480dd27c6cSPeter Brune     flg    = PETSC_FALSE;
2490dd27c6cSPeter Brune     monflg = PETSC_TRUE;
2509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg));
2519566063dSJacob Faibussowitsch     if (flg) PetscCall(SNESFASSetLog(snes,monflg));
252ab8d36c9SPeter Brune   }
253ee78dd50SPeter Brune 
254d0609cedSBarry Smith   PetscOptionsHeadEnd();
255f833ba53SLisandro Dalcin 
2568cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
2571baa6e33SBarry Smith   if (upflg) PetscCall(SNESFASSetNumberSmoothUp(snes,n_up));
2581baa6e33SBarry Smith   if (downflg) PetscCall(SNESFASSetNumberSmoothDown(snes,n_down));
259eff52c0eSPeter Brune 
2609e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
2619e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
2629e764e56SPeter Brune     if (!snes->linesearch) {
2639566063dSJacob Faibussowitsch       PetscCall(SNESGetLineSearch(snes, &linesearch));
2649566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
2659e764e56SPeter Brune     }
2669e764e56SPeter Brune   }
2679e764e56SPeter Brune 
268ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
2699566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
2709566063dSJacob Faibussowitsch   if (next) PetscCall(SNESSetFromOptions(next));
271421d9b32SPeter Brune   PetscFunctionReturn(0);
272421d9b32SPeter Brune }
273421d9b32SPeter Brune 
2749804daf3SBarry Smith #include <petscdraw.h>
27540244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
276421d9b32SPeter Brune {
277421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
278656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
279ab8d36c9SPeter Brune   PetscInt       i;
280ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
281421d9b32SPeter Brune 
282421d9b32SPeter Brune   PetscFunctionBegin;
2839566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
284ab8d36c9SPeter Brune   if (isFine) {
2859566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2869566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
287421d9b32SPeter Brune     if (iascii) {
28863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%" PetscInt_FMT ", cycles=%" PetscInt_FMT "\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles));
289ab8d36c9SPeter Brune       if (fas->galerkin) {
2909566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  Using Galerkin computed coarse grid function evaluation\n"));
291421d9b32SPeter Brune       } else {
2929566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"  Not using Galerkin computed coarse grid function evaluation\n"));
293421d9b32SPeter Brune       }
294ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
2959566063dSJacob Faibussowitsch         PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
2969566063dSJacob Faibussowitsch         PetscCall(SNESFASCycleGetSmootherUp(levelsnes, &smoothu));
2979566063dSJacob Faibussowitsch         PetscCall(SNESFASCycleGetSmootherDown(levelsnes, &smoothd));
298ab8d36c9SPeter Brune         if (!i) {
29963a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"  Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n",i));
300421d9b32SPeter Brune         } else {
30163a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"  Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n",i));
302421d9b32SPeter Brune         }
3039566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
304166b3ea4SJed Brown         if (smoothd) {
3059566063dSJacob Faibussowitsch           PetscCall(SNESView(smoothd,viewer));
306166b3ea4SJed Brown         } else {
3079566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"Not yet available\n"));
308166b3ea4SJed Brown         }
3099566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
310ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
3119566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) same as down solver (pre-smoother)\n"));
312ab8d36c9SPeter Brune         } else if (i) {
31363a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n",i));
3149566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
315166b3ea4SJed Brown           if (smoothu) {
3169566063dSJacob Faibussowitsch             PetscCall(SNESView(smoothu,viewer));
317166b3ea4SJed Brown           } else {
3189566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer,"Not yet available\n"));
319166b3ea4SJed Brown           }
3209566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
321ab8d36c9SPeter Brune         }
322ab8d36c9SPeter Brune       }
323656ede7eSPeter Brune     } else if (isdraw) {
324656ede7eSPeter Brune       PetscDraw draw;
325b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
326656ede7eSPeter Brune       SNES_FAS  *curfas = fas;
3279566063dSJacob Faibussowitsch       PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
3289566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetCurrentPoint(draw,&x,&y));
3299566063dSJacob Faibussowitsch       PetscCall(PetscDrawStringGetSize(draw,&wth,&th));
330656ede7eSPeter Brune       bottom = y - th;
331656ede7eSPeter Brune       while (curfas) {
332b4375e8dSPeter Brune         if (!curfas->smoothu) {
3339566063dSJacob Faibussowitsch           PetscCall(PetscDrawPushCurrentPoint(draw,x,bottom));
3349566063dSJacob Faibussowitsch           if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd,viewer));
3359566063dSJacob Faibussowitsch           PetscCall(PetscDrawPopCurrentPoint(draw));
336b4375e8dSPeter Brune         } else {
337b4375e8dSPeter Brune           w    = 0.5*PetscMin(1.0-x,x);
3389566063dSJacob Faibussowitsch           PetscCall(PetscDrawPushCurrentPoint(draw,x-w,bottom));
3399566063dSJacob Faibussowitsch           if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd,viewer));
3409566063dSJacob Faibussowitsch           PetscCall(PetscDrawPopCurrentPoint(draw));
3419566063dSJacob Faibussowitsch           PetscCall(PetscDrawPushCurrentPoint(draw,x+w,bottom));
3429566063dSJacob Faibussowitsch           if (curfas->smoothu) PetscCall(SNESView(curfas->smoothu,viewer));
3439566063dSJacob Faibussowitsch           PetscCall(PetscDrawPopCurrentPoint(draw));
344b4375e8dSPeter Brune         }
345656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
346656ede7eSPeter Brune         bottom -= 5*th;
3471aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
3480298fd71SBarry Smith         else curfas = NULL;
349656ede7eSPeter Brune       }
350421d9b32SPeter Brune     }
351ab8d36c9SPeter Brune   }
352421d9b32SPeter Brune   PetscFunctionReturn(0);
353421d9b32SPeter Brune }
354421d9b32SPeter Brune 
35539bd7f45SPeter Brune /*
35639bd7f45SPeter Brune Defines the action of the downsmoother
35739bd7f45SPeter Brune  */
35840244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
359b9c2fdf1SPeter Brune {
360742fe5e2SPeter Brune   SNESConvergedReason reason;
361ab8d36c9SPeter Brune   Vec                 FPC;
362ab8d36c9SPeter Brune   SNES                smoothd;
3636cbb2f26SLawrence Mitchell   PetscBool           flg;
3640dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
3656e111a19SKarl Rupp 
366421d9b32SPeter Brune   PetscFunctionBegin;
3679566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetSmootherDown(snes, &smoothd));
3689566063dSJacob Faibussowitsch   PetscCall(SNESSetInitialFunction(smoothd, F));
3699566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve,smoothd,B,X,0));
3709566063dSJacob Faibussowitsch   PetscCall(SNESSolve(smoothd, B, X));
3719566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve,smoothd,B,X,0));
372742fe5e2SPeter Brune   /* check convergence reason for the smoother */
3739566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(smoothd,&reason));
3741ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
375742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
376742fe5e2SPeter Brune     PetscFunctionReturn(0);
377742fe5e2SPeter Brune   }
3786cbb2f26SLawrence Mitchell 
3799566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(smoothd, &FPC, NULL, NULL));
3809566063dSJacob Faibussowitsch   PetscCall(SNESGetAlwaysComputesFinalResidual(smoothd, &flg));
3816cbb2f26SLawrence Mitchell   if (!flg) {
3829566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(smoothd, X, FPC));
3836cbb2f26SLawrence Mitchell   }
3849566063dSJacob Faibussowitsch   PetscCall(VecCopy(FPC, F));
3859566063dSJacob Faibussowitsch   if (fnorm) PetscCall(VecNorm(F,NORM_2,fnorm));
38639bd7f45SPeter Brune   PetscFunctionReturn(0);
38739bd7f45SPeter Brune }
38839bd7f45SPeter Brune 
38939bd7f45SPeter Brune /*
39007144faaSPeter Brune Defines the action of the upsmoother
39139bd7f45SPeter Brune  */
39240244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
3930adebc6cSBarry Smith {
39439bd7f45SPeter Brune   SNESConvergedReason reason;
395ab8d36c9SPeter Brune   Vec                 FPC;
396ab8d36c9SPeter Brune   SNES                smoothu;
3976cbb2f26SLawrence Mitchell   PetscBool           flg;
3980dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
399ab8d36c9SPeter Brune 
4006e111a19SKarl Rupp   PetscFunctionBegin;
4019566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetSmootherUp(snes, &smoothu));
4029566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve,smoothu,0,0,0));
4039566063dSJacob Faibussowitsch   PetscCall(SNESSolve(smoothu, B, X));
4049566063dSJacob Faibussowitsch   if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve,smoothu,0,0,0));
40539bd7f45SPeter Brune   /* check convergence reason for the smoother */
4069566063dSJacob Faibussowitsch   PetscCall(SNESGetConvergedReason(smoothu,&reason));
4071ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
40839bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
40939bd7f45SPeter Brune     PetscFunctionReturn(0);
41039bd7f45SPeter Brune   }
4119566063dSJacob Faibussowitsch   PetscCall(SNESGetFunction(smoothu, &FPC, NULL, NULL));
4129566063dSJacob Faibussowitsch   PetscCall(SNESGetAlwaysComputesFinalResidual(smoothu, &flg));
4136cbb2f26SLawrence Mitchell   if (!flg) {
4149566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(smoothu, X, FPC));
4156cbb2f26SLawrence Mitchell   }
4169566063dSJacob Faibussowitsch   PetscCall(VecCopy(FPC, F));
4179566063dSJacob Faibussowitsch   if (fnorm) PetscCall(VecNorm(F,NORM_2,fnorm));
41839bd7f45SPeter Brune   PetscFunctionReturn(0);
41939bd7f45SPeter Brune }
42039bd7f45SPeter Brune 
421938e4a01SJed Brown /*@
422938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
423938e4a01SJed Brown 
424938e4a01SJed Brown    Collective
425938e4a01SJed Brown 
4264165533cSJose E. Roman    Input Parameter:
427938e4a01SJed Brown .  snes - SNESFAS
428938e4a01SJed Brown 
4294165533cSJose E. Roman    Output Parameter:
430938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
431938e4a01SJed Brown 
432938e4a01SJed Brown    Level: developer
433938e4a01SJed Brown 
434db781477SPatrick Sanan .seealso: `SNESFASSetRestriction()`, `SNESFASRestrict()`
435938e4a01SJed Brown @*/
436938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
437938e4a01SJed Brown {
438f833ba53SLisandro Dalcin   SNES_FAS       *fas;
439938e4a01SJed Brown 
440938e4a01SJed Brown   PetscFunctionBegin;
441f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
442f833ba53SLisandro Dalcin   PetscValidPointer(Xcoarse,2);
443f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
4441aa26658SKarl Rupp   if (fas->rscale) {
4459566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(fas->rscale,Xcoarse));
446f5af7f23SKarl Rupp   } else if (fas->inject) {
4479566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(fas->inject,Xcoarse,NULL));
44813903a91SSatish Balay   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");
449938e4a01SJed Brown   PetscFunctionReturn(0);
450938e4a01SJed Brown }
451938e4a01SJed Brown 
452e9923e8dSJed Brown /*@
453e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
454e9923e8dSJed Brown 
455e9923e8dSJed Brown    Collective
456e9923e8dSJed Brown 
4574165533cSJose E. Roman    Input Parameters:
458e9923e8dSJed Brown +  fine - SNES from which to restrict
459e9923e8dSJed Brown -  Xfine - vector to restrict
460e9923e8dSJed Brown 
4614165533cSJose E. Roman    Output Parameter:
462e9923e8dSJed Brown .  Xcoarse - result of restriction
463e9923e8dSJed Brown 
464e9923e8dSJed Brown    Level: developer
465e9923e8dSJed Brown 
466db781477SPatrick Sanan .seealso: `SNESFASSetRestriction()`, `SNESFASSetInjection()`
467e9923e8dSJed Brown @*/
468e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
469e9923e8dSJed Brown {
470f833ba53SLisandro Dalcin   SNES_FAS       *fas;
471e9923e8dSJed Brown 
472e9923e8dSJed Brown   PetscFunctionBegin;
473f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(fine,SNES_CLASSID,1,SNESFAS);
474e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
475e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
476f833ba53SLisandro Dalcin   fas = (SNES_FAS*)fine->data;
477e9923e8dSJed Brown   if (fas->inject) {
4789566063dSJacob Faibussowitsch     PetscCall(MatRestrict(fas->inject,Xfine,Xcoarse));
479e9923e8dSJed Brown   } else {
4809566063dSJacob Faibussowitsch     PetscCall(MatRestrict(fas->restrct,Xfine,Xcoarse));
4819566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse));
482e9923e8dSJed Brown   }
483e9923e8dSJed Brown   PetscFunctionReturn(0);
484e9923e8dSJed Brown }
485e9923e8dSJed Brown 
48639bd7f45SPeter Brune /*
48739bd7f45SPeter Brune 
4886dfbafefSToby Isaac Performs a variant of FAS using the interpolated total coarse solution
4896dfbafefSToby Isaac 
4906dfbafefSToby Isaac fine problem:   F(x) = b
4916dfbafefSToby Isaac coarse problem: F^c(x^c) = Rb, Initial guess Rx
4926dfbafefSToby Isaac interpolated solution: x^f = I x^c (total solution interpolation
4936dfbafefSToby Isaac 
4946dfbafefSToby Isaac  */
4956dfbafefSToby Isaac static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new)
4966dfbafefSToby Isaac {
4976dfbafefSToby Isaac   Vec                 X_c, B_c;
4986dfbafefSToby Isaac   SNESConvergedReason reason;
4996dfbafefSToby Isaac   SNES                next;
5006dfbafefSToby Isaac   Mat                 restrct, interpolate;
5016dfbafefSToby Isaac   SNES_FAS            *fasc;
5026dfbafefSToby Isaac 
5036dfbafefSToby Isaac   PetscFunctionBegin;
5049566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
5056dfbafefSToby Isaac   if (next) {
5066dfbafefSToby Isaac     fasc = (SNES_FAS*)next->data;
5076dfbafefSToby Isaac 
5089566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
5099566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
5106dfbafefSToby Isaac 
5116dfbafefSToby Isaac     X_c  = next->vec_sol;
5126dfbafefSToby Isaac 
5139566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0));
5146dfbafefSToby Isaac     /* restrict the total solution: Rb */
5159566063dSJacob Faibussowitsch     PetscCall(SNESFASRestrict(snes, X, X_c));
5166dfbafefSToby Isaac     B_c = next->vec_rhs;
5176dfbafefSToby Isaac     if (snes->vec_rhs) {
5186dfbafefSToby Isaac       /* restrict the total rhs defect: Rb */
5199566063dSJacob Faibussowitsch       PetscCall(MatRestrict(restrct, snes->vec_rhs, B_c));
5206dfbafefSToby Isaac     } else {
5219566063dSJacob Faibussowitsch       PetscCall(VecSet(B_c, 0.));
5226dfbafefSToby Isaac     }
5239566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0));
5246dfbafefSToby Isaac 
5259566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next, B_c, X_c));
5269566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next,&reason));
5276dfbafefSToby Isaac     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
5286dfbafefSToby Isaac       snes->reason = SNES_DIVERGED_INNER;
5296dfbafefSToby Isaac       PetscFunctionReturn(0);
5306dfbafefSToby Isaac     }
5316dfbafefSToby Isaac     /* x^f <- Ix^c*/
5326dfbafefSToby Isaac     DM dmc,dmf;
5336dfbafefSToby Isaac 
5349566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(next, &dmc));
5359566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes, &dmf));
5369566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0));
5379566063dSJacob Faibussowitsch     PetscCall(DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new));
5389566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0));
5399566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) X_c, "Coarse solution"));
5409566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view"));
5419566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) X_new, "Updated Fine solution"));
5429566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view"));
5436dfbafefSToby Isaac   }
5446dfbafefSToby Isaac   PetscFunctionReturn(0);
5456dfbafefSToby Isaac }
5466dfbafefSToby Isaac 
5476dfbafefSToby Isaac /*
5486dfbafefSToby Isaac 
54939bd7f45SPeter Brune Performs the FAS coarse correction as:
55039bd7f45SPeter Brune 
551b5527d98SMatthew G. Knepley fine problem:   F(x) = b
552b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c
55339bd7f45SPeter Brune 
554b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
55539bd7f45SPeter Brune 
55639bd7f45SPeter Brune  */
5570adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
5580adebc6cSBarry Smith {
55939bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
56039bd7f45SPeter Brune   SNESConvergedReason reason;
561ab8d36c9SPeter Brune   SNES                next;
562ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
5630dd27c6cSPeter Brune   SNES_FAS            *fasc;
5645fd66863SKarl Rupp 
56539bd7f45SPeter Brune   PetscFunctionBegin;
5669566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
567ab8d36c9SPeter Brune   if (next) {
5680dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
5690dd27c6cSPeter Brune 
5709566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
5719566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
572ab8d36c9SPeter Brune 
573ab8d36c9SPeter Brune     X_c  = next->vec_sol;
574ab8d36c9SPeter Brune     Xo_c = next->work[0];
575ab8d36c9SPeter Brune     F_c  = next->vec_func;
576ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
577efe1f98aSPeter Brune 
5789566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0));
5799566063dSJacob Faibussowitsch     PetscCall(SNESFASRestrict(snes, X, Xo_c));
5805609cbf2SMatthew G. Knepley     /* restrict the defect: R(F(x) - b) */
5819566063dSJacob Faibussowitsch     PetscCall(MatRestrict(restrct, F, B_c));
5829566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0));
5830dd27c6cSPeter Brune 
5849566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual,next,0,0,0));
5855609cbf2SMatthew G. Knepley     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
5869566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(next, Xo_c, F_c));
5879566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual,next,0,0,0));
5880dd27c6cSPeter Brune 
5890dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
5909566063dSJacob Faibussowitsch     PetscCall(VecCopy(B_c, X_c));
5919566063dSJacob Faibussowitsch     PetscCall(VecCopy(F_c, B_c));
5929566063dSJacob Faibussowitsch     PetscCall(VecCopy(X_c, F_c));
593ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
5949566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xo_c, X_c));
595c90fad12SPeter Brune 
596c90fad12SPeter Brune     /* recurse to the next level */
5979566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next, F_c));
5989566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next, B_c, X_c));
5999566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next,&reason));
600742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
601742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
602742fe5e2SPeter Brune       PetscFunctionReturn(0);
603742fe5e2SPeter Brune     }
604fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
6059566063dSJacob Faibussowitsch     PetscCall(VecAXPY(X_c, -1.0, Xo_c));
6060dd27c6cSPeter Brune 
6079566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0));
6089566063dSJacob Faibussowitsch     PetscCall(MatInterpolateAdd(interpolate, X_c, X, X_new));
6099566063dSJacob Faibussowitsch     if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0));
6109566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) X_c, "Coarse correction"));
6119566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view"));
6129566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) X_new, "Updated Fine solution"));
6139566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view"));
614293a7e31SPeter Brune   }
61539bd7f45SPeter Brune   PetscFunctionReturn(0);
61639bd7f45SPeter Brune }
61739bd7f45SPeter Brune 
61839bd7f45SPeter Brune /*
61939bd7f45SPeter Brune 
62039bd7f45SPeter Brune The additive cycle looks like:
62139bd7f45SPeter Brune 
62207144faaSPeter Brune xhat = x
62307144faaSPeter Brune xhat = dS(x, b)
62407144faaSPeter Brune x = coarsecorrection(xhat, b_d)
62507144faaSPeter Brune x = x + nu*(xhat - x);
62639bd7f45SPeter Brune (optional) x = uS(x, b)
62739bd7f45SPeter Brune 
62839bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
62939bd7f45SPeter Brune 
63039bd7f45SPeter Brune  */
63140244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
6320adebc6cSBarry Smith {
63307144faaSPeter Brune   Vec                  F, B, Xhat;
63422c1e704SPeter Brune   Vec                  X_c, Xo_c, F_c, B_c;
63507144faaSPeter Brune   SNESConvergedReason  reason;
63622c1e704SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
637422a814eSBarry Smith   SNESLineSearchReason lsresult;
638ab8d36c9SPeter Brune   SNES                 next;
639ab8d36c9SPeter Brune   Mat                  restrct, interpolate;
6400dd27c6cSPeter Brune   SNES_FAS             *fas = (SNES_FAS*)snes->data,*fasc;
6410dd27c6cSPeter Brune 
64239bd7f45SPeter Brune   PetscFunctionBegin;
6439566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
64439bd7f45SPeter Brune   F    = snes->vec_func;
64539bd7f45SPeter Brune   B    = snes->vec_rhs;
646e7f468e7SPeter Brune   Xhat = snes->work[1];
6479566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, Xhat));
64807144faaSPeter Brune   /* recurse first */
649ab8d36c9SPeter Brune   if (next) {
6500dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
6519566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetRestriction(snes, &restrct));
6529566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate));
6539566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual,snes,0,0,0));
6549566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, Xhat, F));
6559566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual,snes,0,0,0));
6569566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
657ab8d36c9SPeter Brune     X_c  = next->vec_sol;
658ab8d36c9SPeter Brune     Xo_c = next->work[0];
659ab8d36c9SPeter Brune     F_c  = next->vec_func;
660ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
66139bd7f45SPeter Brune 
6629566063dSJacob Faibussowitsch     PetscCall(SNESFASRestrict(snes,Xhat,Xo_c));
66307144faaSPeter Brune     /* restrict the defect */
6649566063dSJacob Faibussowitsch     PetscCall(MatRestrict(restrct, F, B_c));
66507144faaSPeter Brune 
66607144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
6679566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual,next,0,0,0));
6689566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(next, Xo_c, F_c));
6699566063dSJacob Faibussowitsch     if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual,next,0,0,0));
6709566063dSJacob Faibussowitsch     PetscCall(VecCopy(B_c, X_c));
6719566063dSJacob Faibussowitsch     PetscCall(VecCopy(F_c, B_c));
6729566063dSJacob Faibussowitsch     PetscCall(VecCopy(X_c, F_c));
67307144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
6749566063dSJacob Faibussowitsch     PetscCall(VecCopy(Xo_c, X_c));
67507144faaSPeter Brune 
67607144faaSPeter Brune     /* recurse */
6779566063dSJacob Faibussowitsch     PetscCall(SNESSetInitialFunction(next, F_c));
6789566063dSJacob Faibussowitsch     PetscCall(SNESSolve(next, B_c, X_c));
67907144faaSPeter Brune 
68007144faaSPeter Brune     /* smooth on this level */
6819566063dSJacob Faibussowitsch     PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &fnorm));
68207144faaSPeter Brune 
6839566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(next,&reason));
68407144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
68507144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
68607144faaSPeter Brune       PetscFunctionReturn(0);
68707144faaSPeter Brune     }
68807144faaSPeter Brune 
68907144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
6909566063dSJacob Faibussowitsch     PetscCall(VecAYPX(X_c, -1.0, Xo_c));
6919566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(interpolate, X_c, Xhat));
69207144faaSPeter Brune 
693ddebd997SPeter Brune     /* additive correction of the coarse direction*/
6949566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat));
6959566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
6969566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm));
697422a814eSBarry Smith     if (lsresult) {
6989e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6999e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
7009e764e56SPeter Brune         PetscFunctionReturn(0);
7019e764e56SPeter Brune       }
7029e764e56SPeter Brune     }
70307144faaSPeter Brune   } else {
7049566063dSJacob Faibussowitsch     PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
70507144faaSPeter Brune   }
70639bd7f45SPeter Brune   PetscFunctionReturn(0);
70739bd7f45SPeter Brune }
70839bd7f45SPeter Brune 
70939bd7f45SPeter Brune /*
71039bd7f45SPeter Brune 
71139bd7f45SPeter Brune Defines the FAS cycle as:
71239bd7f45SPeter Brune 
713b5527d98SMatthew G. Knepley fine problem: F(x) = b
71439bd7f45SPeter Brune coarse problem: F^c(x) = b^c
71539bd7f45SPeter Brune 
716b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
71739bd7f45SPeter Brune 
71839bd7f45SPeter Brune correction:
71939bd7f45SPeter Brune 
72039bd7f45SPeter Brune x = x + I(x^c - Rx)
72139bd7f45SPeter Brune 
72239bd7f45SPeter Brune  */
72340244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
7240adebc6cSBarry Smith {
72539bd7f45SPeter Brune 
72639bd7f45SPeter Brune   Vec            F,B;
72734d65b3cSPeter Brune   SNES           next;
72839bd7f45SPeter Brune 
72939bd7f45SPeter Brune   PetscFunctionBegin;
73039bd7f45SPeter Brune   F = snes->vec_func;
73139bd7f45SPeter Brune   B = snes->vec_rhs;
73239bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7339566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes, &next));
7349566063dSJacob Faibussowitsch   PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm));
73534d65b3cSPeter Brune   if (next) {
7369566063dSJacob Faibussowitsch     PetscCall(SNESFASCoarseCorrection(snes, X, F, X));
7379566063dSJacob Faibussowitsch     PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm));
738fe6f9142SPeter Brune   }
739fa9694d7SPeter Brune   PetscFunctionReturn(0);
740421d9b32SPeter Brune }
741421d9b32SPeter Brune 
74240244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
7438c94862eSPeter Brune {
7448c94862eSPeter Brune   SNES           next;
7458c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7468c94862eSPeter Brune   PetscBool      isFine;
7478c94862eSPeter Brune 
7488c94862eSPeter Brune   PetscFunctionBegin;
7498c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7509566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes,&isFine));
7519566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes,&next));
7528c94862eSPeter Brune   fas->full_stage = 0;
7539566063dSJacob Faibussowitsch   if (next) PetscCall(SNESFASCycleSetupPhase_Full(next));
7548c94862eSPeter Brune   PetscFunctionReturn(0);
7558c94862eSPeter Brune }
7568c94862eSPeter Brune 
75740244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
7588c94862eSPeter Brune {
7598c94862eSPeter Brune   Vec            F,B;
7608c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7618c94862eSPeter Brune   PetscBool      isFine;
7628c94862eSPeter Brune   SNES           next;
7638c94862eSPeter Brune 
7648c94862eSPeter Brune   PetscFunctionBegin;
7658c94862eSPeter Brune   F = snes->vec_func;
7668c94862eSPeter Brune   B = snes->vec_rhs;
7679566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes,&isFine));
7689566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes,&next));
7698c94862eSPeter Brune 
7701baa6e33SBarry Smith   if (isFine) PetscCall(SNESFASCycleSetupPhase_Full(snes));
7718c94862eSPeter Brune 
7728c94862eSPeter Brune   if (fas->full_stage == 0) {
773928e959bSPeter Brune     /* downsweep */
7748c94862eSPeter Brune     if (next) {
7758c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
7769566063dSJacob Faibussowitsch       if (fas->full_downsweep) PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm));
777e140b65aSPatrick Farrell       fas->full_downsweep = PETSC_TRUE;
7789566063dSJacob Faibussowitsch       if (fas->full_total) PetscCall(SNESFASInterpolatedCoarseSolution(snes,X,X));
7799566063dSJacob Faibussowitsch       else                 PetscCall(SNESFASCoarseCorrection(snes,X,F,X));
7806dfbafefSToby Isaac       fas->full_total = PETSC_FALSE;
7819566063dSJacob Faibussowitsch       PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm));
7828c94862eSPeter Brune       if (fas->level != 1) next->max_its -= 1;
7838c94862eSPeter Brune     } else {
784a3a80b83SMatthew G. Knepley       /* The smoother on the coarse level is the coarse solver */
7859566063dSJacob Faibussowitsch       PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm));
7868c94862eSPeter Brune     }
7878c94862eSPeter Brune     fas->full_stage = 1;
7888c94862eSPeter Brune   } else if (fas->full_stage == 1) {
7899566063dSJacob Faibussowitsch     if (snes->iter == 0) PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm));
7908c94862eSPeter Brune     if (next) {
7919566063dSJacob Faibussowitsch       PetscCall(SNESFASCoarseCorrection(snes,X,F,X));
7929566063dSJacob Faibussowitsch       PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm));
7938c94862eSPeter Brune     }
7948c94862eSPeter Brune   }
7958c94862eSPeter Brune   /* final v-cycle */
7968c94862eSPeter Brune   if (isFine) {
7978c94862eSPeter Brune     if (next) {
7989566063dSJacob Faibussowitsch       PetscCall(SNESFASCoarseCorrection(snes,X,F,X));
7999566063dSJacob Faibussowitsch       PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm));
8008c94862eSPeter Brune     }
8018c94862eSPeter Brune   }
8028c94862eSPeter Brune   PetscFunctionReturn(0);
8038c94862eSPeter Brune }
8048c94862eSPeter Brune 
80540244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
80634d65b3cSPeter Brune {
80734d65b3cSPeter Brune   Vec            F,B;
80834d65b3cSPeter Brune   SNES           next;
80934d65b3cSPeter Brune 
81034d65b3cSPeter Brune   PetscFunctionBegin;
81134d65b3cSPeter Brune   F = snes->vec_func;
81234d65b3cSPeter Brune   B = snes->vec_rhs;
8139566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleGetCorrection(snes,&next));
81434d65b3cSPeter Brune   if (next) {
8159566063dSJacob Faibussowitsch     PetscCall(SNESFASCoarseCorrection(snes,X,F,X));
8169566063dSJacob Faibussowitsch     PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm));
81734d65b3cSPeter Brune   } else {
8189566063dSJacob Faibussowitsch     PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm));
81934d65b3cSPeter Brune   }
82034d65b3cSPeter Brune   PetscFunctionReturn(0);
82134d65b3cSPeter Brune }
82234d65b3cSPeter Brune 
823fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE;
824fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
825fffbeea8SBarry Smith                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
826fffbeea8SBarry Smith                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
827fffbeea8SBarry Smith                             "  year = 2013,\n"
828fffbeea8SBarry Smith                             "  type = Preprint,\n"
829fffbeea8SBarry Smith                             "  number = {ANL/MCS-P2010-0112},\n"
830fffbeea8SBarry Smith                             "  institution = {Argonne National Laboratory}\n}\n";
831fffbeea8SBarry Smith 
83240244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes)
833421d9b32SPeter Brune {
834f833ba53SLisandro Dalcin   PetscInt       i;
835ddb5aff1SPeter Brune   Vec            X, F;
836fe6f9142SPeter Brune   PetscReal      fnorm;
837b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
838b17ce1afSJed Brown   DM             dm;
839e70c42e5SPeter Brune   PetscBool      isFine;
840b17ce1afSJed Brown 
841421d9b32SPeter Brune   PetscFunctionBegin;
8420b121fc5SBarry 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);
843c579b300SPatrick Farrell 
8449566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite));
845fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
846fa9694d7SPeter Brune   X            = snes->vec_sol;
847f5a6d4f9SBarry Smith   F            = snes->vec_func;
848293a7e31SPeter Brune 
8499566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
850293a7e31SPeter Brune   /* norm setup */
8519566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
852fe6f9142SPeter Brune   snes->iter = 0;
853f833ba53SLisandro Dalcin   snes->norm = 0;
8549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
855e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
8569566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual,snes,0,0,0));
8579566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes,X,F));
8589566063dSJacob Faibussowitsch     if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual,snes,0,0,0));
8591aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
860e4ed7901SPeter Brune 
8619566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
862422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
8639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
864fe6f9142SPeter Brune   snes->norm = fnorm;
8659566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
8669566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes,fnorm,0));
8679566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes,snes->iter,fnorm));
868fe6f9142SPeter Brune 
869fe6f9142SPeter Brune   /* test convergence */
870*dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes,converged ,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
871fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
872e4ed7901SPeter Brune 
873b9c2fdf1SPeter Brune   if (isFine) {
874b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
8759566063dSJacob Faibussowitsch     PetscCall(SNESGetDM(snes,&dm));
876b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
877b17ce1afSJed Brown       DM dmcoarse;
8789566063dSJacob Faibussowitsch       PetscCall(SNESGetDM(ffas->next,&dmcoarse));
8799566063dSJacob Faibussowitsch       PetscCall(DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse));
880b17ce1afSJed Brown       dm   = dmcoarse;
881b17ce1afSJed Brown     }
882b9c2fdf1SPeter Brune   }
883b17ce1afSJed Brown 
884f833ba53SLisandro Dalcin   for (i = 0; i < snes->max_its; i++) {
885fe6f9142SPeter Brune     /* Call general purpose update function */
886*dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes,update, snes->iter);
887f833ba53SLisandro Dalcin 
88807144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
8899566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Multiplicative(snes, X));
8908c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
8919566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Additive(snes, X));
8928c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
8939566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Full(snes, X));
89434d65b3cSPeter Brune     } else if (fas->fastype == SNES_FAS_KASKADE) {
8959566063dSJacob Faibussowitsch       PetscCall(SNESFASCycle_Kaskade(snes, X));
8966c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");
897742fe5e2SPeter Brune 
898742fe5e2SPeter Brune     /* check for FAS cycle divergence */
8991aa26658SKarl Rupp     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
900b9c2fdf1SPeter Brune 
901c90fad12SPeter Brune     /* Monitor convergence */
9029566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
903c90fad12SPeter Brune     snes->iter = i+1;
9049566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
9059566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
9069566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
907c90fad12SPeter Brune     /* Test for convergence */
90866585501SPeter Brune     if (isFine) {
909*dbbe0bcdSBarry Smith       PetscUseTypeMethod(snes,converged ,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);
910c90fad12SPeter Brune       if (snes->reason) break;
911fe6f9142SPeter Brune     }
91266585501SPeter Brune   }
913f833ba53SLisandro Dalcin   if (i == snes->max_its) {
91463a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", i));
915fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
916fe6f9142SPeter Brune   }
917421d9b32SPeter Brune   PetscFunctionReturn(0);
918421d9b32SPeter Brune }
91940244768SBarry Smith 
92040244768SBarry Smith /*MC
92140244768SBarry Smith 
92240244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
92340244768SBarry Smith 
92440244768SBarry Smith    The nonlinear problem is solved by correction using coarse versions
92540244768SBarry Smith    of the nonlinear problem.  This problem is perturbed so that a projected
92640244768SBarry Smith    solution of the fine problem elicits no correction from the coarse problem.
92740244768SBarry Smith 
92840244768SBarry Smith Options Database:
92940244768SBarry Smith +   -snes_fas_levels -  The number of levels
93040244768SBarry Smith .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
93140244768SBarry Smith .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  Additive or multiplicative cycle
93240244768SBarry Smith .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
93340244768SBarry Smith .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
93440244768SBarry Smith .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
93540244768SBarry Smith .   -snes_fas_monitor -  Monitor progress of all of the levels
93640244768SBarry Smith .   -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
93740244768SBarry Smith .   -fas_levels_snes_ -  SNES options for all smoothers
93840244768SBarry Smith .   -fas_levels_cycle_snes_ -  SNES options for all cycles
93940244768SBarry Smith .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
94040244768SBarry Smith .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
94140244768SBarry Smith -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
94240244768SBarry Smith 
94340244768SBarry Smith Notes:
94440244768SBarry Smith    The organization of the FAS solver is slightly different from the organization of PCMG
94540244768SBarry Smith    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
94640244768SBarry Smith    The cycle SNES instance may be used for monitoring convergence on a particular level.
94740244768SBarry Smith 
94840244768SBarry Smith Level: beginner
94940244768SBarry Smith 
95040244768SBarry Smith    References:
951606c0280SSatish Balay .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
95240244768SBarry Smith    SIAM Review, 57(4), 2015
95340244768SBarry Smith 
954db781477SPatrick Sanan .seealso: `PCMG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`
95540244768SBarry Smith M*/
95640244768SBarry Smith 
95740244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
95840244768SBarry Smith {
95940244768SBarry Smith   SNES_FAS       *fas;
96040244768SBarry Smith 
96140244768SBarry Smith   PetscFunctionBegin;
96240244768SBarry Smith   snes->ops->destroy        = SNESDestroy_FAS;
96340244768SBarry Smith   snes->ops->setup          = SNESSetUp_FAS;
96440244768SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
96540244768SBarry Smith   snes->ops->view           = SNESView_FAS;
96640244768SBarry Smith   snes->ops->solve          = SNESSolve_FAS;
96740244768SBarry Smith   snes->ops->reset          = SNESReset_FAS;
96840244768SBarry Smith 
96940244768SBarry Smith   snes->usesksp = PETSC_FALSE;
970efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
97140244768SBarry Smith 
97240244768SBarry Smith   if (!snes->tolerancesset) {
97340244768SBarry Smith     snes->max_funcs = 30000;
97440244768SBarry Smith     snes->max_its   = 10000;
97540244768SBarry Smith   }
97640244768SBarry Smith 
9774fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
9784fc747eaSLawrence Mitchell 
9799566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(snes,&fas));
98040244768SBarry Smith 
98140244768SBarry Smith   snes->data                  = (void*) fas;
98240244768SBarry Smith   fas->level                  = 0;
98340244768SBarry Smith   fas->levels                 = 1;
98440244768SBarry Smith   fas->n_cycles               = 1;
98540244768SBarry Smith   fas->max_up_it              = 1;
98640244768SBarry Smith   fas->max_down_it            = 1;
98740244768SBarry Smith   fas->smoothu                = NULL;
98840244768SBarry Smith   fas->smoothd                = NULL;
98940244768SBarry Smith   fas->next                   = NULL;
99040244768SBarry Smith   fas->previous               = NULL;
99140244768SBarry Smith   fas->fine                   = snes;
99240244768SBarry Smith   fas->interpolate            = NULL;
99340244768SBarry Smith   fas->restrct                = NULL;
99440244768SBarry Smith   fas->inject                 = NULL;
99540244768SBarry Smith   fas->usedmfornumberoflevels = PETSC_FALSE;
99640244768SBarry Smith   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
99740244768SBarry Smith   fas->full_downsweep         = PETSC_FALSE;
9986dfbafefSToby Isaac   fas->full_total             = PETSC_FALSE;
99940244768SBarry Smith 
100040244768SBarry Smith   fas->eventsmoothsetup    = 0;
100140244768SBarry Smith   fas->eventsmoothsolve    = 0;
100240244768SBarry Smith   fas->eventresidual       = 0;
100340244768SBarry Smith   fas->eventinterprestrict = 0;
100440244768SBarry Smith   PetscFunctionReturn(0);
100540244768SBarry Smith }
1006