xref: /petsc/src/snes/impls/fas/fas.c (revision 534ebe21143b28d6375c549c032c2822f1bc913f)
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 
19d3bc2379SPeter Brune    The nonlinear problem is solved by correction using coarse versions
20d3bc2379SPeter Brune    of the nonlinear problem.  This problem is perturbed so that a projected
21d3bc2379SPeter Brune    solution of the fine problem elicits no correction from the coarse problem.
22d3bc2379SPeter Brune 
23d3bc2379SPeter Brune Options Database:
24d3bc2379SPeter Brune +   -snes_fas_levels -  The number of levels
25d3bc2379SPeter Brune .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
26d3bc2379SPeter Brune .   -snes_fas_type<additive, multiplicative>  -  Additive or multiplicative cycle
27d3bc2379SPeter Brune .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
28d3bc2379SPeter Brune .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
29d3bc2379SPeter Brune .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
30d3bc2379SPeter Brune .   -snes_fas_monitor -  Monitor progress of all of the levels
31d3bc2379SPeter Brune .   -fas_levels_snes_ -  SNES options for all smoothers
32d3bc2379SPeter Brune .   -fas_levels_cycle_snes -  SNES options for all cycles
33d3bc2379SPeter Brune .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
34d3bc2379SPeter Brune .   -fas_levels_i_cycle_snes - SNES options for the cycle on level i
35d3bc2379SPeter Brune -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
36d3bc2379SPeter Brune 
37d3bc2379SPeter Brune Notes:
38d3bc2379SPeter Brune    The organization of the FAS solver is slightly different from the organization of PCMG
39d3bc2379SPeter Brune    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
40d3bc2379SPeter Brune    The cycle SNES instance may be used for monitoring convergence on a particular level.
411fbfccc6SJed Brown 
421fbfccc6SJed Brown Level: advanced
431fbfccc6SJed Brown 
44d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
451fbfccc6SJed Brown M*/
46421d9b32SPeter Brune 
47421d9b32SPeter Brune #undef __FUNCT__
48421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
491fbfccc6SJed Brown PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes)
50421d9b32SPeter Brune {
51421d9b32SPeter Brune   SNES_FAS * fas;
52421d9b32SPeter Brune   PetscErrorCode ierr;
53421d9b32SPeter Brune 
54421d9b32SPeter Brune   PetscFunctionBegin;
55421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
56421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
57421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
58421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
59421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
60421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
61421d9b32SPeter Brune 
62ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
63ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
64ed020824SBarry Smith 
6588976e71SPeter Brune   if (!snes->tolerancesset) {
660e444f03SPeter Brune     snes->max_funcs = 30000;
670e444f03SPeter Brune     snes->max_its   = 10000;
6888976e71SPeter Brune   }
690e444f03SPeter Brune 
70421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
71421d9b32SPeter Brune   snes->data                  = (void*) fas;
72421d9b32SPeter Brune   fas->level                  = 0;
73293a7e31SPeter Brune   fas->levels                 = 1;
74ee78dd50SPeter Brune   fas->n_cycles               = 1;
75ee78dd50SPeter Brune   fas->max_up_it              = 1;
76ee78dd50SPeter Brune   fas->max_down_it            = 1;
77ab8d36c9SPeter Brune   fas->smoothu                = PETSC_NULL;
78ab8d36c9SPeter Brune   fas->smoothd                = PETSC_NULL;
79421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
806273346dSPeter Brune   fas->previous               = PETSC_NULL;
81ab8d36c9SPeter Brune   fas->fine                   = snes;
82421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
83421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
84efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
85cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
86cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
87ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
88efe1f98aSPeter Brune   PetscFunctionReturn(0);
89efe1f98aSPeter Brune }
90efe1f98aSPeter Brune 
91421d9b32SPeter Brune #undef __FUNCT__
92421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
93421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
94421d9b32SPeter Brune {
9577df8cc4SPeter Brune   PetscErrorCode ierr = 0;
96421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
97421d9b32SPeter Brune 
98421d9b32SPeter Brune   PetscFunctionBegin;
99ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
100ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
1013dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
102bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
103bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
104bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
105742fe5e2SPeter Brune   if (fas->next)      ierr = SNESReset(fas->next);CHKERRQ(ierr);
10622c1e704SPeter Brune 
107421d9b32SPeter Brune   PetscFunctionReturn(0);
108421d9b32SPeter Brune }
109421d9b32SPeter Brune 
110421d9b32SPeter Brune #undef __FUNCT__
111421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
112421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
113421d9b32SPeter Brune {
114421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
115742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
116421d9b32SPeter Brune 
117421d9b32SPeter Brune   PetscFunctionBegin;
118421d9b32SPeter Brune   /* recursively resets and then destroys */
11979d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
120421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
121421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
122ee1fd11aSPeter Brune 
123421d9b32SPeter Brune   PetscFunctionReturn(0);
124421d9b32SPeter Brune }
125421d9b32SPeter Brune 
126421d9b32SPeter Brune #undef __FUNCT__
127421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
128421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
129421d9b32SPeter Brune {
13048bfdf8aSPeter Brune   SNES_FAS                *fas = (SNES_FAS *) snes->data;
131421d9b32SPeter Brune   PetscErrorCode          ierr;
132efe1f98aSPeter Brune   VecScatter              injscatter;
133d1adcc6fSPeter Brune   PetscInt                dm_levels;
1343dccd265SPeter Brune   Vec                     vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
135ab8d36c9SPeter Brune   SNES                    next;
136ab8d36c9SPeter Brune   PetscBool               isFine;
137421d9b32SPeter Brune   PetscFunctionBegin;
138eff52c0eSPeter Brune 
139ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
140ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
141d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
142d1adcc6fSPeter Brune     dm_levels++;
143cc05f883SPeter Brune     if (dm_levels > fas->levels) {
1442e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
1453dccd265SPeter Brune       vec_sol = snes->vec_sol;
1463dccd265SPeter Brune       vec_func = snes->vec_func;
1473dccd265SPeter Brune       vec_sol_update = snes->vec_sol_update;
1483dccd265SPeter Brune       vec_rhs = snes->vec_rhs;
1493dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
1503dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
1513dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
1523dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
1533dccd265SPeter Brune 
1543dccd265SPeter Brune       /* reset the number of levels */
155d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
156cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1573dccd265SPeter Brune 
1583dccd265SPeter Brune       snes->vec_sol = vec_sol;
1593dccd265SPeter Brune       snes->vec_func = vec_func;
1603dccd265SPeter Brune       snes->vec_rhs = vec_rhs;
1613dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
162d1adcc6fSPeter Brune     }
163d1adcc6fSPeter Brune   }
164ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
165ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
1663dccd265SPeter Brune 
16707144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
168e4ed7901SPeter Brune     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
16907144faaSPeter Brune   } else {
170e7f468e7SPeter Brune     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
17107144faaSPeter Brune   }
172cc05f883SPeter Brune 
173ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
174ab8d36c9SPeter Brune   if (!fas->smoothd) {
175ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
176ab8d36c9SPeter Brune   }
177*534ebe21SPeter Brune   if (!fas->smoothu && fas->level != fas->levels - 1) {
178*534ebe21SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
179*534ebe21SPeter Brune   }
180ab8d36c9SPeter Brune 
18179d9a41aSPeter Brune   if (snes->dm) {
182ab8d36c9SPeter Brune     /* set the smoother DMs properly */
183ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
184ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
18579d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
186ab8d36c9SPeter Brune     if (next) {
18779d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
188ab8d36c9SPeter Brune       if (!next->dm) {
189ab8d36c9SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
190ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
19179d9a41aSPeter Brune       }
19279d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
19379d9a41aSPeter Brune       if (!fas->interpolate) {
194ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
195bccf9bb3SJed Brown         if (!fas->restrct) {
196bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
19779d9a41aSPeter Brune           fas->restrct = fas->interpolate;
19879d9a41aSPeter Brune         }
199bccf9bb3SJed Brown       }
20079d9a41aSPeter Brune       /* set the injection from the DM */
20179d9a41aSPeter Brune       if (!fas->inject) {
202ab8d36c9SPeter Brune         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
20379d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
20479d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
20579d9a41aSPeter Brune       }
20679d9a41aSPeter Brune     }
20779d9a41aSPeter Brune   }
20879d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
20979d9a41aSPeter Brune   if (fas->galerkin) {
210ab8d36c9SPeter Brune     if (next)
211ab8d36c9SPeter Brune       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
212ab8d36c9SPeter Brune     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
213ab8d36c9SPeter Brune     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
21479d9a41aSPeter Brune   }
21579d9a41aSPeter Brune 
216*534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
217*534ebe21SPeter Brune   if(fas->smoothd){
218*534ebe21SPeter Brune     if (fas->level == 0) {
219*534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
220*534ebe21SPeter Brune      } else {
221*534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
222*534ebe21SPeter Brune     }
223*534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
224*534ebe21SPeter Brune   }
225*534ebe21SPeter Brune 
226*534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
227*534ebe21SPeter Brune   if(fas->smoothu){
228*534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
229*534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
230*534ebe21SPeter Brune     } else {
231*534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
232*534ebe21SPeter Brune     }
233*534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
234*534ebe21SPeter Brune   }
235d06165b7SPeter Brune 
236ab8d36c9SPeter Brune   if (next) {
23779d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
238ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
239ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
240ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
24179d9a41aSPeter Brune   }
2426273346dSPeter Brune   /* setup FAS work vectors */
2436273346dSPeter Brune   if (fas->galerkin) {
2446273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
2456273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
2466273346dSPeter Brune   }
247421d9b32SPeter Brune   PetscFunctionReturn(0);
248421d9b32SPeter Brune }
249421d9b32SPeter Brune 
250421d9b32SPeter Brune #undef __FUNCT__
251421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
252421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
253421d9b32SPeter Brune {
254ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
255ee78dd50SPeter Brune   PetscInt       levels = 1;
2564d26bfa5SPeter Brune   PetscBool      flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE;
257421d9b32SPeter Brune   PetscErrorCode ierr;
258ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
25907144faaSPeter Brune   SNESFASType    fastype;
260fde0ff24SPeter Brune   const char     *optionsprefix;
261f1c6b773SPeter Brune   SNESLineSearch linesearch;
26266585501SPeter Brune   PetscInt       m, n_up, n_down;
263ab8d36c9SPeter Brune   SNES           next;
264ab8d36c9SPeter Brune   PetscBool      isFine;
265421d9b32SPeter Brune 
266421d9b32SPeter Brune   PetscFunctionBegin;
267ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
268c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
269ee78dd50SPeter Brune 
270ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
271ab8d36c9SPeter Brune   if (isFine) {
272ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
273c732cbdbSBarry Smith     if (!flg && snes->dm) {
274c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
275c732cbdbSBarry Smith       levels++;
276d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
277c732cbdbSBarry Smith     }
278ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
27907144faaSPeter Brune     fastype = fas->fastype;
28007144faaSPeter Brune     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
28107144faaSPeter Brune     if (flg) {
28207144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
28307144faaSPeter Brune     }
284ee78dd50SPeter Brune 
285fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
286ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
287ab8d36c9SPeter Brune     if (flg) {
288ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
289fde0ff24SPeter Brune     }
290fde0ff24SPeter Brune 
291ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
292ab8d36c9SPeter Brune     if (flg) {
293ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
294ab8d36c9SPeter Brune     }
295ee78dd50SPeter Brune 
29666585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
297162d76ddSPeter Brune 
29866585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
299162d76ddSPeter Brune 
300c8c899caSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
301c8c899caSPeter Brune     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
302ab8d36c9SPeter Brune   }
303ee78dd50SPeter Brune 
304421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
3058cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
306162d76ddSPeter Brune   if (upflg) {
30766585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
308162d76ddSPeter Brune   }
309162d76ddSPeter Brune   if (downflg) {
31066585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
311162d76ddSPeter Brune   }
312eff52c0eSPeter Brune 
3139e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3149e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3159e764e56SPeter Brune     if (!snes->linesearch) {
316f1c6b773SPeter Brune       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
3171a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
3189e764e56SPeter Brune     }
3199e764e56SPeter Brune   }
3209e764e56SPeter Brune 
321ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
322ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
323ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
324421d9b32SPeter Brune   PetscFunctionReturn(0);
325421d9b32SPeter Brune }
326421d9b32SPeter Brune 
327421d9b32SPeter Brune #undef __FUNCT__
328421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
329421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
330421d9b32SPeter Brune {
331421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
332ab8d36c9SPeter Brune   PetscBool      isFine, iascii;
333ab8d36c9SPeter Brune   PetscInt       i;
334421d9b32SPeter Brune   PetscErrorCode ierr;
335ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
336421d9b32SPeter Brune 
337421d9b32SPeter Brune   PetscFunctionBegin;
338ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
339ab8d36c9SPeter Brune   if (isFine) {
340421d9b32SPeter Brune     ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
341421d9b32SPeter Brune     if (iascii) {
342ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
343ab8d36c9SPeter Brune       if (fas->galerkin) {
344ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
345421d9b32SPeter Brune       } else {
346ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
347421d9b32SPeter Brune       }
348ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
349ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
350ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
351ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
352ab8d36c9SPeter Brune         if (!i) {
353ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
354421d9b32SPeter Brune         } else {
355ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
356421d9b32SPeter Brune         }
357ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
358ab8d36c9SPeter Brune         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
359ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
360ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
361ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
362ab8d36c9SPeter Brune         } else if (i) {
363ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
364ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
365ab8d36c9SPeter Brune           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
366ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
367ab8d36c9SPeter Brune         }
368ab8d36c9SPeter Brune       }
369421d9b32SPeter Brune     } else {
370421d9b32SPeter Brune       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
371421d9b32SPeter Brune     }
372ab8d36c9SPeter Brune   }
373421d9b32SPeter Brune   PetscFunctionReturn(0);
374421d9b32SPeter Brune }
375421d9b32SPeter Brune 
376421d9b32SPeter Brune #undef __FUNCT__
37791f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private"
37839bd7f45SPeter Brune /*
37939bd7f45SPeter Brune Defines the action of the downsmoother
38039bd7f45SPeter Brune  */
38191f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
382b9c2fdf1SPeter Brune {
38339bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
384742fe5e2SPeter Brune   SNESConvergedReason reason;
385ab8d36c9SPeter Brune   Vec                 FPC;
386ab8d36c9SPeter Brune   SNES                smoothd;
387421d9b32SPeter Brune   PetscFunctionBegin;
388ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
389e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
390e4ed7901SPeter Brune   ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr);
391ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
392742fe5e2SPeter Brune   /* check convergence reason for the smoother */
393ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
394e70c42e5SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
395742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
396742fe5e2SPeter Brune     PetscFunctionReturn(0);
397742fe5e2SPeter Brune   }
398ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
3994b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
400b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
40139bd7f45SPeter Brune   PetscFunctionReturn(0);
40239bd7f45SPeter Brune }
40339bd7f45SPeter Brune 
40439bd7f45SPeter Brune 
40539bd7f45SPeter Brune #undef __FUNCT__
40691f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private"
40739bd7f45SPeter Brune /*
40807144faaSPeter Brune Defines the action of the upsmoother
40939bd7f45SPeter Brune  */
41091f99d7cSPeter Brune PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) {
41139bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
41239bd7f45SPeter Brune   SNESConvergedReason reason;
413ab8d36c9SPeter Brune   Vec                 FPC;
414ab8d36c9SPeter Brune   SNES                smoothu;
41539bd7f45SPeter Brune   PetscFunctionBegin;
416ab8d36c9SPeter Brune 
417ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
418ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
41939bd7f45SPeter Brune   /* check convergence reason for the smoother */
420ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
42139bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
42239bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
42339bd7f45SPeter Brune     PetscFunctionReturn(0);
42439bd7f45SPeter Brune   }
425ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
4264b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
427b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
42839bd7f45SPeter Brune   PetscFunctionReturn(0);
42939bd7f45SPeter Brune }
43039bd7f45SPeter Brune 
43139bd7f45SPeter Brune #undef __FUNCT__
432938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
433938e4a01SJed Brown /*@
434938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
435938e4a01SJed Brown 
436938e4a01SJed Brown    Collective
437938e4a01SJed Brown 
438938e4a01SJed Brown    Input Arguments:
439938e4a01SJed Brown .  snes - SNESFAS
440938e4a01SJed Brown 
441938e4a01SJed Brown    Output Arguments:
442938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
443938e4a01SJed Brown 
444938e4a01SJed Brown    Level: developer
445938e4a01SJed Brown 
446938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
447938e4a01SJed Brown @*/
448938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
449938e4a01SJed Brown {
450938e4a01SJed Brown   PetscErrorCode ierr;
451938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
452938e4a01SJed Brown 
453938e4a01SJed Brown   PetscFunctionBegin;
454938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
455938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
456938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
457938e4a01SJed Brown   PetscFunctionReturn(0);
458938e4a01SJed Brown }
459938e4a01SJed Brown 
460e9923e8dSJed Brown #undef __FUNCT__
461e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
462e9923e8dSJed Brown /*@
463e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
464e9923e8dSJed Brown 
465e9923e8dSJed Brown    Collective
466e9923e8dSJed Brown 
467e9923e8dSJed Brown    Input Arguments:
468e9923e8dSJed Brown +  fine - SNES from which to restrict
469e9923e8dSJed Brown -  Xfine - vector to restrict
470e9923e8dSJed Brown 
471e9923e8dSJed Brown    Output Arguments:
472e9923e8dSJed Brown .  Xcoarse - result of restriction
473e9923e8dSJed Brown 
474e9923e8dSJed Brown    Level: developer
475e9923e8dSJed Brown 
476e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
477e9923e8dSJed Brown @*/
478e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
479e9923e8dSJed Brown {
480e9923e8dSJed Brown   PetscErrorCode ierr;
481e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
482e9923e8dSJed Brown 
483e9923e8dSJed Brown   PetscFunctionBegin;
484e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
485e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
486e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
487e9923e8dSJed Brown   if (fas->inject) {
488e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
489e9923e8dSJed Brown   } else {
490e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
491e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
492e9923e8dSJed Brown   }
493e9923e8dSJed Brown   PetscFunctionReturn(0);
494e9923e8dSJed Brown }
495e9923e8dSJed Brown 
496e9923e8dSJed Brown #undef __FUNCT__
4978c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection"
49839bd7f45SPeter Brune /*
49939bd7f45SPeter Brune 
50039bd7f45SPeter Brune Performs the FAS coarse correction as:
50139bd7f45SPeter Brune 
50239bd7f45SPeter Brune fine problem: F(x) = 0
50339bd7f45SPeter Brune coarse problem: F^c(x) = b^c
50439bd7f45SPeter Brune 
50539bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
50639bd7f45SPeter Brune 
50739bd7f45SPeter Brune  */
5088c40d5fbSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
50939bd7f45SPeter Brune   PetscErrorCode      ierr;
51039bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
51139bd7f45SPeter Brune   SNESConvergedReason reason;
512ab8d36c9SPeter Brune   SNES                next;
513ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
51439bd7f45SPeter Brune   PetscFunctionBegin;
515ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
516ab8d36c9SPeter Brune   if (next) {
517ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
518ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
519ab8d36c9SPeter Brune 
520ab8d36c9SPeter Brune     X_c  = next->vec_sol;
521ab8d36c9SPeter Brune     Xo_c = next->work[0];
522ab8d36c9SPeter Brune     F_c  = next->vec_func;
523ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
524efe1f98aSPeter Brune 
525938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
526293a7e31SPeter Brune     /* restrict the defect */
527ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
528b9c2fdf1SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
529ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
530e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
531b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
532e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
533ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
534ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
535c90fad12SPeter Brune 
536c90fad12SPeter Brune     /* recurse to the next level */
537e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
538ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
539ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
540742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
541742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
542742fe5e2SPeter Brune       PetscFunctionReturn(0);
543742fe5e2SPeter Brune     }
544fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
545fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
546ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
547293a7e31SPeter Brune   }
54839bd7f45SPeter Brune   PetscFunctionReturn(0);
54939bd7f45SPeter Brune }
55039bd7f45SPeter Brune 
55139bd7f45SPeter Brune #undef __FUNCT__
5522cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive"
55339bd7f45SPeter Brune /*
55439bd7f45SPeter Brune 
55539bd7f45SPeter Brune The additive cycle looks like:
55639bd7f45SPeter Brune 
55707144faaSPeter Brune xhat = x
55807144faaSPeter Brune xhat = dS(x, b)
55907144faaSPeter Brune x = coarsecorrection(xhat, b_d)
56007144faaSPeter Brune x = x + nu*(xhat - x);
56139bd7f45SPeter Brune (optional) x = uS(x, b)
56239bd7f45SPeter Brune 
56339bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
56439bd7f45SPeter Brune 
56539bd7f45SPeter Brune  */
56691f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) {
56707144faaSPeter Brune   Vec                 F, B, Xhat;
56822c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
56939bd7f45SPeter Brune   PetscErrorCode      ierr;
57007144faaSPeter Brune   SNESConvergedReason reason;
57122c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
57222c1e704SPeter Brune   PetscBool           lssuccess;
573ab8d36c9SPeter Brune   SNES                next;
574ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
57539bd7f45SPeter Brune   PetscFunctionBegin;
576ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
57739bd7f45SPeter Brune   F = snes->vec_func;
57839bd7f45SPeter Brune   B = snes->vec_rhs;
579e7f468e7SPeter Brune   Xhat = snes->work[1];
58007144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
58107144faaSPeter Brune   /* recurse first */
582ab8d36c9SPeter Brune   if (next) {
583ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
584ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
58507144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
586c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
587ab8d36c9SPeter Brune     X_c  = next->vec_sol;
588ab8d36c9SPeter Brune     Xo_c = next->work[0];
589ab8d36c9SPeter Brune     F_c  = next->vec_func;
590ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
59139bd7f45SPeter Brune 
592938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
59307144faaSPeter Brune     /* restrict the defect */
594ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
59507144faaSPeter Brune 
59607144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
597ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
598e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
599b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
600e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
60107144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
60207144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
60307144faaSPeter Brune 
60407144faaSPeter Brune     /* recurse */
605e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
606ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
60707144faaSPeter Brune 
60807144faaSPeter Brune     /* smooth on this level */
60991f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
61007144faaSPeter Brune 
611ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
61207144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
61307144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
61407144faaSPeter Brune       PetscFunctionReturn(0);
61507144faaSPeter Brune     }
61607144faaSPeter Brune 
61707144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
618c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
619ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
62007144faaSPeter Brune 
621ddebd997SPeter Brune     /* additive correction of the coarse direction*/
622f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
623f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
6249e764e56SPeter Brune     if (!lssuccess) {
6259e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6269e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6279e764e56SPeter Brune         PetscFunctionReturn(0);
6289e764e56SPeter Brune       }
6299e764e56SPeter Brune     }
630b9c2fdf1SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
63107144faaSPeter Brune   } else {
63291f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
63307144faaSPeter Brune   }
63439bd7f45SPeter Brune   PetscFunctionReturn(0);
63539bd7f45SPeter Brune }
63639bd7f45SPeter Brune 
63739bd7f45SPeter Brune #undef __FUNCT__
6382cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative"
63939bd7f45SPeter Brune /*
64039bd7f45SPeter Brune 
64139bd7f45SPeter Brune Defines the FAS cycle as:
64239bd7f45SPeter Brune 
64339bd7f45SPeter Brune fine problem: F(x) = 0
64439bd7f45SPeter Brune coarse problem: F^c(x) = b^c
64539bd7f45SPeter Brune 
64639bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
64739bd7f45SPeter Brune 
64839bd7f45SPeter Brune correction:
64939bd7f45SPeter Brune 
65039bd7f45SPeter Brune x = x + I(x^c - Rx)
65139bd7f45SPeter Brune 
65239bd7f45SPeter Brune  */
65391f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) {
65439bd7f45SPeter Brune 
65539bd7f45SPeter Brune   PetscErrorCode      ierr;
65639bd7f45SPeter Brune   Vec                 F,B;
65739bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
65839bd7f45SPeter Brune 
65939bd7f45SPeter Brune   PetscFunctionBegin;
66039bd7f45SPeter Brune   F = snes->vec_func;
66139bd7f45SPeter Brune   B = snes->vec_rhs;
66239bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
66391f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
664c90fad12SPeter Brune   if (fas->level != 0) {
6658c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
66691f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
667fe6f9142SPeter Brune   }
668fa9694d7SPeter Brune   PetscFunctionReturn(0);
669421d9b32SPeter Brune }
670421d9b32SPeter Brune 
671421d9b32SPeter Brune #undef __FUNCT__
672421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
673421d9b32SPeter Brune 
674421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
675421d9b32SPeter Brune {
676fa9694d7SPeter Brune   PetscErrorCode ierr;
677fe6f9142SPeter Brune   PetscInt       i, maxits;
678ddb5aff1SPeter Brune   Vec            X, F;
679fe6f9142SPeter Brune   PetscReal      fnorm;
680b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
681b17ce1afSJed Brown   DM             dm;
682e70c42e5SPeter Brune   PetscBool      isFine;
683b17ce1afSJed Brown 
684421d9b32SPeter Brune   PetscFunctionBegin;
685fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
686fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
687fa9694d7SPeter Brune   X = snes->vec_sol;
688f5a6d4f9SBarry Smith   F = snes->vec_func;
689293a7e31SPeter Brune 
69018a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
691293a7e31SPeter Brune   /*norm setup */
692fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
693fe6f9142SPeter Brune   snes->iter = 0;
694fe6f9142SPeter Brune   snes->norm = 0.;
695fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
696e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
697fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
698fe6f9142SPeter Brune     if (snes->domainerror) {
699fe6f9142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
700fe6f9142SPeter Brune       PetscFunctionReturn(0);
701fe6f9142SPeter Brune     }
702e4ed7901SPeter Brune   } else {
703e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
704e4ed7901SPeter Brune   }
705e4ed7901SPeter Brune 
706e4ed7901SPeter Brune   if (!snes->norm_init_set) {
707fe6f9142SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
708fe6f9142SPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
709fe6f9142SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
710e4ed7901SPeter Brune   } else {
711e4ed7901SPeter Brune     fnorm = snes->norm_init;
712e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
713e4ed7901SPeter Brune   }
714e4ed7901SPeter Brune 
715fe6f9142SPeter Brune   snes->norm = fnorm;
716fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
717fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
718fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
719fe6f9142SPeter Brune 
720fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
721fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
722fe6f9142SPeter Brune   /* test convergence */
723fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
724fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
725e4ed7901SPeter Brune 
726b17ce1afSJed Brown 
727b9c2fdf1SPeter Brune   if (isFine) {
728b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
729b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
730b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
731b17ce1afSJed Brown       DM dmcoarse;
732b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
733b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
734b17ce1afSJed Brown       dm = dmcoarse;
735b17ce1afSJed Brown     }
736b9c2fdf1SPeter Brune   }
737b17ce1afSJed Brown 
738fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
739fe6f9142SPeter Brune     /* Call general purpose update function */
740646217ecSPeter Brune 
741fe6f9142SPeter Brune     if (snes->ops->update) {
742fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
743fe6f9142SPeter Brune     }
74407144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
74591f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
74607144faaSPeter Brune     } else {
74791f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
74807144faaSPeter Brune     }
749742fe5e2SPeter Brune 
750742fe5e2SPeter Brune     /* check for FAS cycle divergence */
751742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
752742fe5e2SPeter Brune       PetscFunctionReturn(0);
753742fe5e2SPeter Brune     }
754b9c2fdf1SPeter Brune 
755c90fad12SPeter Brune     /* Monitor convergence */
756c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
757c90fad12SPeter Brune     snes->iter = i+1;
758c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
759c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
760c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
761c90fad12SPeter Brune     /* Test for convergence */
76266585501SPeter Brune     if (isFine) {
763b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
764c90fad12SPeter Brune       if (snes->reason) break;
765fe6f9142SPeter Brune     }
76666585501SPeter Brune   }
767fe6f9142SPeter Brune   if (i == maxits) {
768fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
769fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
770fe6f9142SPeter Brune   }
771421d9b32SPeter Brune   PetscFunctionReturn(0);
772421d9b32SPeter Brune }
773