xref: /petsc/src/snes/impls/fas/fas.c (revision 60d4fc614ff537ecd43db92193c85cc85809f93f)
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;
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);}
2456dfbafefSToby Isaac       ierr   = 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);CHKERRQ(ierr);
2466dfbafefSToby Isaac       if (flg) {ierr = SNESFASFullSetTotal(snes,fas->full_total);CHKERRQ(ierr);}
247928e959bSPeter Brune     }
248928e959bSPeter Brune 
24966585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
250162d76ddSPeter Brune 
25166585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
252162d76ddSPeter Brune 
253d142ab34SLawrence Mitchell     {
254d142ab34SLawrence Mitchell       PetscViewer viewer;
255d142ab34SLawrence Mitchell       PetscViewerFormat format;
25616413a6aSBarry Smith       ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_fas_monitor",&viewer,&format,&monflg);CHKERRQ(ierr);
257d142ab34SLawrence Mitchell       if (monflg) {
258d142ab34SLawrence Mitchell         PetscViewerAndFormat *vf;
259d142ab34SLawrence Mitchell         ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
260d142ab34SLawrence Mitchell         ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
261d142ab34SLawrence Mitchell         ierr = SNESFASSetMonitor(snes,vf,PETSC_TRUE);CHKERRQ(ierr);
262d142ab34SLawrence Mitchell       }
263d142ab34SLawrence Mitchell     }
2640dd27c6cSPeter Brune     flg    = PETSC_FALSE;
2650dd27c6cSPeter Brune     monflg = PETSC_TRUE;
2660dd27c6cSPeter Brune     ierr   = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
2670dd27c6cSPeter Brune     if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
268ab8d36c9SPeter Brune   }
269ee78dd50SPeter Brune 
270421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
271f833ba53SLisandro Dalcin 
2728cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
273162d76ddSPeter Brune   if (upflg) {
27466585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
275162d76ddSPeter Brune   }
276162d76ddSPeter Brune   if (downflg) {
27766585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
278162d76ddSPeter Brune   }
279eff52c0eSPeter Brune 
2809e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
2819e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
2829e764e56SPeter Brune     if (!snes->linesearch) {
2837601faf0SJed Brown       ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
2841a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
2859e764e56SPeter Brune     }
2869e764e56SPeter Brune   }
2879e764e56SPeter Brune 
288ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
289f833ba53SLisandro Dalcin   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
290ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
291421d9b32SPeter Brune   PetscFunctionReturn(0);
292421d9b32SPeter Brune }
293421d9b32SPeter Brune 
2949804daf3SBarry Smith #include <petscdraw.h>
29540244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
296421d9b32SPeter Brune {
297421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
298656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
299ab8d36c9SPeter Brune   PetscInt       i;
300421d9b32SPeter Brune   PetscErrorCode ierr;
301ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
302421d9b32SPeter Brune 
303421d9b32SPeter Brune   PetscFunctionBegin;
304ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
305ab8d36c9SPeter Brune   if (isFine) {
306251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
307656ede7eSPeter Brune     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
308421d9b32SPeter Brune     if (iascii) {
309efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
310ab8d36c9SPeter Brune       if (fas->galerkin) {
311ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
312421d9b32SPeter Brune       } else {
313ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
314421d9b32SPeter Brune       }
315ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
316ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
317ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
318ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
319ab8d36c9SPeter Brune         if (!i) {
320ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
321421d9b32SPeter Brune         } else {
322ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
323421d9b32SPeter Brune         }
324ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
325166b3ea4SJed Brown         if (smoothd) {
326ab8d36c9SPeter Brune           ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
327166b3ea4SJed Brown         } else {
328166b3ea4SJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
329166b3ea4SJed Brown         }
330ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
331ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
332ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
333ab8d36c9SPeter Brune         } else if (i) {
334ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
335ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
336166b3ea4SJed Brown           if (smoothu) {
337ab8d36c9SPeter Brune             ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
338166b3ea4SJed Brown           } else {
339166b3ea4SJed Brown             ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
340166b3ea4SJed Brown           }
341ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
342ab8d36c9SPeter Brune         }
343ab8d36c9SPeter Brune       }
344656ede7eSPeter Brune     } else if (isdraw) {
345656ede7eSPeter Brune       PetscDraw draw;
346b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
347656ede7eSPeter Brune       SNES_FAS  *curfas = fas;
348656ede7eSPeter Brune       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
349656ede7eSPeter Brune       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
350656ede7eSPeter Brune       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
351656ede7eSPeter Brune       bottom = y - th;
352656ede7eSPeter Brune       while (curfas) {
353b4375e8dSPeter Brune         if (!curfas->smoothu) {
354656ede7eSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
355*60d4fc61SSatish Balay           if (curfas->smoothd) {ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);}
356656ede7eSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
357b4375e8dSPeter Brune         } else {
358b4375e8dSPeter Brune           w    = 0.5*PetscMin(1.0-x,x);
359b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
360*60d4fc61SSatish Balay           if (curfas->smoothd) {ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);}
361b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
362b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
363*60d4fc61SSatish Balay           if (curfas->smoothu) {ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);}
364b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
365b4375e8dSPeter Brune         }
366656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
367656ede7eSPeter Brune         bottom -= 5*th;
3681aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
3690298fd71SBarry Smith         else curfas = NULL;
370656ede7eSPeter Brune       }
371421d9b32SPeter Brune     }
372ab8d36c9SPeter Brune   }
373421d9b32SPeter Brune   PetscFunctionReturn(0);
374421d9b32SPeter Brune }
375421d9b32SPeter Brune 
37639bd7f45SPeter Brune /*
37739bd7f45SPeter Brune Defines the action of the downsmoother
37839bd7f45SPeter Brune  */
37940244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
380b9c2fdf1SPeter Brune {
381f833ba53SLisandro Dalcin   PetscErrorCode      ierr;
382742fe5e2SPeter Brune   SNESConvergedReason reason;
383ab8d36c9SPeter Brune   Vec                 FPC;
384ab8d36c9SPeter Brune   SNES                smoothd;
3856cbb2f26SLawrence Mitchell   PetscBool           flg;
3860dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
3876e111a19SKarl Rupp 
388421d9b32SPeter Brune   PetscFunctionBegin;
389ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
390e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
391217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);}
392ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
393217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);}
394742fe5e2SPeter Brune   /* check convergence reason for the smoother */
395ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
3961ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
397742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
398742fe5e2SPeter Brune     PetscFunctionReturn(0);
399742fe5e2SPeter Brune   }
4006cbb2f26SLawrence Mitchell 
4010298fd71SBarry Smith   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
4026cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothd, &flg);CHKERRQ(ierr);
4036cbb2f26SLawrence Mitchell   if (!flg) {
4046cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothd, X, FPC);CHKERRQ(ierr);
4056cbb2f26SLawrence Mitchell   }
4064b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
40771dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
40839bd7f45SPeter Brune   PetscFunctionReturn(0);
40939bd7f45SPeter Brune }
41039bd7f45SPeter Brune 
41139bd7f45SPeter Brune /*
41207144faaSPeter Brune Defines the action of the upsmoother
41339bd7f45SPeter Brune  */
41440244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
4150adebc6cSBarry Smith {
416f833ba53SLisandro Dalcin   PetscErrorCode      ierr;
41739bd7f45SPeter Brune   SNESConvergedReason reason;
418ab8d36c9SPeter Brune   Vec                 FPC;
419ab8d36c9SPeter Brune   SNES                smoothu;
4206cbb2f26SLawrence Mitchell   PetscBool           flg;
4210dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
422ab8d36c9SPeter Brune 
4236e111a19SKarl Rupp   PetscFunctionBegin;
424ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
425217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);}
426ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
427217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);}
42839bd7f45SPeter Brune   /* check convergence reason for the smoother */
429ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
4301ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
43139bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
43239bd7f45SPeter Brune     PetscFunctionReturn(0);
43339bd7f45SPeter Brune   }
4340298fd71SBarry Smith   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
4356cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothu, &flg);CHKERRQ(ierr);
4366cbb2f26SLawrence Mitchell   if (!flg) {
4376cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothu, X, FPC);CHKERRQ(ierr);
4386cbb2f26SLawrence Mitchell   }
4394b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
44071dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
44139bd7f45SPeter Brune   PetscFunctionReturn(0);
44239bd7f45SPeter Brune }
44339bd7f45SPeter Brune 
444938e4a01SJed Brown /*@
445938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
446938e4a01SJed Brown 
447938e4a01SJed Brown    Collective
448938e4a01SJed Brown 
449938e4a01SJed Brown    Input Arguments:
450938e4a01SJed Brown .  snes - SNESFAS
451938e4a01SJed Brown 
452938e4a01SJed Brown    Output Arguments:
453938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
454938e4a01SJed Brown 
455938e4a01SJed Brown    Level: developer
456938e4a01SJed Brown 
457938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
458938e4a01SJed Brown @*/
459938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
460938e4a01SJed Brown {
461938e4a01SJed Brown   PetscErrorCode ierr;
462f833ba53SLisandro Dalcin   SNES_FAS       *fas;
463938e4a01SJed Brown 
464938e4a01SJed Brown   PetscFunctionBegin;
465f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
466f833ba53SLisandro Dalcin   PetscValidPointer(Xcoarse,2);
467f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
4681aa26658SKarl Rupp   if (fas->rscale) {
4691aa26658SKarl Rupp     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
470f5af7f23SKarl Rupp   } else if (fas->inject) {
4712a7a6963SBarry Smith     ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
47213903a91SSatish Balay   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");
473938e4a01SJed Brown   PetscFunctionReturn(0);
474938e4a01SJed Brown }
475938e4a01SJed Brown 
476e9923e8dSJed Brown /*@
477e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
478e9923e8dSJed Brown 
479e9923e8dSJed Brown    Collective
480e9923e8dSJed Brown 
481e9923e8dSJed Brown    Input Arguments:
482e9923e8dSJed Brown +  fine - SNES from which to restrict
483e9923e8dSJed Brown -  Xfine - vector to restrict
484e9923e8dSJed Brown 
485e9923e8dSJed Brown    Output Arguments:
486e9923e8dSJed Brown .  Xcoarse - result of restriction
487e9923e8dSJed Brown 
488e9923e8dSJed Brown    Level: developer
489e9923e8dSJed Brown 
490e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
491e9923e8dSJed Brown @*/
492e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
493e9923e8dSJed Brown {
494e9923e8dSJed Brown   PetscErrorCode ierr;
495f833ba53SLisandro Dalcin   SNES_FAS       *fas;
496e9923e8dSJed Brown 
497e9923e8dSJed Brown   PetscFunctionBegin;
498f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(fine,SNES_CLASSID,1,SNESFAS);
499e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
500e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
501f833ba53SLisandro Dalcin   fas = (SNES_FAS*)fine->data;
502e9923e8dSJed Brown   if (fas->inject) {
503e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
504e9923e8dSJed Brown   } else {
505e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
506e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
507e9923e8dSJed Brown   }
508e9923e8dSJed Brown   PetscFunctionReturn(0);
509e9923e8dSJed Brown }
510e9923e8dSJed Brown 
51139bd7f45SPeter Brune /*
51239bd7f45SPeter Brune 
5136dfbafefSToby Isaac Performs a variant of FAS using the interpolated total coarse solution
5146dfbafefSToby Isaac 
5156dfbafefSToby Isaac fine problem:   F(x) = b
5166dfbafefSToby Isaac coarse problem: F^c(x^c) = Rb, Initial guess Rx
5176dfbafefSToby Isaac interpolated solution: x^f = I x^c (total solution interpolation
5186dfbafefSToby Isaac 
5196dfbafefSToby Isaac  */
5206dfbafefSToby Isaac static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new)
5216dfbafefSToby Isaac {
5226dfbafefSToby Isaac   PetscErrorCode      ierr;
5236dfbafefSToby Isaac   Vec                 X_c, B_c;
5246dfbafefSToby Isaac   SNESConvergedReason reason;
5256dfbafefSToby Isaac   SNES                next;
5266dfbafefSToby Isaac   Mat                 restrct, interpolate;
5276dfbafefSToby Isaac   SNES_FAS            *fasc;
5286dfbafefSToby Isaac 
5296dfbafefSToby Isaac   PetscFunctionBegin;
5306dfbafefSToby Isaac   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
5316dfbafefSToby Isaac   if (next) {
5326dfbafefSToby Isaac     fasc = (SNES_FAS*)next->data;
5336dfbafefSToby Isaac 
5346dfbafefSToby Isaac     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
5356dfbafefSToby Isaac     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
5366dfbafefSToby Isaac 
5376dfbafefSToby Isaac     X_c  = next->vec_sol;
5386dfbafefSToby Isaac 
5396dfbafefSToby Isaac     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
5406dfbafefSToby Isaac     /* restrict the total solution: Rb */
5416dfbafefSToby Isaac     ierr = SNESFASRestrict(snes, X, X_c);CHKERRQ(ierr);
5426dfbafefSToby Isaac     B_c = next->vec_rhs;
5436dfbafefSToby Isaac     if (snes->vec_rhs) {
5446dfbafefSToby Isaac       /* restrict the total rhs defect: Rb */
5456dfbafefSToby Isaac       ierr = MatRestrict(restrct, snes->vec_rhs, B_c);CHKERRQ(ierr);
5466dfbafefSToby Isaac     } else {
5476dfbafefSToby Isaac       ierr = VecSet(B_c, 0.);CHKERRQ(ierr);
5486dfbafefSToby Isaac     }
5496dfbafefSToby Isaac     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
5506dfbafefSToby Isaac 
5516dfbafefSToby Isaac     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
5526dfbafefSToby Isaac     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
5536dfbafefSToby Isaac     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
5546dfbafefSToby Isaac       snes->reason = SNES_DIVERGED_INNER;
5556dfbafefSToby Isaac       PetscFunctionReturn(0);
5566dfbafefSToby Isaac     }
5576dfbafefSToby Isaac     /* x^f <- Ix^c*/
5586dfbafefSToby Isaac     DM dmc,dmf;
5596dfbafefSToby Isaac 
5606dfbafefSToby Isaac     ierr = SNESGetDM(next, &dmc);CHKERRQ(ierr);
5616dfbafefSToby Isaac     ierr = SNESGetDM(snes, &dmf);CHKERRQ(ierr);
5626dfbafefSToby Isaac     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
5636dfbafefSToby Isaac     ierr = DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new);CHKERRQ(ierr);
5646dfbafefSToby Isaac     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
5656dfbafefSToby Isaac     ierr = PetscObjectSetName((PetscObject) X_c, "Coarse solution");CHKERRQ(ierr);
5666dfbafefSToby Isaac     ierr = VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view");CHKERRQ(ierr);
5676dfbafefSToby Isaac     ierr = PetscObjectSetName((PetscObject) X_new, "Updated Fine solution");CHKERRQ(ierr);
5686dfbafefSToby Isaac     ierr = VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view");CHKERRQ(ierr);
5696dfbafefSToby Isaac   }
5706dfbafefSToby Isaac   PetscFunctionReturn(0);
5716dfbafefSToby Isaac }
5726dfbafefSToby Isaac 
5736dfbafefSToby Isaac /*
5746dfbafefSToby Isaac 
57539bd7f45SPeter Brune Performs the FAS coarse correction as:
57639bd7f45SPeter Brune 
577b5527d98SMatthew G. Knepley fine problem:   F(x) = b
578b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c
57939bd7f45SPeter Brune 
580b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
58139bd7f45SPeter Brune 
58239bd7f45SPeter Brune  */
5830adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
5840adebc6cSBarry Smith {
58539bd7f45SPeter Brune   PetscErrorCode      ierr;
58639bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
58739bd7f45SPeter Brune   SNESConvergedReason reason;
588ab8d36c9SPeter Brune   SNES                next;
589ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
5900dd27c6cSPeter Brune   SNES_FAS            *fasc;
5915fd66863SKarl Rupp 
59239bd7f45SPeter Brune   PetscFunctionBegin;
593ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
594ab8d36c9SPeter Brune   if (next) {
5950dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
5960dd27c6cSPeter Brune 
597ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
598ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
599ab8d36c9SPeter Brune 
600ab8d36c9SPeter Brune     X_c  = next->vec_sol;
601ab8d36c9SPeter Brune     Xo_c = next->work[0];
602ab8d36c9SPeter Brune     F_c  = next->vec_func;
603ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
604efe1f98aSPeter Brune 
605217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
606938e4a01SJed Brown     ierr = SNESFASRestrict(snes, X, Xo_c);CHKERRQ(ierr);
6075609cbf2SMatthew G. Knepley     /* restrict the defect: R(F(x) - b) */
608ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
609217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
6100dd27c6cSPeter Brune 
611217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
6125609cbf2SMatthew G. Knepley     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
613ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
614217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
6150dd27c6cSPeter Brune 
6160dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
617e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
618b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
619e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
620ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
621ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
622c90fad12SPeter Brune 
623c90fad12SPeter Brune     /* recurse to the next level */
624e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
625ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
626ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
627742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
628742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
629742fe5e2SPeter Brune       PetscFunctionReturn(0);
630742fe5e2SPeter Brune     }
631fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
632fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
6330dd27c6cSPeter Brune 
634217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
635ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
6366dfbafefSToby Isaac     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
637fb8dc43dSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) X_c, "Coarse correction");CHKERRQ(ierr);
638fb8dc43dSMatthew G. Knepley     ierr = VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view");CHKERRQ(ierr);
639fb8dc43dSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) X_new, "Updated Fine solution");CHKERRQ(ierr);
640fb8dc43dSMatthew G. Knepley     ierr = VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view");CHKERRQ(ierr);
641293a7e31SPeter Brune   }
64239bd7f45SPeter Brune   PetscFunctionReturn(0);
64339bd7f45SPeter Brune }
64439bd7f45SPeter Brune 
64539bd7f45SPeter Brune /*
64639bd7f45SPeter Brune 
64739bd7f45SPeter Brune The additive cycle looks like:
64839bd7f45SPeter Brune 
64907144faaSPeter Brune xhat = x
65007144faaSPeter Brune xhat = dS(x, b)
65107144faaSPeter Brune x = coarsecorrection(xhat, b_d)
65207144faaSPeter Brune x = x + nu*(xhat - x);
65339bd7f45SPeter Brune (optional) x = uS(x, b)
65439bd7f45SPeter Brune 
65539bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
65639bd7f45SPeter Brune 
65739bd7f45SPeter Brune  */
65840244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
6590adebc6cSBarry Smith {
66007144faaSPeter Brune   Vec                  F, B, Xhat;
66122c1e704SPeter Brune   Vec                  X_c, Xo_c, F_c, B_c;
66239bd7f45SPeter Brune   PetscErrorCode       ierr;
66307144faaSPeter Brune   SNESConvergedReason  reason;
66422c1e704SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
665422a814eSBarry Smith   SNESLineSearchReason lsresult;
666ab8d36c9SPeter Brune   SNES                 next;
667ab8d36c9SPeter Brune   Mat                  restrct, interpolate;
6680dd27c6cSPeter Brune   SNES_FAS             *fas = (SNES_FAS*)snes->data,*fasc;
6690dd27c6cSPeter Brune 
67039bd7f45SPeter Brune   PetscFunctionBegin;
671ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
67239bd7f45SPeter Brune   F    = snes->vec_func;
67339bd7f45SPeter Brune   B    = snes->vec_rhs;
674e7f468e7SPeter Brune   Xhat = snes->work[1];
67507144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
67607144faaSPeter Brune   /* recurse first */
677ab8d36c9SPeter Brune   if (next) {
6780dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
679ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
680ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
681217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
68207144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
683217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
684c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
685ab8d36c9SPeter Brune     X_c  = next->vec_sol;
686ab8d36c9SPeter Brune     Xo_c = next->work[0];
687ab8d36c9SPeter Brune     F_c  = next->vec_func;
688ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
68939bd7f45SPeter Brune 
690938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
69107144faaSPeter Brune     /* restrict the defect */
692ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
69307144faaSPeter Brune 
69407144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
695217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
696ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
697217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
698e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
699b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
700e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
70107144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
70207144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
70307144faaSPeter Brune 
70407144faaSPeter Brune     /* recurse */
705e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
706ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
70707144faaSPeter Brune 
70807144faaSPeter Brune     /* smooth on this level */
70991f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
71007144faaSPeter Brune 
711ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
71207144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
71307144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
71407144faaSPeter Brune       PetscFunctionReturn(0);
71507144faaSPeter Brune     }
71607144faaSPeter Brune 
71707144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
718c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
719ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
72007144faaSPeter Brune 
721ddebd997SPeter Brune     /* additive correction of the coarse direction*/
722f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
723422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
724422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
725422a814eSBarry Smith     if (lsresult) {
7269e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
7279e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
7289e764e56SPeter Brune         PetscFunctionReturn(0);
7299e764e56SPeter Brune       }
7309e764e56SPeter Brune     }
73107144faaSPeter Brune   } else {
73291f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
73307144faaSPeter Brune   }
73439bd7f45SPeter Brune   PetscFunctionReturn(0);
73539bd7f45SPeter Brune }
73639bd7f45SPeter Brune 
73739bd7f45SPeter Brune /*
73839bd7f45SPeter Brune 
73939bd7f45SPeter Brune Defines the FAS cycle as:
74039bd7f45SPeter Brune 
741b5527d98SMatthew G. Knepley fine problem: F(x) = b
74239bd7f45SPeter Brune coarse problem: F^c(x) = b^c
74339bd7f45SPeter Brune 
744b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
74539bd7f45SPeter Brune 
74639bd7f45SPeter Brune correction:
74739bd7f45SPeter Brune 
74839bd7f45SPeter Brune x = x + I(x^c - Rx)
74939bd7f45SPeter Brune 
75039bd7f45SPeter Brune  */
75140244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
7520adebc6cSBarry Smith {
75339bd7f45SPeter Brune 
75439bd7f45SPeter Brune   PetscErrorCode ierr;
75539bd7f45SPeter Brune   Vec            F,B;
75634d65b3cSPeter Brune   SNES           next;
75739bd7f45SPeter Brune 
75839bd7f45SPeter Brune   PetscFunctionBegin;
75939bd7f45SPeter Brune   F = snes->vec_func;
76039bd7f45SPeter Brune   B = snes->vec_rhs;
76139bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
76234d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
76391f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
76434d65b3cSPeter Brune   if (next) {
7658c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
76691f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
767fe6f9142SPeter Brune   }
768fa9694d7SPeter Brune   PetscFunctionReturn(0);
769421d9b32SPeter Brune }
770421d9b32SPeter Brune 
77140244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
7728c94862eSPeter Brune {
7738c94862eSPeter Brune   SNES           next;
7748c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7758c94862eSPeter Brune   PetscBool      isFine;
7768c94862eSPeter Brune   PetscErrorCode ierr;
7778c94862eSPeter Brune 
7788c94862eSPeter Brune   PetscFunctionBegin;
7798c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7808c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7818c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
7828c94862eSPeter Brune   fas->full_stage = 0;
7838c94862eSPeter Brune   if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);}
7848c94862eSPeter Brune   PetscFunctionReturn(0);
7858c94862eSPeter Brune }
7868c94862eSPeter Brune 
78740244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
7888c94862eSPeter Brune {
7898c94862eSPeter Brune   PetscErrorCode ierr;
7908c94862eSPeter Brune   Vec            F,B;
7918c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7928c94862eSPeter Brune   PetscBool      isFine;
7938c94862eSPeter Brune   SNES           next;
7948c94862eSPeter Brune 
7958c94862eSPeter Brune   PetscFunctionBegin;
7968c94862eSPeter Brune   F = snes->vec_func;
7978c94862eSPeter Brune   B = snes->vec_rhs;
7988c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7998c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
8008c94862eSPeter Brune 
8018c94862eSPeter Brune   if (isFine) {
8028c94862eSPeter Brune     ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr);
8038c94862eSPeter Brune   }
8048c94862eSPeter Brune 
8058c94862eSPeter Brune   if (fas->full_stage == 0) {
806928e959bSPeter Brune     /* downsweep */
8078c94862eSPeter Brune     if (next) {
8088c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
809e140b65aSPatrick Farrell       if (fas->full_downsweep) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
810e140b65aSPatrick Farrell       fas->full_downsweep = PETSC_TRUE;
8116dfbafefSToby Isaac       if (fas->full_total) {ierr = SNESFASInterpolatedCoarseSolution(snes,X,X);CHKERRQ(ierr);}
8126dfbafefSToby Isaac       else                 {ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);}
8136dfbafefSToby Isaac       fas->full_total = PETSC_FALSE;
8148c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8158c94862eSPeter Brune       if (fas->level != 1) next->max_its -= 1;
8168c94862eSPeter Brune     } else {
817a3a80b83SMatthew G. Knepley       /* The smoother on the coarse level is the coarse solver */
8188c94862eSPeter Brune       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8198c94862eSPeter Brune     }
8208c94862eSPeter Brune     fas->full_stage = 1;
8218c94862eSPeter Brune   } else if (fas->full_stage == 1) {
8228c94862eSPeter Brune     if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
8238c94862eSPeter Brune     if (next) {
8248c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8258c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8268c94862eSPeter Brune     }
8278c94862eSPeter Brune   }
8288c94862eSPeter Brune   /* final v-cycle */
8298c94862eSPeter Brune   if (isFine) {
8308c94862eSPeter Brune     if (next) {
8318c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8328c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8338c94862eSPeter Brune     }
8348c94862eSPeter Brune   }
8358c94862eSPeter Brune   PetscFunctionReturn(0);
8368c94862eSPeter Brune }
8378c94862eSPeter Brune 
83840244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
83934d65b3cSPeter Brune {
84034d65b3cSPeter Brune   PetscErrorCode ierr;
84134d65b3cSPeter Brune   Vec            F,B;
84234d65b3cSPeter Brune   SNES           next;
84334d65b3cSPeter Brune 
84434d65b3cSPeter Brune   PetscFunctionBegin;
84534d65b3cSPeter Brune   F = snes->vec_func;
84634d65b3cSPeter Brune   B = snes->vec_rhs;
84734d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
84834d65b3cSPeter Brune   if (next) {
84934d65b3cSPeter Brune     ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
85034d65b3cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
85134d65b3cSPeter Brune   } else {
85234d65b3cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
85334d65b3cSPeter Brune   }
85434d65b3cSPeter Brune   PetscFunctionReturn(0);
85534d65b3cSPeter Brune }
85634d65b3cSPeter Brune 
857fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE;
858fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
859fffbeea8SBarry Smith                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
860fffbeea8SBarry Smith                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
861fffbeea8SBarry Smith                             "  year = 2013,\n"
862fffbeea8SBarry Smith                             "  type = Preprint,\n"
863fffbeea8SBarry Smith                             "  number = {ANL/MCS-P2010-0112},\n"
864fffbeea8SBarry Smith                             "  institution = {Argonne National Laboratory}\n}\n";
865fffbeea8SBarry Smith 
86640244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes)
867421d9b32SPeter Brune {
868fa9694d7SPeter Brune   PetscErrorCode ierr;
869f833ba53SLisandro Dalcin   PetscInt       i;
870ddb5aff1SPeter Brune   Vec            X, F;
871fe6f9142SPeter Brune   PetscReal      fnorm;
872b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
873b17ce1afSJed Brown   DM             dm;
874e70c42e5SPeter Brune   PetscBool      isFine;
875b17ce1afSJed Brown 
876421d9b32SPeter Brune   PetscFunctionBegin;
8776c4ed002SBarry 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);
878c579b300SPatrick Farrell 
879fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
880fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
881fa9694d7SPeter Brune   X            = snes->vec_sol;
882f5a6d4f9SBarry Smith   F            = snes->vec_func;
883293a7e31SPeter Brune 
88418a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
885293a7e31SPeter Brune   /* norm setup */
886e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
887fe6f9142SPeter Brune   snes->iter = 0;
888f833ba53SLisandro Dalcin   snes->norm = 0;
889e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
890e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
891217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
892fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
893217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
8941aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
895e4ed7901SPeter Brune 
896fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
897422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
898e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
899fe6f9142SPeter Brune   snes->norm = fnorm;
900e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
901a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
902f833ba53SLisandro Dalcin   ierr       = SNESMonitor(snes,snes->iter,fnorm);CHKERRQ(ierr);
903fe6f9142SPeter Brune 
904fe6f9142SPeter Brune   /* test convergence */
905fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
906fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
907e4ed7901SPeter Brune 
908b9c2fdf1SPeter Brune   if (isFine) {
909b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
910b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
911b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
912b17ce1afSJed Brown       DM dmcoarse;
913b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
914b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
915b17ce1afSJed Brown       dm   = dmcoarse;
916b17ce1afSJed Brown     }
917b9c2fdf1SPeter Brune   }
918b17ce1afSJed Brown 
919f833ba53SLisandro Dalcin   for (i = 0; i < snes->max_its; i++) {
920fe6f9142SPeter Brune     /* Call general purpose update function */
921fe6f9142SPeter Brune     if (snes->ops->update) {
922fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
923fe6f9142SPeter Brune     }
924f833ba53SLisandro Dalcin 
92507144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
92691f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
9278c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
92891f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
9298c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
9308c94862eSPeter Brune       ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr);
93134d65b3cSPeter Brune     } else if (fas->fastype == SNES_FAS_KASKADE) {
93234d65b3cSPeter Brune       ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr);
9336c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");
934742fe5e2SPeter Brune 
935742fe5e2SPeter Brune     /* check for FAS cycle divergence */
9361aa26658SKarl Rupp     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
937b9c2fdf1SPeter Brune 
938c90fad12SPeter Brune     /* Monitor convergence */
939e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
940c90fad12SPeter Brune     snes->iter = i+1;
941e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
942a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
943c90fad12SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
944c90fad12SPeter Brune     /* Test for convergence */
94566585501SPeter Brune     if (isFine) {
946b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
947c90fad12SPeter Brune       if (snes->reason) break;
948fe6f9142SPeter Brune     }
94966585501SPeter Brune   }
950f833ba53SLisandro Dalcin   if (i == snes->max_its) {
951f833ba53SLisandro Dalcin     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", i);CHKERRQ(ierr);
952fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
953fe6f9142SPeter Brune   }
954421d9b32SPeter Brune   PetscFunctionReturn(0);
955421d9b32SPeter Brune }
95640244768SBarry Smith 
95740244768SBarry Smith /*MC
95840244768SBarry Smith 
95940244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
96040244768SBarry Smith 
96140244768SBarry Smith    The nonlinear problem is solved by correction using coarse versions
96240244768SBarry Smith    of the nonlinear problem.  This problem is perturbed so that a projected
96340244768SBarry Smith    solution of the fine problem elicits no correction from the coarse problem.
96440244768SBarry Smith 
96540244768SBarry Smith Options Database:
96640244768SBarry Smith +   -snes_fas_levels -  The number of levels
96740244768SBarry Smith .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
96840244768SBarry Smith .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  Additive or multiplicative cycle
96940244768SBarry Smith .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
97040244768SBarry Smith .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
97140244768SBarry Smith .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
97240244768SBarry Smith .   -snes_fas_monitor -  Monitor progress of all of the levels
97340244768SBarry Smith .   -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
97440244768SBarry Smith .   -fas_levels_snes_ -  SNES options for all smoothers
97540244768SBarry Smith .   -fas_levels_cycle_snes_ -  SNES options for all cycles
97640244768SBarry Smith .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
97740244768SBarry Smith .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
97840244768SBarry Smith -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
97940244768SBarry Smith 
98040244768SBarry Smith Notes:
98140244768SBarry Smith    The organization of the FAS solver is slightly different from the organization of PCMG
98240244768SBarry Smith    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
98340244768SBarry Smith    The cycle SNES instance may be used for monitoring convergence on a particular level.
98440244768SBarry Smith 
98540244768SBarry Smith Level: beginner
98640244768SBarry Smith 
98740244768SBarry Smith    References:
98840244768SBarry Smith . 1. -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
98940244768SBarry Smith    SIAM Review, 57(4), 2015
99040244768SBarry Smith 
99140244768SBarry Smith .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
99240244768SBarry Smith M*/
99340244768SBarry Smith 
99440244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
99540244768SBarry Smith {
99640244768SBarry Smith   SNES_FAS       *fas;
99740244768SBarry Smith   PetscErrorCode ierr;
99840244768SBarry Smith 
99940244768SBarry Smith   PetscFunctionBegin;
100040244768SBarry Smith   snes->ops->destroy        = SNESDestroy_FAS;
100140244768SBarry Smith   snes->ops->setup          = SNESSetUp_FAS;
100240244768SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
100340244768SBarry Smith   snes->ops->view           = SNESView_FAS;
100440244768SBarry Smith   snes->ops->solve          = SNESSolve_FAS;
100540244768SBarry Smith   snes->ops->reset          = SNESReset_FAS;
100640244768SBarry Smith 
100740244768SBarry Smith   snes->usesksp = PETSC_FALSE;
1008efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
100940244768SBarry Smith 
101040244768SBarry Smith   if (!snes->tolerancesset) {
101140244768SBarry Smith     snes->max_funcs = 30000;
101240244768SBarry Smith     snes->max_its   = 10000;
101340244768SBarry Smith   }
101440244768SBarry Smith 
10154fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
10164fc747eaSLawrence Mitchell 
101740244768SBarry Smith   ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr);
101840244768SBarry Smith 
101940244768SBarry Smith   snes->data                  = (void*) fas;
102040244768SBarry Smith   fas->level                  = 0;
102140244768SBarry Smith   fas->levels                 = 1;
102240244768SBarry Smith   fas->n_cycles               = 1;
102340244768SBarry Smith   fas->max_up_it              = 1;
102440244768SBarry Smith   fas->max_down_it            = 1;
102540244768SBarry Smith   fas->smoothu                = NULL;
102640244768SBarry Smith   fas->smoothd                = NULL;
102740244768SBarry Smith   fas->next                   = NULL;
102840244768SBarry Smith   fas->previous               = NULL;
102940244768SBarry Smith   fas->fine                   = snes;
103040244768SBarry Smith   fas->interpolate            = NULL;
103140244768SBarry Smith   fas->restrct                = NULL;
103240244768SBarry Smith   fas->inject                 = NULL;
103340244768SBarry Smith   fas->usedmfornumberoflevels = PETSC_FALSE;
103440244768SBarry Smith   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
103540244768SBarry Smith   fas->full_downsweep         = PETSC_FALSE;
10366dfbafefSToby Isaac   fas->full_total             = PETSC_FALSE;
103740244768SBarry Smith 
103840244768SBarry Smith   fas->eventsmoothsetup    = 0;
103940244768SBarry Smith   fas->eventsmoothsolve    = 0;
104040244768SBarry Smith   fas->eventresidual       = 0;
104140244768SBarry Smith   fas->eventinterprestrict = 0;
104240244768SBarry Smith   PetscFunctionReturn(0);
104340244768SBarry Smith }
1044