xref: /petsc/src/snes/impls/fas/fas.c (revision 9e5d08929599755f0fc0004c3575ee121dac326a)
1421d9b32SPeter Brune /* Defines the basic SNES object */
2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnes.h"  I*/
3421d9b32SPeter Brune 
4*9e5d0892SLisandro 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;
9f833ba53SLisandro Dalcin   PetscErrorCode ierr;
10421d9b32SPeter Brune 
11421d9b32SPeter Brune   PetscFunctionBegin;
12ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
13ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
143dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
15bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
16bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
17bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
187cd33a7bSPeter Brune   ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr);
197cd33a7bSPeter Brune   ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr);
20d477d801SPeter Brune   if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);}
21421d9b32SPeter Brune   PetscFunctionReturn(0);
22421d9b32SPeter Brune }
23421d9b32SPeter Brune 
2440244768SBarry Smith static PetscErrorCode SNESDestroy_FAS(SNES snes)
25421d9b32SPeter Brune {
26421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
27f833ba53SLisandro Dalcin   PetscErrorCode ierr;
28421d9b32SPeter Brune 
29421d9b32SPeter Brune   PetscFunctionBegin;
30421d9b32SPeter Brune   /* recursively resets and then destroys */
31f833ba53SLisandro Dalcin   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
321aa26658SKarl Rupp   ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
33421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
34421d9b32SPeter Brune   PetscFunctionReturn(0);
35421d9b32SPeter Brune }
36421d9b32SPeter Brune 
37f833ba53SLisandro Dalcin static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth)
38f833ba53SLisandro Dalcin {
39f833ba53SLisandro Dalcin   SNESLineSearch linesearch;
40f833ba53SLisandro Dalcin   SNESLineSearch slinesearch;
41f833ba53SLisandro Dalcin   void           *lsprectx,*lspostctx;
42f833ba53SLisandro Dalcin   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
43f833ba53SLisandro Dalcin   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
44f833ba53SLisandro Dalcin   PetscErrorCode ierr;
45f833ba53SLisandro Dalcin 
46f833ba53SLisandro Dalcin   PetscFunctionBegin;
47f833ba53SLisandro Dalcin   if (!snes->linesearch) PetscFunctionReturn(0);
48f833ba53SLisandro Dalcin   ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
49f833ba53SLisandro Dalcin   ierr = SNESGetLineSearch(smooth,&slinesearch);CHKERRQ(ierr);
50f833ba53SLisandro Dalcin   ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
51f833ba53SLisandro Dalcin   ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
52f833ba53SLisandro Dalcin   ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
53f833ba53SLisandro Dalcin   ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
54f833ba53SLisandro Dalcin   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
55f833ba53SLisandro Dalcin   PetscFunctionReturn(0);
56f833ba53SLisandro Dalcin }
57f833ba53SLisandro Dalcin 
58f833ba53SLisandro Dalcin static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth)
59f833ba53SLisandro Dalcin {
60f833ba53SLisandro Dalcin   SNES_FAS       *fas = (SNES_FAS*) snes->data;
61f833ba53SLisandro Dalcin   PetscErrorCode ierr;
62f833ba53SLisandro Dalcin 
63f833ba53SLisandro Dalcin   PetscFunctionBegin;
64f833ba53SLisandro Dalcin   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth);CHKERRQ(ierr);
65f833ba53SLisandro Dalcin   ierr = SNESSetFromOptions(smooth);CHKERRQ(ierr);
66f833ba53SLisandro Dalcin   ierr = SNESFASSetUpLineSearch_Private(snes, smooth);CHKERRQ(ierr);
67f833ba53SLisandro Dalcin 
68f833ba53SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
69f833ba53SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
70f833ba53SLisandro Dalcin   ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
71f833ba53SLisandro Dalcin   smooth->vec_sol        = snes->vec_sol;
72f833ba53SLisandro Dalcin   smooth->vec_sol_update = snes->vec_sol_update;
73f833ba53SLisandro Dalcin   smooth->vec_func       = snes->vec_func;
74f833ba53SLisandro Dalcin 
75f833ba53SLisandro Dalcin   if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,smooth,0,0,0);CHKERRQ(ierr);}
76f833ba53SLisandro Dalcin   ierr = SNESSetUp(smooth);CHKERRQ(ierr);
77f833ba53SLisandro Dalcin   if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,smooth,0,0,0);CHKERRQ(ierr);}
78f833ba53SLisandro Dalcin   PetscFunctionReturn(0);
79f833ba53SLisandro Dalcin }
80f833ba53SLisandro Dalcin 
8140244768SBarry Smith static PetscErrorCode SNESSetUp_FAS(SNES snes)
82421d9b32SPeter Brune {
8348bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
84421d9b32SPeter Brune   PetscErrorCode ierr;
85d1adcc6fSPeter Brune   PetscInt       dm_levels;
86ab8d36c9SPeter Brune   SNES           next;
87f833ba53SLisandro Dalcin   PetscBool      isFine, hasCreateRestriction, hasCreateInjection;
88eff52c0eSPeter Brune 
896b2b7091SBarry Smith   PetscFunctionBegin;
90ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
91ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
92d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
93d1adcc6fSPeter Brune     dm_levels++;
94cc05f883SPeter Brune     if (dm_levels > fas->levels) {
953dccd265SPeter Brune       /* reset the number of levels */
960298fd71SBarry Smith       ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr);
97cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
98d1adcc6fSPeter Brune     }
99d1adcc6fSPeter Brune   }
100ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
101ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
1023dccd265SPeter Brune 
103fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
104cc05f883SPeter Brune 
105ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
106ab8d36c9SPeter Brune   if (!fas->smoothd) {
107ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
108ab8d36c9SPeter Brune   }
109ab8d36c9SPeter Brune 
11079d9a41aSPeter Brune   if (snes->dm) {
111ab8d36c9SPeter Brune     /* set the smoother DMs properly */
112f833ba53SLisandro Dalcin     if (fas->smoothu) {ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);}
113ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
11479d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
115ab8d36c9SPeter Brune     if (next) {
11679d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
117ab8d36c9SPeter Brune       if (!next->dm) {
118ce94432eSBarry Smith         ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr);
119ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
12079d9a41aSPeter Brune       }
12179d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
12279d9a41aSPeter Brune       if (!fas->interpolate) {
123ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
124bccf9bb3SJed Brown         if (!fas->restrct) {
1250a7266b2SPatrick Farrell           ierr = DMHasCreateRestriction(next->dm, &hasCreateRestriction);CHKERRQ(ierr);
1260a7266b2SPatrick Farrell           /* DM can create restrictions, use that */
1270a7266b2SPatrick Farrell           if (hasCreateRestriction) {
1280a7266b2SPatrick Farrell             ierr = DMCreateRestriction(next->dm, snes->dm, &fas->restrct);CHKERRQ(ierr);
1290a7266b2SPatrick Farrell           } else {
130bccf9bb3SJed Brown             ierr         = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
13179d9a41aSPeter Brune             fas->restrct = fas->interpolate;
13279d9a41aSPeter Brune           }
133bccf9bb3SJed Brown         }
1340a7266b2SPatrick Farrell       }
13579d9a41aSPeter Brune       /* set the injection from the DM */
13679d9a41aSPeter Brune       if (!fas->inject) {
137f833ba53SLisandro Dalcin         ierr = DMHasCreateInjection(next->dm, &hasCreateInjection);CHKERRQ(ierr);
138f833ba53SLisandro Dalcin         if (hasCreateInjection) {
1396dbf9973SLawrence Mitchell           ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr);
14079d9a41aSPeter Brune         }
14179d9a41aSPeter Brune       }
14279d9a41aSPeter Brune     }
14323e68893SLawrence Mitchell   }
144f833ba53SLisandro Dalcin 
14579d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
14679d9a41aSPeter Brune   if (fas->galerkin) {
1471aa26658SKarl Rupp     if (next) {
14825acbd8eSLisandro Dalcin       ierr = SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next);CHKERRQ(ierr);
1491aa26658SKarl Rupp     }
1501aa26658SKarl Rupp     if (fas->smoothd && fas->level != fas->levels - 1) {
15125acbd8eSLisandro Dalcin       ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes);CHKERRQ(ierr);
1521aa26658SKarl Rupp     }
1531aa26658SKarl Rupp     if (fas->smoothu && fas->level != fas->levels - 1) {
15425acbd8eSLisandro Dalcin       ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes);CHKERRQ(ierr);
1551aa26658SKarl Rupp     }
15679d9a41aSPeter Brune   }
15779d9a41aSPeter Brune 
158534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
159534ebe21SPeter Brune   if (fas->smoothd) {
160bc3f2f05SPeter Brune     if (fas->level == 0 && fas->levels != 1) {
161365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
162534ebe21SPeter Brune     } else {
163365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
164534ebe21SPeter Brune     }
165f833ba53SLisandro Dalcin     ierr = SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd);CHKERRQ(ierr);
166534ebe21SPeter Brune   }
167534ebe21SPeter Brune 
168534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
169534ebe21SPeter Brune   if (fas->smoothu) {
170534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
171365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
172534ebe21SPeter Brune     } else {
173365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
174534ebe21SPeter Brune     }
175f833ba53SLisandro Dalcin     ierr = SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu);CHKERRQ(ierr);
176534ebe21SPeter Brune   }
177d06165b7SPeter Brune 
178ab8d36c9SPeter Brune   if (next) {
17979d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
180ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
181ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
1827fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr);
183f833ba53SLisandro Dalcin     ierr = SNESFASSetUpLineSearch_Private(snes, next);CHKERRQ(ierr);
184ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
18579d9a41aSPeter Brune   }
186f833ba53SLisandro Dalcin 
1876273346dSPeter Brune   /* setup FAS work vectors */
1886273346dSPeter Brune   if (fas->galerkin) {
1896273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
1906273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
1916273346dSPeter Brune   }
192421d9b32SPeter Brune   PetscFunctionReturn(0);
193421d9b32SPeter Brune }
194421d9b32SPeter Brune 
19540244768SBarry Smith static PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes)
196421d9b32SPeter Brune {
197ee78dd50SPeter Brune   SNES_FAS       *fas   = (SNES_FAS*) snes->data;
198ee78dd50SPeter Brune   PetscInt       levels = 1;
19987f44e3fSPeter Brune   PetscBool      flg    = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE;
200421d9b32SPeter Brune   PetscErrorCode ierr;
20107144faaSPeter Brune   SNESFASType    fastype;
202fde0ff24SPeter Brune   const char     *optionsprefix;
203f1c6b773SPeter Brune   SNESLineSearch linesearch;
20466585501SPeter Brune   PetscInt       m, n_up, n_down;
205ab8d36c9SPeter Brune   SNES           next;
206ab8d36c9SPeter Brune   PetscBool      isFine;
207421d9b32SPeter Brune 
208421d9b32SPeter Brune   PetscFunctionBegin;
209ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
210e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");CHKERRQ(ierr);
211ee78dd50SPeter Brune 
212ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
213ab8d36c9SPeter Brune   if (isFine) {
214ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
215c732cbdbSBarry Smith     if (!flg && snes->dm) {
216c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
217c732cbdbSBarry Smith       levels++;
218d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
219c732cbdbSBarry Smith     }
2200298fd71SBarry Smith     ierr    = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr);
22107144faaSPeter Brune     fastype = fas->fastype;
22207144faaSPeter Brune     ierr    = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
22307144faaSPeter Brune     if (flg) {
22407144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
22507144faaSPeter Brune     }
226ee78dd50SPeter Brune 
227fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
228ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
229ab8d36c9SPeter Brune     if (flg) {
230ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
231fde0ff24SPeter Brune     }
23287f44e3fSPeter Brune     ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr);
23387f44e3fSPeter Brune     if (flg) {
23487f44e3fSPeter Brune       ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr);
23587f44e3fSPeter Brune     }
236fde0ff24SPeter Brune 
237ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
238ab8d36c9SPeter Brune     if (flg) {
239ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
240ab8d36c9SPeter Brune     }
241ee78dd50SPeter Brune 
242928e959bSPeter Brune     if (fas->fastype == SNES_FAS_FULL) {
243ba672821SBarry Smith       ierr   = PetscOptionsBool("-snes_fas_full_downsweep","Smooth on the initial down sweep for full FAS cycles","SNESFASFullSetDownSweep",fas->full_downsweep,&fas->full_downsweep,&flg);CHKERRQ(ierr);
244ba672821SBarry Smith       if (flg) {ierr = SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);}
245928e959bSPeter Brune     }
246928e959bSPeter Brune 
24766585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
248162d76ddSPeter Brune 
24966585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
250162d76ddSPeter Brune 
251d142ab34SLawrence Mitchell     {
252d142ab34SLawrence Mitchell       PetscViewer viewer;
253d142ab34SLawrence Mitchell       PetscViewerFormat format;
25416413a6aSBarry Smith       ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_fas_monitor",&viewer,&format,&monflg);CHKERRQ(ierr);
255d142ab34SLawrence Mitchell       if (monflg) {
256d142ab34SLawrence Mitchell         PetscViewerAndFormat *vf;
257d142ab34SLawrence Mitchell         ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
258d142ab34SLawrence Mitchell         ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
259d142ab34SLawrence Mitchell         ierr = SNESFASSetMonitor(snes,vf,PETSC_TRUE);CHKERRQ(ierr);
260d142ab34SLawrence Mitchell       }
261d142ab34SLawrence Mitchell     }
2620dd27c6cSPeter Brune     flg    = PETSC_FALSE;
2630dd27c6cSPeter Brune     monflg = PETSC_TRUE;
2640dd27c6cSPeter Brune     ierr   = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
2650dd27c6cSPeter Brune     if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
266ab8d36c9SPeter Brune   }
267ee78dd50SPeter Brune 
268421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
269f833ba53SLisandro Dalcin 
2708cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
271162d76ddSPeter Brune   if (upflg) {
27266585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
273162d76ddSPeter Brune   }
274162d76ddSPeter Brune   if (downflg) {
27566585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
276162d76ddSPeter Brune   }
277eff52c0eSPeter Brune 
2789e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
2799e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
2809e764e56SPeter Brune     if (!snes->linesearch) {
2817601faf0SJed Brown       ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
2821a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
2839e764e56SPeter Brune     }
2849e764e56SPeter Brune   }
2859e764e56SPeter Brune 
286ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
287f833ba53SLisandro Dalcin   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
288ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
289421d9b32SPeter Brune   PetscFunctionReturn(0);
290421d9b32SPeter Brune }
291421d9b32SPeter Brune 
2929804daf3SBarry Smith #include <petscdraw.h>
29340244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
294421d9b32SPeter Brune {
295421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
296656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
297ab8d36c9SPeter Brune   PetscInt       i;
298421d9b32SPeter Brune   PetscErrorCode ierr;
299ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
300421d9b32SPeter Brune 
301421d9b32SPeter Brune   PetscFunctionBegin;
302ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
303ab8d36c9SPeter Brune   if (isFine) {
304251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
305656ede7eSPeter Brune     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
306421d9b32SPeter Brune     if (iascii) {
307efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
308ab8d36c9SPeter Brune       if (fas->galerkin) {
309ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
310421d9b32SPeter Brune       } else {
311ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
312421d9b32SPeter Brune       }
313ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
314ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
315ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
316ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
317ab8d36c9SPeter Brune         if (!i) {
318ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
319421d9b32SPeter Brune         } else {
320ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
321421d9b32SPeter Brune         }
322ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
323166b3ea4SJed Brown         if (smoothd) {
324ab8d36c9SPeter Brune           ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
325166b3ea4SJed Brown         } else {
326166b3ea4SJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
327166b3ea4SJed Brown         }
328ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
329ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
330ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
331ab8d36c9SPeter Brune         } else if (i) {
332ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
333ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
334166b3ea4SJed Brown           if (smoothu) {
335ab8d36c9SPeter Brune             ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
336166b3ea4SJed Brown           } else {
337166b3ea4SJed Brown             ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
338166b3ea4SJed Brown           }
339ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
340ab8d36c9SPeter Brune         }
341ab8d36c9SPeter Brune       }
342656ede7eSPeter Brune     } else if (isdraw) {
343656ede7eSPeter Brune       PetscDraw draw;
344b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
345656ede7eSPeter Brune       SNES_FAS  *curfas = fas;
346656ede7eSPeter Brune       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
347656ede7eSPeter Brune       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
348656ede7eSPeter Brune       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
349656ede7eSPeter Brune       bottom = y - th;
350656ede7eSPeter Brune       while (curfas) {
351b4375e8dSPeter Brune         if (!curfas->smoothu) {
352656ede7eSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
353656ede7eSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
354656ede7eSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
355b4375e8dSPeter Brune         } else {
356b4375e8dSPeter Brune           w    = 0.5*PetscMin(1.0-x,x);
357b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
358b4375e8dSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
359b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
360b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
361b4375e8dSPeter Brune           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
362b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
363b4375e8dSPeter Brune         }
364656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
365656ede7eSPeter Brune         bottom -= 5*th;
3661aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
3670298fd71SBarry Smith         else curfas = NULL;
368656ede7eSPeter Brune       }
369421d9b32SPeter Brune     }
370ab8d36c9SPeter Brune   }
371421d9b32SPeter Brune   PetscFunctionReturn(0);
372421d9b32SPeter Brune }
373421d9b32SPeter Brune 
37439bd7f45SPeter Brune /*
37539bd7f45SPeter Brune Defines the action of the downsmoother
37639bd7f45SPeter Brune  */
37740244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
378b9c2fdf1SPeter Brune {
379f833ba53SLisandro Dalcin   PetscErrorCode      ierr;
380742fe5e2SPeter Brune   SNESConvergedReason reason;
381ab8d36c9SPeter Brune   Vec                 FPC;
382ab8d36c9SPeter Brune   SNES                smoothd;
3836cbb2f26SLawrence Mitchell   PetscBool           flg;
3840dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
3856e111a19SKarl Rupp 
386421d9b32SPeter Brune   PetscFunctionBegin;
387ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
388e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
389217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);}
390ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
391217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);}
392742fe5e2SPeter Brune   /* check convergence reason for the smoother */
393ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
3941ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
395742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
396742fe5e2SPeter Brune     PetscFunctionReturn(0);
397742fe5e2SPeter Brune   }
3986cbb2f26SLawrence Mitchell 
3990298fd71SBarry Smith   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
4006cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothd, &flg);CHKERRQ(ierr);
4016cbb2f26SLawrence Mitchell   if (!flg) {
4026cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothd, X, FPC);CHKERRQ(ierr);
4036cbb2f26SLawrence Mitchell   }
4044b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
40571dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
40639bd7f45SPeter Brune   PetscFunctionReturn(0);
40739bd7f45SPeter Brune }
40839bd7f45SPeter Brune 
40939bd7f45SPeter Brune 
41039bd7f45SPeter Brune /*
41107144faaSPeter Brune Defines the action of the upsmoother
41239bd7f45SPeter Brune  */
41340244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
4140adebc6cSBarry Smith {
415f833ba53SLisandro Dalcin   PetscErrorCode      ierr;
41639bd7f45SPeter Brune   SNESConvergedReason reason;
417ab8d36c9SPeter Brune   Vec                 FPC;
418ab8d36c9SPeter Brune   SNES                smoothu;
4196cbb2f26SLawrence Mitchell   PetscBool           flg;
4200dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
421ab8d36c9SPeter Brune 
4226e111a19SKarl Rupp   PetscFunctionBegin;
423ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
424217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);}
425ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
426217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);}
42739bd7f45SPeter Brune   /* check convergence reason for the smoother */
428ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
4291ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
43039bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
43139bd7f45SPeter Brune     PetscFunctionReturn(0);
43239bd7f45SPeter Brune   }
4330298fd71SBarry Smith   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
4346cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothu, &flg);CHKERRQ(ierr);
4356cbb2f26SLawrence Mitchell   if (!flg) {
4366cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothu, X, FPC);CHKERRQ(ierr);
4376cbb2f26SLawrence Mitchell   }
4384b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
43971dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
44039bd7f45SPeter Brune   PetscFunctionReturn(0);
44139bd7f45SPeter Brune }
44239bd7f45SPeter Brune 
443938e4a01SJed Brown /*@
444938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
445938e4a01SJed Brown 
446938e4a01SJed Brown    Collective
447938e4a01SJed Brown 
448938e4a01SJed Brown    Input Arguments:
449938e4a01SJed Brown .  snes - SNESFAS
450938e4a01SJed Brown 
451938e4a01SJed Brown    Output Arguments:
452938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
453938e4a01SJed Brown 
454938e4a01SJed Brown    Level: developer
455938e4a01SJed Brown 
456938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
457938e4a01SJed Brown @*/
458938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
459938e4a01SJed Brown {
460938e4a01SJed Brown   PetscErrorCode ierr;
461f833ba53SLisandro Dalcin   SNES_FAS       *fas;
462938e4a01SJed Brown 
463938e4a01SJed Brown   PetscFunctionBegin;
464f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
465f833ba53SLisandro Dalcin   PetscValidPointer(Xcoarse,2);
466f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
4671aa26658SKarl Rupp   if (fas->rscale) {
4681aa26658SKarl Rupp     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
469f5af7f23SKarl Rupp   } else if (fas->inject) {
4702a7a6963SBarry Smith     ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
47113903a91SSatish Balay   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");
472938e4a01SJed Brown   PetscFunctionReturn(0);
473938e4a01SJed Brown }
474938e4a01SJed Brown 
475e9923e8dSJed Brown /*@
476e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
477e9923e8dSJed Brown 
478e9923e8dSJed Brown    Collective
479e9923e8dSJed Brown 
480e9923e8dSJed Brown    Input Arguments:
481e9923e8dSJed Brown +  fine - SNES from which to restrict
482e9923e8dSJed Brown -  Xfine - vector to restrict
483e9923e8dSJed Brown 
484e9923e8dSJed Brown    Output Arguments:
485e9923e8dSJed Brown .  Xcoarse - result of restriction
486e9923e8dSJed Brown 
487e9923e8dSJed Brown    Level: developer
488e9923e8dSJed Brown 
489e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
490e9923e8dSJed Brown @*/
491e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
492e9923e8dSJed Brown {
493e9923e8dSJed Brown   PetscErrorCode ierr;
494f833ba53SLisandro Dalcin   SNES_FAS       *fas;
495e9923e8dSJed Brown 
496e9923e8dSJed Brown   PetscFunctionBegin;
497f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(fine,SNES_CLASSID,1,SNESFAS);
498e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
499e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
500f833ba53SLisandro Dalcin   fas = (SNES_FAS*)fine->data;
501e9923e8dSJed Brown   if (fas->inject) {
502e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
503e9923e8dSJed Brown   } else {
504e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
505e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
506e9923e8dSJed Brown   }
507e9923e8dSJed Brown   PetscFunctionReturn(0);
508e9923e8dSJed Brown }
509e9923e8dSJed Brown 
51039bd7f45SPeter Brune /*
51139bd7f45SPeter Brune 
51239bd7f45SPeter Brune Performs the FAS coarse correction as:
51339bd7f45SPeter Brune 
514b5527d98SMatthew G. Knepley fine problem:   F(x) = b
515b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c
51639bd7f45SPeter Brune 
517b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
51839bd7f45SPeter Brune 
51939bd7f45SPeter Brune  */
5200adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
5210adebc6cSBarry Smith {
52239bd7f45SPeter Brune   PetscErrorCode      ierr;
52339bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
52439bd7f45SPeter Brune   SNESConvergedReason reason;
525ab8d36c9SPeter Brune   SNES                next;
526ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
5270dd27c6cSPeter Brune   SNES_FAS            *fasc;
5285fd66863SKarl Rupp 
52939bd7f45SPeter Brune   PetscFunctionBegin;
530ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
531ab8d36c9SPeter Brune   if (next) {
5320dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
5330dd27c6cSPeter Brune 
534ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
535ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
536ab8d36c9SPeter Brune 
537ab8d36c9SPeter Brune     X_c  = next->vec_sol;
538ab8d36c9SPeter Brune     Xo_c = next->work[0];
539ab8d36c9SPeter Brune     F_c  = next->vec_func;
540ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
541efe1f98aSPeter Brune 
542217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
543938e4a01SJed Brown     ierr = SNESFASRestrict(snes, X, Xo_c);CHKERRQ(ierr);
5445609cbf2SMatthew G. Knepley     /* restrict the defect: R(F(x) - b) */
545ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
546217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
5470dd27c6cSPeter Brune 
548217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
5495609cbf2SMatthew G. Knepley     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
550ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
551217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
5520dd27c6cSPeter Brune 
5530dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
554e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
555b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
556e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
557ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
558ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
559c90fad12SPeter Brune 
560c90fad12SPeter Brune     /* recurse to the next level */
561e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
562ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
563ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
564742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
565742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
566742fe5e2SPeter Brune       PetscFunctionReturn(0);
567742fe5e2SPeter Brune     }
568fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
569fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
5700dd27c6cSPeter Brune 
571217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
572ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
573fb8dc43dSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) X_c, "Coarse correction");CHKERRQ(ierr);
574fb8dc43dSMatthew G. Knepley     ierr = VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view");CHKERRQ(ierr);
575fb8dc43dSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) X_new, "Updated Fine solution");CHKERRQ(ierr);
576fb8dc43dSMatthew G. Knepley     ierr = VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view");CHKERRQ(ierr);
577217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
578293a7e31SPeter Brune   }
57939bd7f45SPeter Brune   PetscFunctionReturn(0);
58039bd7f45SPeter Brune }
58139bd7f45SPeter Brune 
58239bd7f45SPeter Brune /*
58339bd7f45SPeter Brune 
58439bd7f45SPeter Brune The additive cycle looks like:
58539bd7f45SPeter Brune 
58607144faaSPeter Brune xhat = x
58707144faaSPeter Brune xhat = dS(x, b)
58807144faaSPeter Brune x = coarsecorrection(xhat, b_d)
58907144faaSPeter Brune x = x + nu*(xhat - x);
59039bd7f45SPeter Brune (optional) x = uS(x, b)
59139bd7f45SPeter Brune 
59239bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
59339bd7f45SPeter Brune 
59439bd7f45SPeter Brune  */
59540244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
5960adebc6cSBarry Smith {
59707144faaSPeter Brune   Vec                  F, B, Xhat;
59822c1e704SPeter Brune   Vec                  X_c, Xo_c, F_c, B_c;
59939bd7f45SPeter Brune   PetscErrorCode       ierr;
60007144faaSPeter Brune   SNESConvergedReason  reason;
60122c1e704SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
602422a814eSBarry Smith   SNESLineSearchReason lsresult;
603ab8d36c9SPeter Brune   SNES                 next;
604ab8d36c9SPeter Brune   Mat                  restrct, interpolate;
6050dd27c6cSPeter Brune   SNES_FAS             *fas = (SNES_FAS*)snes->data,*fasc;
6060dd27c6cSPeter Brune 
60739bd7f45SPeter Brune   PetscFunctionBegin;
608ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
60939bd7f45SPeter Brune   F    = snes->vec_func;
61039bd7f45SPeter Brune   B    = snes->vec_rhs;
611e7f468e7SPeter Brune   Xhat = snes->work[1];
61207144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
61307144faaSPeter Brune   /* recurse first */
614ab8d36c9SPeter Brune   if (next) {
6150dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
616ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
617ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
618217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
61907144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
620217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
621c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
622ab8d36c9SPeter Brune     X_c  = next->vec_sol;
623ab8d36c9SPeter Brune     Xo_c = next->work[0];
624ab8d36c9SPeter Brune     F_c  = next->vec_func;
625ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
62639bd7f45SPeter Brune 
627938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
62807144faaSPeter Brune     /* restrict the defect */
629ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
63007144faaSPeter Brune 
63107144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
632217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
633ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
634217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
635e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
636b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
637e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
63807144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
63907144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
64007144faaSPeter Brune 
64107144faaSPeter Brune     /* recurse */
642e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
643ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
64407144faaSPeter Brune 
64507144faaSPeter Brune     /* smooth on this level */
64691f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
64707144faaSPeter Brune 
648ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
64907144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
65007144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
65107144faaSPeter Brune       PetscFunctionReturn(0);
65207144faaSPeter Brune     }
65307144faaSPeter Brune 
65407144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
655c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
656ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
65707144faaSPeter Brune 
658ddebd997SPeter Brune     /* additive correction of the coarse direction*/
659f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
660422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
661422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
662422a814eSBarry Smith     if (lsresult) {
6639e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6649e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6659e764e56SPeter Brune         PetscFunctionReturn(0);
6669e764e56SPeter Brune       }
6679e764e56SPeter Brune     }
66807144faaSPeter Brune   } else {
66991f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
67007144faaSPeter Brune   }
67139bd7f45SPeter Brune   PetscFunctionReturn(0);
67239bd7f45SPeter Brune }
67339bd7f45SPeter Brune 
67439bd7f45SPeter Brune /*
67539bd7f45SPeter Brune 
67639bd7f45SPeter Brune Defines the FAS cycle as:
67739bd7f45SPeter Brune 
678b5527d98SMatthew G. Knepley fine problem: F(x) = b
67939bd7f45SPeter Brune coarse problem: F^c(x) = b^c
68039bd7f45SPeter Brune 
681b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
68239bd7f45SPeter Brune 
68339bd7f45SPeter Brune correction:
68439bd7f45SPeter Brune 
68539bd7f45SPeter Brune x = x + I(x^c - Rx)
68639bd7f45SPeter Brune 
68739bd7f45SPeter Brune  */
68840244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
6890adebc6cSBarry Smith {
69039bd7f45SPeter Brune 
69139bd7f45SPeter Brune   PetscErrorCode ierr;
69239bd7f45SPeter Brune   Vec            F,B;
69334d65b3cSPeter Brune   SNES           next;
69439bd7f45SPeter Brune 
69539bd7f45SPeter Brune   PetscFunctionBegin;
69639bd7f45SPeter Brune   F = snes->vec_func;
69739bd7f45SPeter Brune   B = snes->vec_rhs;
69839bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
69934d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
70091f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
70134d65b3cSPeter Brune   if (next) {
7028c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
70391f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
704fe6f9142SPeter Brune   }
705fa9694d7SPeter Brune   PetscFunctionReturn(0);
706421d9b32SPeter Brune }
707421d9b32SPeter Brune 
70840244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
7098c94862eSPeter Brune {
7108c94862eSPeter Brune   SNES           next;
7118c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7128c94862eSPeter Brune   PetscBool      isFine;
7138c94862eSPeter Brune   PetscErrorCode ierr;
7148c94862eSPeter Brune 
7158c94862eSPeter Brune   PetscFunctionBegin;
7168c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7178c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7188c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
7198c94862eSPeter Brune   fas->full_stage = 0;
7208c94862eSPeter Brune   if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);}
7218c94862eSPeter Brune   PetscFunctionReturn(0);
7228c94862eSPeter Brune }
7238c94862eSPeter Brune 
72440244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
7258c94862eSPeter Brune {
7268c94862eSPeter Brune   PetscErrorCode ierr;
7278c94862eSPeter Brune   Vec            F,B;
7288c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7298c94862eSPeter Brune   PetscBool      isFine;
7308c94862eSPeter Brune   SNES           next;
7318c94862eSPeter Brune 
7328c94862eSPeter Brune   PetscFunctionBegin;
7338c94862eSPeter Brune   F = snes->vec_func;
7348c94862eSPeter Brune   B = snes->vec_rhs;
7358c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7368c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
7378c94862eSPeter Brune 
7388c94862eSPeter Brune   if (isFine) {
7398c94862eSPeter Brune     ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr);
7408c94862eSPeter Brune   }
7418c94862eSPeter Brune 
7428c94862eSPeter Brune   if (fas->full_stage == 0) {
743928e959bSPeter Brune     /* downsweep */
7448c94862eSPeter Brune     if (next) {
7458c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
746e140b65aSPatrick Farrell       if (fas->full_downsweep) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
747e140b65aSPatrick Farrell       fas->full_downsweep = PETSC_TRUE;
7488c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7498c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7508c94862eSPeter Brune       if (fas->level != 1) next->max_its -= 1;
7518c94862eSPeter Brune     } else {
752a3a80b83SMatthew G. Knepley       /* The smoother on the coarse level is the coarse solver */
7538c94862eSPeter Brune       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7548c94862eSPeter Brune     }
7558c94862eSPeter Brune     fas->full_stage = 1;
7568c94862eSPeter Brune   } else if (fas->full_stage == 1) {
7578c94862eSPeter Brune     if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
7588c94862eSPeter Brune     if (next) {
7598c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7608c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7618c94862eSPeter Brune     }
7628c94862eSPeter Brune   }
7638c94862eSPeter Brune   /* final v-cycle */
7648c94862eSPeter Brune   if (isFine) {
7658c94862eSPeter Brune     if (next) {
7668c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7678c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7688c94862eSPeter Brune     }
7698c94862eSPeter Brune   }
7708c94862eSPeter Brune   PetscFunctionReturn(0);
7718c94862eSPeter Brune }
7728c94862eSPeter Brune 
77340244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
77434d65b3cSPeter Brune {
77534d65b3cSPeter Brune   PetscErrorCode ierr;
77634d65b3cSPeter Brune   Vec            F,B;
77734d65b3cSPeter Brune   SNES           next;
77834d65b3cSPeter Brune 
77934d65b3cSPeter Brune   PetscFunctionBegin;
78034d65b3cSPeter Brune   F = snes->vec_func;
78134d65b3cSPeter Brune   B = snes->vec_rhs;
78234d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
78334d65b3cSPeter Brune   if (next) {
78434d65b3cSPeter Brune     ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
78534d65b3cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
78634d65b3cSPeter Brune   } else {
78734d65b3cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
78834d65b3cSPeter Brune   }
78934d65b3cSPeter Brune   PetscFunctionReturn(0);
79034d65b3cSPeter Brune }
79134d65b3cSPeter Brune 
792fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE;
793fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
794fffbeea8SBarry Smith                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
795fffbeea8SBarry Smith                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
796fffbeea8SBarry Smith                             "  year = 2013,\n"
797fffbeea8SBarry Smith                             "  type = Preprint,\n"
798fffbeea8SBarry Smith                             "  number = {ANL/MCS-P2010-0112},\n"
799fffbeea8SBarry Smith                             "  institution = {Argonne National Laboratory}\n}\n";
800fffbeea8SBarry Smith 
80140244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes)
802421d9b32SPeter Brune {
803fa9694d7SPeter Brune   PetscErrorCode ierr;
804f833ba53SLisandro Dalcin   PetscInt       i;
805ddb5aff1SPeter Brune   Vec            X, F;
806fe6f9142SPeter Brune   PetscReal      fnorm;
807b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
808b17ce1afSJed Brown   DM             dm;
809e70c42e5SPeter Brune   PetscBool      isFine;
810b17ce1afSJed Brown 
811421d9b32SPeter Brune   PetscFunctionBegin;
8126c4ed002SBarry Smith   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
813c579b300SPatrick Farrell 
814fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
815fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
816fa9694d7SPeter Brune   X            = snes->vec_sol;
817f5a6d4f9SBarry Smith   F            = snes->vec_func;
818293a7e31SPeter Brune 
81918a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
820293a7e31SPeter Brune   /* norm setup */
821e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
822fe6f9142SPeter Brune   snes->iter = 0;
823f833ba53SLisandro Dalcin   snes->norm = 0;
824e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
825e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
826217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
827fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
828217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
8291aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
830e4ed7901SPeter Brune 
831fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
832422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
833e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
834fe6f9142SPeter Brune   snes->norm = fnorm;
835e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
836a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
837f833ba53SLisandro Dalcin   ierr       = SNESMonitor(snes,snes->iter,fnorm);CHKERRQ(ierr);
838fe6f9142SPeter Brune 
839fe6f9142SPeter Brune   /* test convergence */
840fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
841fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
842e4ed7901SPeter Brune 
843b17ce1afSJed Brown 
844b9c2fdf1SPeter Brune   if (isFine) {
845b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
846b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
847b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
848b17ce1afSJed Brown       DM dmcoarse;
849b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
850b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
851b17ce1afSJed Brown       dm   = dmcoarse;
852b17ce1afSJed Brown     }
853b9c2fdf1SPeter Brune   }
854b17ce1afSJed Brown 
855f833ba53SLisandro Dalcin   for (i = 0; i < snes->max_its; i++) {
856fe6f9142SPeter Brune     /* Call general purpose update function */
857fe6f9142SPeter Brune     if (snes->ops->update) {
858fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
859fe6f9142SPeter Brune     }
860f833ba53SLisandro Dalcin 
86107144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
86291f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
8638c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
86491f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
8658c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
8668c94862eSPeter Brune       ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr);
86734d65b3cSPeter Brune     } else if (fas->fastype == SNES_FAS_KASKADE) {
86834d65b3cSPeter Brune       ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr);
8696c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");
870742fe5e2SPeter Brune 
871742fe5e2SPeter Brune     /* check for FAS cycle divergence */
8721aa26658SKarl Rupp     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
873b9c2fdf1SPeter Brune 
874c90fad12SPeter Brune     /* Monitor convergence */
875e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
876c90fad12SPeter Brune     snes->iter = i+1;
877e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
878a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
879c90fad12SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
880c90fad12SPeter Brune     /* Test for convergence */
88166585501SPeter Brune     if (isFine) {
882b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
883c90fad12SPeter Brune       if (snes->reason) break;
884fe6f9142SPeter Brune     }
88566585501SPeter Brune   }
886f833ba53SLisandro Dalcin   if (i == snes->max_its) {
887f833ba53SLisandro Dalcin     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", i);CHKERRQ(ierr);
888fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
889fe6f9142SPeter Brune   }
890421d9b32SPeter Brune   PetscFunctionReturn(0);
891421d9b32SPeter Brune }
89240244768SBarry Smith 
89340244768SBarry Smith /*MC
89440244768SBarry Smith 
89540244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
89640244768SBarry Smith 
89740244768SBarry Smith    The nonlinear problem is solved by correction using coarse versions
89840244768SBarry Smith    of the nonlinear problem.  This problem is perturbed so that a projected
89940244768SBarry Smith    solution of the fine problem elicits no correction from the coarse problem.
90040244768SBarry Smith 
90140244768SBarry Smith Options Database:
90240244768SBarry Smith +   -snes_fas_levels -  The number of levels
90340244768SBarry Smith .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
90440244768SBarry Smith .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  Additive or multiplicative cycle
90540244768SBarry Smith .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
90640244768SBarry Smith .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
90740244768SBarry Smith .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
90840244768SBarry Smith .   -snes_fas_monitor -  Monitor progress of all of the levels
90940244768SBarry Smith .   -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
91040244768SBarry Smith .   -fas_levels_snes_ -  SNES options for all smoothers
91140244768SBarry Smith .   -fas_levels_cycle_snes_ -  SNES options for all cycles
91240244768SBarry Smith .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
91340244768SBarry Smith .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
91440244768SBarry Smith -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
91540244768SBarry Smith 
91640244768SBarry Smith Notes:
91740244768SBarry Smith    The organization of the FAS solver is slightly different from the organization of PCMG
91840244768SBarry Smith    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
91940244768SBarry Smith    The cycle SNES instance may be used for monitoring convergence on a particular level.
92040244768SBarry Smith 
92140244768SBarry Smith Level: beginner
92240244768SBarry Smith 
92340244768SBarry Smith    References:
92440244768SBarry Smith . 1. -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
92540244768SBarry Smith    SIAM Review, 57(4), 2015
92640244768SBarry Smith 
92740244768SBarry Smith .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
92840244768SBarry Smith M*/
92940244768SBarry Smith 
93040244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
93140244768SBarry Smith {
93240244768SBarry Smith   SNES_FAS       *fas;
93340244768SBarry Smith   PetscErrorCode ierr;
93440244768SBarry Smith 
93540244768SBarry Smith   PetscFunctionBegin;
93640244768SBarry Smith   snes->ops->destroy        = SNESDestroy_FAS;
93740244768SBarry Smith   snes->ops->setup          = SNESSetUp_FAS;
93840244768SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
93940244768SBarry Smith   snes->ops->view           = SNESView_FAS;
94040244768SBarry Smith   snes->ops->solve          = SNESSolve_FAS;
94140244768SBarry Smith   snes->ops->reset          = SNESReset_FAS;
94240244768SBarry Smith 
94340244768SBarry Smith   snes->usesksp = PETSC_FALSE;
944efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
94540244768SBarry Smith 
94640244768SBarry Smith   if (!snes->tolerancesset) {
94740244768SBarry Smith     snes->max_funcs = 30000;
94840244768SBarry Smith     snes->max_its   = 10000;
94940244768SBarry Smith   }
95040244768SBarry Smith 
9514fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
9524fc747eaSLawrence Mitchell 
95340244768SBarry Smith   ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr);
95440244768SBarry Smith 
95540244768SBarry Smith   snes->data                  = (void*) fas;
95640244768SBarry Smith   fas->level                  = 0;
95740244768SBarry Smith   fas->levels                 = 1;
95840244768SBarry Smith   fas->n_cycles               = 1;
95940244768SBarry Smith   fas->max_up_it              = 1;
96040244768SBarry Smith   fas->max_down_it            = 1;
96140244768SBarry Smith   fas->smoothu                = NULL;
96240244768SBarry Smith   fas->smoothd                = NULL;
96340244768SBarry Smith   fas->next                   = NULL;
96440244768SBarry Smith   fas->previous               = NULL;
96540244768SBarry Smith   fas->fine                   = snes;
96640244768SBarry Smith   fas->interpolate            = NULL;
96740244768SBarry Smith   fas->restrct                = NULL;
96840244768SBarry Smith   fas->inject                 = NULL;
96940244768SBarry Smith   fas->usedmfornumberoflevels = PETSC_FALSE;
97040244768SBarry Smith   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
97140244768SBarry Smith   fas->full_downsweep         = PETSC_FALSE;
97240244768SBarry Smith 
97340244768SBarry Smith   fas->eventsmoothsetup    = 0;
97440244768SBarry Smith   fas->eventsmoothsolve    = 0;
97540244768SBarry Smith   fas->eventresidual       = 0;
97640244768SBarry Smith   fas->eventinterprestrict = 0;
97740244768SBarry Smith   PetscFunctionReturn(0);
97840244768SBarry Smith }
979