xref: /petsc/src/snes/impls/fas/fas.c (revision c8c899ca14d28eb40f2b95e79d22697c89496ce1)
1421d9b32SPeter Brune /* Defines the basic SNES object */
26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3421d9b32SPeter Brune 
407144faaSPeter Brune const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0};
507144faaSPeter Brune 
6421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
7421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
8421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
9421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
10421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
11421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
126273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
13ab8d36c9SPeter Brune extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *);
14421d9b32SPeter Brune 
151fbfccc6SJed Brown /*MC
161fbfccc6SJed Brown 
171fbfccc6SJed Brown SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
181fbfccc6SJed Brown 
191fbfccc6SJed Brown The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
201fbfccc6SJed Brown 
211fbfccc6SJed Brown Level: advanced
221fbfccc6SJed Brown 
231fbfccc6SJed Brown .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
241fbfccc6SJed Brown M*/
25421d9b32SPeter Brune 
26421d9b32SPeter Brune #undef __FUNCT__
27421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
281fbfccc6SJed Brown PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes)
29421d9b32SPeter Brune {
30421d9b32SPeter Brune   SNES_FAS * fas;
31421d9b32SPeter Brune   PetscErrorCode ierr;
32421d9b32SPeter Brune 
33421d9b32SPeter Brune   PetscFunctionBegin;
34421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
35421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
36421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
37421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
38421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
39421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
40421d9b32SPeter Brune 
41ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
42ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
43ed020824SBarry Smith 
4488976e71SPeter Brune   if (!snes->tolerancesset) {
450e444f03SPeter Brune     snes->max_funcs = 30000;
460e444f03SPeter Brune     snes->max_its   = 10000;
4788976e71SPeter Brune   }
480e444f03SPeter Brune 
49421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
50421d9b32SPeter Brune   snes->data                  = (void*) fas;
51421d9b32SPeter Brune   fas->level                  = 0;
52293a7e31SPeter Brune   fas->levels                 = 1;
53ee78dd50SPeter Brune   fas->n_cycles               = 1;
54ee78dd50SPeter Brune   fas->max_up_it              = 1;
55ee78dd50SPeter Brune   fas->max_down_it            = 1;
56ab8d36c9SPeter Brune   fas->smoothu                = PETSC_NULL;
57ab8d36c9SPeter Brune   fas->smoothd                = PETSC_NULL;
58421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
596273346dSPeter Brune   fas->previous               = PETSC_NULL;
60ab8d36c9SPeter Brune   fas->fine                   = snes;
61421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
62421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
63efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
64cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
65cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
66ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
67efe1f98aSPeter Brune   PetscFunctionReturn(0);
68efe1f98aSPeter Brune }
69efe1f98aSPeter Brune 
70421d9b32SPeter Brune #undef __FUNCT__
71421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
72421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
73421d9b32SPeter Brune {
7477df8cc4SPeter Brune   PetscErrorCode ierr = 0;
75421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
76421d9b32SPeter Brune 
77421d9b32SPeter Brune   PetscFunctionBegin;
78ab8d36c9SPeter Brune   if (fas->smoothu != fas->smoothd) {
79ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
80ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
81ab8d36c9SPeter Brune   } else {
82ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
83ab8d36c9SPeter Brune     fas->smoothu = PETSC_NULL;
84ab8d36c9SPeter Brune   }
853dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
86bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
87bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
88bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
89742fe5e2SPeter Brune   if (fas->next)      ierr = SNESReset(fas->next);CHKERRQ(ierr);
9022c1e704SPeter Brune 
91421d9b32SPeter Brune   PetscFunctionReturn(0);
92421d9b32SPeter Brune }
93421d9b32SPeter Brune 
94421d9b32SPeter Brune #undef __FUNCT__
95421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
96421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
97421d9b32SPeter Brune {
98421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
99742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
100421d9b32SPeter Brune 
101421d9b32SPeter Brune   PetscFunctionBegin;
102421d9b32SPeter Brune   /* recursively resets and then destroys */
10379d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
104421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
105421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
106ee1fd11aSPeter Brune 
107421d9b32SPeter Brune   PetscFunctionReturn(0);
108421d9b32SPeter Brune }
109421d9b32SPeter Brune 
110421d9b32SPeter Brune #undef __FUNCT__
111421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
112421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
113421d9b32SPeter Brune {
11448bfdf8aSPeter Brune   SNES_FAS                *fas = (SNES_FAS *) snes->data;
115421d9b32SPeter Brune   PetscErrorCode          ierr;
116efe1f98aSPeter Brune   VecScatter              injscatter;
117d1adcc6fSPeter Brune   PetscInt                dm_levels;
1183dccd265SPeter Brune   Vec                     vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
119ab8d36c9SPeter Brune   SNES                    next;
120ab8d36c9SPeter Brune   PetscBool               isFine;
121421d9b32SPeter Brune   PetscFunctionBegin;
122eff52c0eSPeter Brune 
123ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
124ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
125d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
126d1adcc6fSPeter Brune     dm_levels++;
127cc05f883SPeter Brune     if (dm_levels > fas->levels) {
1282e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
1293dccd265SPeter Brune       vec_sol = snes->vec_sol;
1303dccd265SPeter Brune       vec_func = snes->vec_func;
1313dccd265SPeter Brune       vec_sol_update = snes->vec_sol_update;
1323dccd265SPeter Brune       vec_rhs = snes->vec_rhs;
1333dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
1343dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
1353dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
1363dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
1373dccd265SPeter Brune 
1383dccd265SPeter Brune       /* reset the number of levels */
139d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
140cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1413dccd265SPeter Brune 
1423dccd265SPeter Brune       snes->vec_sol = vec_sol;
1433dccd265SPeter Brune       snes->vec_func = vec_func;
1443dccd265SPeter Brune       snes->vec_rhs = vec_rhs;
1453dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
146d1adcc6fSPeter Brune     }
147d1adcc6fSPeter Brune   }
148ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
149ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
1503dccd265SPeter Brune 
15107144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
152fde0ff24SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
15307144faaSPeter Brune   } else {
154ddebd997SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
15507144faaSPeter Brune   }
156cc05f883SPeter Brune 
157ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
158ab8d36c9SPeter Brune   if (!fas->smoothd) {
159ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
160ab8d36c9SPeter Brune   }
161ab8d36c9SPeter Brune 
16279d9a41aSPeter Brune   if (snes->dm) {
163ab8d36c9SPeter Brune     /* set the smoother DMs properly */
164ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
165ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
16679d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
167ab8d36c9SPeter Brune     if (next) {
16879d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
169ab8d36c9SPeter Brune       if (!next->dm) {
170ab8d36c9SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
171ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
17279d9a41aSPeter Brune       }
17379d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
17479d9a41aSPeter Brune       if (!fas->interpolate) {
175ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
176bccf9bb3SJed Brown         if (!fas->restrct) {
177bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
17879d9a41aSPeter Brune           fas->restrct = fas->interpolate;
17979d9a41aSPeter Brune         }
180bccf9bb3SJed Brown       }
18179d9a41aSPeter Brune       /* set the injection from the DM */
18279d9a41aSPeter Brune       if (!fas->inject) {
183ab8d36c9SPeter Brune         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
18479d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
18579d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
18679d9a41aSPeter Brune       }
18779d9a41aSPeter Brune     }
18879d9a41aSPeter Brune   }
18979d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
19079d9a41aSPeter Brune 
19179d9a41aSPeter Brune   if (fas->galerkin) {
192ab8d36c9SPeter Brune     if (next)
193ab8d36c9SPeter Brune       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
194ab8d36c9SPeter Brune     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
195ab8d36c9SPeter Brune     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
19679d9a41aSPeter Brune   }
19779d9a41aSPeter Brune 
198ab8d36c9SPeter Brune   if (next) {
19979d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
200ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
201ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
202ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
20379d9a41aSPeter Brune   }
2046273346dSPeter Brune   /* setup FAS work vectors */
2056273346dSPeter Brune   if (fas->galerkin) {
2066273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
2076273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
2086273346dSPeter Brune   }
209421d9b32SPeter Brune   PetscFunctionReturn(0);
210421d9b32SPeter Brune }
211421d9b32SPeter Brune 
212421d9b32SPeter Brune #undef __FUNCT__
213421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
214421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
215421d9b32SPeter Brune {
216ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
217ee78dd50SPeter Brune   PetscInt       levels = 1;
218ab8d36c9SPeter Brune   PetscBool      flg, monflg, galerkinflg;
219421d9b32SPeter Brune   PetscErrorCode ierr;
220ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
22107144faaSPeter Brune   SNESFASType    fastype;
222fde0ff24SPeter Brune   const char     *optionsprefix;
223f1c6b773SPeter Brune   SNESLineSearch linesearch;
224ab8d36c9SPeter Brune   PetscInt       m;
225ab8d36c9SPeter Brune   SNES           next;
226ab8d36c9SPeter Brune   PetscBool      isFine;
227421d9b32SPeter Brune 
228421d9b32SPeter Brune   PetscFunctionBegin;
229ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
230c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
231ee78dd50SPeter Brune 
232ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
233ab8d36c9SPeter Brune   if (isFine) {
234ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
235c732cbdbSBarry Smith     if (!flg && snes->dm) {
236c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
237c732cbdbSBarry Smith       levels++;
238d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
239c732cbdbSBarry Smith     }
240ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
24107144faaSPeter Brune     fastype = fas->fastype;
24207144faaSPeter Brune     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
24307144faaSPeter Brune     if (flg) {
24407144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
24507144faaSPeter Brune     }
246ee78dd50SPeter Brune 
247fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
248ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
249ab8d36c9SPeter Brune     if (flg) {
250ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
251fde0ff24SPeter Brune     }
252fde0ff24SPeter Brune 
253ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
254ab8d36c9SPeter Brune     if (flg) {
255ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
256ab8d36c9SPeter Brune     }
257ee78dd50SPeter Brune 
258ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
259ab8d36c9SPeter Brune     if (flg) {
260ab8d36c9SPeter Brune       ierr = SNESFASSetNumberSmoothUp(snes,m);CHKERRQ(ierr);
261ab8d36c9SPeter Brune     }
262ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
263ab8d36c9SPeter Brune     if (flg) {
264ab8d36c9SPeter Brune       ierr = SNESFASSetNumberSmoothDown(snes,m);CHKERRQ(ierr);
265ab8d36c9SPeter Brune     }
266*c8c899caSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
267*c8c899caSPeter Brune     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
268ab8d36c9SPeter Brune   }
269ee78dd50SPeter Brune 
270421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2718cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
272eff52c0eSPeter Brune 
2739e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
2749e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
2759e764e56SPeter Brune     if (!snes->linesearch) {
276f1c6b773SPeter Brune       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
277f1c6b773SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
2789e764e56SPeter Brune     }
2799e764e56SPeter Brune   }
2809e764e56SPeter Brune 
281ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
282ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
283ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
284421d9b32SPeter Brune   PetscFunctionReturn(0);
285421d9b32SPeter Brune }
286421d9b32SPeter Brune 
287421d9b32SPeter Brune #undef __FUNCT__
288421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
289421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
290421d9b32SPeter Brune {
291421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
292ab8d36c9SPeter Brune   PetscBool      isFine, iascii;
293ab8d36c9SPeter Brune   PetscInt       i;
294421d9b32SPeter Brune   PetscErrorCode ierr;
295ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
296421d9b32SPeter Brune 
297421d9b32SPeter Brune   PetscFunctionBegin;
298ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
299ab8d36c9SPeter Brune   if (isFine) {
300421d9b32SPeter Brune     ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
301421d9b32SPeter Brune     if (iascii) {
302ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
303ab8d36c9SPeter Brune       if (fas->galerkin) {
304ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
305421d9b32SPeter Brune       } else {
306ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
307421d9b32SPeter Brune       }
308ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
309ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
310ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
311ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
312ab8d36c9SPeter Brune         if (!i) {
313ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
314421d9b32SPeter Brune         } else {
315ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
316421d9b32SPeter Brune         }
317ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
318ab8d36c9SPeter Brune         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
319ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
320ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
321ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
322ab8d36c9SPeter Brune         } else if (i) {
323ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
324ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
325ab8d36c9SPeter Brune           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
326ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
327ab8d36c9SPeter Brune         }
328ab8d36c9SPeter Brune       }
329421d9b32SPeter Brune     } else {
330421d9b32SPeter Brune       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
331421d9b32SPeter Brune     }
332ab8d36c9SPeter Brune   }
333421d9b32SPeter Brune   PetscFunctionReturn(0);
334421d9b32SPeter Brune }
335421d9b32SPeter Brune 
336421d9b32SPeter Brune #undef __FUNCT__
33739bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth"
33839bd7f45SPeter Brune /*
33939bd7f45SPeter Brune Defines the action of the downsmoother
34039bd7f45SPeter Brune  */
34139bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
34239bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
343742fe5e2SPeter Brune   SNESConvergedReason reason;
344ab8d36c9SPeter Brune   Vec                 FPC;
345ab8d36c9SPeter Brune   SNES                smoothd;
346421d9b32SPeter Brune   PetscFunctionBegin;
347ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
348ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
349742fe5e2SPeter Brune   /* check convergence reason for the smoother */
350ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
351742fe5e2SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
352742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
353742fe5e2SPeter Brune     PetscFunctionReturn(0);
354742fe5e2SPeter Brune   }
355ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
3564b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
35739bd7f45SPeter Brune   PetscFunctionReturn(0);
35839bd7f45SPeter Brune }
35939bd7f45SPeter Brune 
36039bd7f45SPeter Brune 
36139bd7f45SPeter Brune #undef __FUNCT__
36239bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth"
36339bd7f45SPeter Brune /*
36407144faaSPeter Brune Defines the action of the upsmoother
36539bd7f45SPeter Brune  */
36639bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
36739bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
36839bd7f45SPeter Brune   SNESConvergedReason reason;
369ab8d36c9SPeter Brune   Vec                 FPC;
370ab8d36c9SPeter Brune   SNES                smoothu;
37139bd7f45SPeter Brune   PetscFunctionBegin;
372ab8d36c9SPeter Brune 
373ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
374ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
37539bd7f45SPeter Brune   /* check convergence reason for the smoother */
376ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
37739bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
37839bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
37939bd7f45SPeter Brune     PetscFunctionReturn(0);
38039bd7f45SPeter Brune   }
381ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
3824b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
38339bd7f45SPeter Brune   PetscFunctionReturn(0);
38439bd7f45SPeter Brune }
38539bd7f45SPeter Brune 
38639bd7f45SPeter Brune #undef __FUNCT__
387938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
388938e4a01SJed Brown /*@
389938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
390938e4a01SJed Brown 
391938e4a01SJed Brown    Collective
392938e4a01SJed Brown 
393938e4a01SJed Brown    Input Arguments:
394938e4a01SJed Brown .  snes - SNESFAS
395938e4a01SJed Brown 
396938e4a01SJed Brown    Output Arguments:
397938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
398938e4a01SJed Brown 
399938e4a01SJed Brown    Level: developer
400938e4a01SJed Brown 
401938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
402938e4a01SJed Brown @*/
403938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
404938e4a01SJed Brown {
405938e4a01SJed Brown   PetscErrorCode ierr;
406938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
407938e4a01SJed Brown 
408938e4a01SJed Brown   PetscFunctionBegin;
409938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
410938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
411938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
412938e4a01SJed Brown   PetscFunctionReturn(0);
413938e4a01SJed Brown }
414938e4a01SJed Brown 
415e9923e8dSJed Brown #undef __FUNCT__
416e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
417e9923e8dSJed Brown /*@
418e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
419e9923e8dSJed Brown 
420e9923e8dSJed Brown    Collective
421e9923e8dSJed Brown 
422e9923e8dSJed Brown    Input Arguments:
423e9923e8dSJed Brown +  fine - SNES from which to restrict
424e9923e8dSJed Brown -  Xfine - vector to restrict
425e9923e8dSJed Brown 
426e9923e8dSJed Brown    Output Arguments:
427e9923e8dSJed Brown .  Xcoarse - result of restriction
428e9923e8dSJed Brown 
429e9923e8dSJed Brown    Level: developer
430e9923e8dSJed Brown 
431e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
432e9923e8dSJed Brown @*/
433e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
434e9923e8dSJed Brown {
435e9923e8dSJed Brown   PetscErrorCode ierr;
436e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
437e9923e8dSJed Brown 
438e9923e8dSJed Brown   PetscFunctionBegin;
439e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
440e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
441e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
442e9923e8dSJed Brown   if (fas->inject) {
443e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
444e9923e8dSJed Brown   } else {
445e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
446e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
447e9923e8dSJed Brown   }
448e9923e8dSJed Brown   PetscFunctionReturn(0);
449e9923e8dSJed Brown }
450e9923e8dSJed Brown 
451e9923e8dSJed Brown #undef __FUNCT__
45239bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection"
45339bd7f45SPeter Brune /*
45439bd7f45SPeter Brune 
45539bd7f45SPeter Brune Performs the FAS coarse correction as:
45639bd7f45SPeter Brune 
45739bd7f45SPeter Brune fine problem: F(x) = 0
45839bd7f45SPeter Brune coarse problem: F^c(x) = b^c
45939bd7f45SPeter Brune 
46039bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
46139bd7f45SPeter Brune 
46239bd7f45SPeter Brune  */
46339a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
46439bd7f45SPeter Brune   PetscErrorCode      ierr;
46539bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
46639bd7f45SPeter Brune   SNESConvergedReason reason;
467ab8d36c9SPeter Brune   SNES                next;
468ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
46939bd7f45SPeter Brune   PetscFunctionBegin;
470ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
471ab8d36c9SPeter Brune   if (next) {
472ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
473ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
474ab8d36c9SPeter Brune 
475ab8d36c9SPeter Brune     X_c  = next->vec_sol;
476ab8d36c9SPeter Brune     Xo_c = next->work[0];
477ab8d36c9SPeter Brune     F_c  = next->vec_func;
478ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
479efe1f98aSPeter Brune 
480938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
481293a7e31SPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
482293a7e31SPeter Brune 
483293a7e31SPeter Brune     /* restrict the defect */
484ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
485293a7e31SPeter Brune 
486c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
487ab8d36c9SPeter Brune     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
488ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
489742fe5e2SPeter Brune 
490293a7e31SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
491c90fad12SPeter Brune 
492ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
493ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
494c90fad12SPeter Brune 
495c90fad12SPeter Brune     /* recurse to the next level */
496ab8d36c9SPeter Brune     next->vec_rhs = B_c;
497ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
498ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
499742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
500742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
501742fe5e2SPeter Brune       PetscFunctionReturn(0);
502742fe5e2SPeter Brune     }
503ee78dd50SPeter Brune 
504fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
505fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
506ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
507293a7e31SPeter Brune   }
50839bd7f45SPeter Brune   PetscFunctionReturn(0);
50939bd7f45SPeter Brune }
51039bd7f45SPeter Brune 
51139bd7f45SPeter Brune #undef __FUNCT__
51239bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive"
51339bd7f45SPeter Brune /*
51439bd7f45SPeter Brune 
51539bd7f45SPeter Brune The additive cycle looks like:
51639bd7f45SPeter Brune 
51707144faaSPeter Brune xhat = x
51807144faaSPeter Brune xhat = dS(x, b)
51907144faaSPeter Brune x = coarsecorrection(xhat, b_d)
52007144faaSPeter Brune x = x + nu*(xhat - x);
52139bd7f45SPeter Brune (optional) x = uS(x, b)
52239bd7f45SPeter Brune 
52339bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
52439bd7f45SPeter Brune 
52539bd7f45SPeter Brune  */
52639bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
52707144faaSPeter Brune   Vec                 F, B, Xhat;
52822c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
52939bd7f45SPeter Brune   PetscErrorCode      ierr;
53007144faaSPeter Brune   SNESConvergedReason reason;
53122c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
53222c1e704SPeter Brune   PetscBool           lssuccess;
533ab8d36c9SPeter Brune   SNES                next;
534ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
53539bd7f45SPeter Brune   PetscFunctionBegin;
536ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
53739bd7f45SPeter Brune   F = snes->vec_func;
53839bd7f45SPeter Brune   B = snes->vec_rhs;
539fde0ff24SPeter Brune   Xhat = snes->work[3];
54007144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
54107144faaSPeter Brune   /* recurse first */
542ab8d36c9SPeter Brune   if (next) {
543ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
544ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
54507144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
54639bd7f45SPeter Brune 
547ab8d36c9SPeter Brune     X_c  = next->vec_sol;
548ab8d36c9SPeter Brune     Xo_c = next->work[0];
549ab8d36c9SPeter Brune     F_c  = next->vec_func;
550ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
55139bd7f45SPeter Brune 
552938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
55307144faaSPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
55407144faaSPeter Brune 
55507144faaSPeter Brune     /* restrict the defect */
556ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
55707144faaSPeter Brune 
55807144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
559ab8d36c9SPeter Brune     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
560ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
56107144faaSPeter Brune 
56207144faaSPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
56307144faaSPeter Brune 
56407144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
56507144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
56607144faaSPeter Brune 
56707144faaSPeter Brune     /* recurse */
568ab8d36c9SPeter Brune     next->vec_rhs = B_c;
569ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
57007144faaSPeter Brune 
57107144faaSPeter Brune     /* smooth on this level */
57207144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
57307144faaSPeter Brune 
574ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
57507144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
57607144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
57707144faaSPeter Brune       PetscFunctionReturn(0);
57807144faaSPeter Brune     }
57907144faaSPeter Brune 
58007144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
581c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
582ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
58307144faaSPeter Brune 
584ddebd997SPeter Brune     /* additive correction of the coarse direction*/
585ddebd997SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
586ddebd997SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
587f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
588f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
5899e764e56SPeter Brune     if (!lssuccess) {
5909e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
5919e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
5929e764e56SPeter Brune         PetscFunctionReturn(0);
5939e764e56SPeter Brune       }
5949e764e56SPeter Brune     }
595f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
59607144faaSPeter Brune   } else {
59707144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
59807144faaSPeter Brune   }
59939bd7f45SPeter Brune   PetscFunctionReturn(0);
60039bd7f45SPeter Brune }
60139bd7f45SPeter Brune 
60239bd7f45SPeter Brune #undef __FUNCT__
60339bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative"
60439bd7f45SPeter Brune /*
60539bd7f45SPeter Brune 
60639bd7f45SPeter Brune Defines the FAS cycle as:
60739bd7f45SPeter Brune 
60839bd7f45SPeter Brune fine problem: F(x) = 0
60939bd7f45SPeter Brune coarse problem: F^c(x) = b^c
61039bd7f45SPeter Brune 
61139bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
61239bd7f45SPeter Brune 
61339bd7f45SPeter Brune correction:
61439bd7f45SPeter Brune 
61539bd7f45SPeter Brune x = x + I(x^c - Rx)
61639bd7f45SPeter Brune 
61739bd7f45SPeter Brune  */
61839bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
61939bd7f45SPeter Brune 
62039bd7f45SPeter Brune   PetscErrorCode      ierr;
62139bd7f45SPeter Brune   Vec                 F,B;
62239bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
62339bd7f45SPeter Brune 
62439bd7f45SPeter Brune   PetscFunctionBegin;
62539bd7f45SPeter Brune   F = snes->vec_func;
62639bd7f45SPeter Brune   B = snes->vec_rhs;
62739bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
62839bd7f45SPeter Brune   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
62939bd7f45SPeter Brune 
63039a4b097SPeter Brune   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
63139bd7f45SPeter Brune 
632c90fad12SPeter Brune   if (fas->level != 0) {
63339bd7f45SPeter Brune     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
634fe6f9142SPeter Brune   }
635fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
636fa9694d7SPeter Brune 
637fa9694d7SPeter Brune   PetscFunctionReturn(0);
638421d9b32SPeter Brune }
639421d9b32SPeter Brune 
640421d9b32SPeter Brune #undef __FUNCT__
641421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
642421d9b32SPeter Brune 
643421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
644421d9b32SPeter Brune {
645fa9694d7SPeter Brune   PetscErrorCode ierr;
646fe6f9142SPeter Brune   PetscInt       i, maxits;
647ddb5aff1SPeter Brune   Vec            X, F;
648fe6f9142SPeter Brune   PetscReal      fnorm;
649b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
650b17ce1afSJed Brown   DM             dm;
651b17ce1afSJed Brown 
652421d9b32SPeter Brune   PetscFunctionBegin;
653fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
654fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
655fa9694d7SPeter Brune   X = snes->vec_sol;
656f5a6d4f9SBarry Smith   F = snes->vec_func;
657293a7e31SPeter Brune 
658293a7e31SPeter Brune   /*norm setup */
659fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
660fe6f9142SPeter Brune   snes->iter = 0;
661fe6f9142SPeter Brune   snes->norm = 0.;
662fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
663fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
664fe6f9142SPeter Brune   if (snes->domainerror) {
665fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
666fe6f9142SPeter Brune     PetscFunctionReturn(0);
667fe6f9142SPeter Brune   }
668fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
669fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
670fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
671fe6f9142SPeter Brune   snes->norm = fnorm;
672fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
673fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
674fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
675fe6f9142SPeter Brune 
676fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
677fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
678fe6f9142SPeter Brune   /* test convergence */
679fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
680fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
681b17ce1afSJed Brown 
682b17ce1afSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
683b17ce1afSJed Brown   for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
684b17ce1afSJed Brown     DM dmcoarse;
685b17ce1afSJed Brown     ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
686b17ce1afSJed Brown     ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
687b17ce1afSJed Brown     dm = dmcoarse;
688b17ce1afSJed Brown   }
689b17ce1afSJed Brown 
690fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
691fe6f9142SPeter Brune     /* Call general purpose update function */
692646217ecSPeter Brune 
693fe6f9142SPeter Brune     if (snes->ops->update) {
694fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
695fe6f9142SPeter Brune     }
69607144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
69739bd7f45SPeter Brune       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
69807144faaSPeter Brune     } else {
69907144faaSPeter Brune       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
70007144faaSPeter Brune     }
701742fe5e2SPeter Brune 
702742fe5e2SPeter Brune     /* check for FAS cycle divergence */
703742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
704742fe5e2SPeter Brune       PetscFunctionReturn(0);
705742fe5e2SPeter Brune     }
706c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
707c90fad12SPeter Brune     /* Monitor convergence */
708c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
709c90fad12SPeter Brune     snes->iter = i+1;
710c90fad12SPeter Brune     snes->norm = fnorm;
711c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
712c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
713c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
714c90fad12SPeter Brune     /* Test for convergence */
715c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
716c90fad12SPeter Brune     if (snes->reason) break;
717fe6f9142SPeter Brune   }
718fe6f9142SPeter Brune   if (i == maxits) {
719fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
720fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
721fe6f9142SPeter Brune   }
722421d9b32SPeter Brune   PetscFunctionReturn(0);
723421d9b32SPeter Brune }
724