xref: /petsc/src/snes/impls/fas/fas.c (revision f89ba88e4820770605f92f5e73984a1592b2e9a1)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3 
4 const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0};
5 
6 extern PetscErrorCode SNESDestroy_FAS(SNES snes);
7 extern PetscErrorCode SNESSetUp_FAS(SNES snes);
8 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
9 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
10 extern PetscErrorCode SNESSolve_FAS(SNES snes);
11 extern PetscErrorCode SNESReset_FAS(SNES snes);
12 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
13 extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *);
14 
15 /*MC
16 
17 SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
18 
19    The nonlinear problem is solved by correction using coarse versions
20    of the nonlinear problem.  This problem is perturbed so that a projected
21    solution of the fine problem elicits no correction from the coarse problem.
22 
23 Options Database:
24 +   -snes_fas_levels -  The number of levels
25 .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
26 .   -snes_fas_type<additive, multiplicative>  -  Additive or multiplicative cycle
27 .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
28 .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
29 .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
30 .   -snes_fas_monitor -  Monitor progress of all of the levels
31 .   -fas_levels_snes_ -  SNES options for all smoothers
32 .   -fas_levels_cycle_snes_ -  SNES options for all cycles
33 .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
34 .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
35 -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
36 
37 Notes:
38    The organization of the FAS solver is slightly different from the organization of PCMG
39    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
40    The cycle SNES instance may be used for monitoring convergence on a particular level.
41 
42 Level: beginner
43 
44 .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
45 M*/
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "SNESCreate_FAS"
49 PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes)
50 {
51   SNES_FAS * fas;
52   PetscErrorCode ierr;
53 
54   PetscFunctionBegin;
55   snes->ops->destroy        = SNESDestroy_FAS;
56   snes->ops->setup          = SNESSetUp_FAS;
57   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
58   snes->ops->view           = SNESView_FAS;
59   snes->ops->solve          = SNESSolve_FAS;
60   snes->ops->reset          = SNESReset_FAS;
61 
62   snes->usesksp             = PETSC_FALSE;
63   snes->usespc              = PETSC_FALSE;
64 
65   if (!snes->tolerancesset) {
66     snes->max_funcs = 30000;
67     snes->max_its   = 10000;
68   }
69 
70   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
71   snes->data                  = (void*) fas;
72   fas->level                  = 0;
73   fas->levels                 = 1;
74   fas->n_cycles               = 1;
75   fas->max_up_it              = 1;
76   fas->max_down_it            = 1;
77   fas->smoothu                = PETSC_NULL;
78   fas->smoothd                = PETSC_NULL;
79   fas->next                   = PETSC_NULL;
80   fas->previous               = PETSC_NULL;
81   fas->fine                   = snes;
82   fas->interpolate            = PETSC_NULL;
83   fas->restrct                = PETSC_NULL;
84   fas->inject                 = PETSC_NULL;
85   fas->monitor                = PETSC_NULL;
86   fas->usedmfornumberoflevels = PETSC_FALSE;
87   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "SNESReset_FAS"
93 PetscErrorCode SNESReset_FAS(SNES snes)
94 {
95   PetscErrorCode ierr = 0;
96   SNES_FAS * fas = (SNES_FAS *)snes->data;
97 
98   PetscFunctionBegin;
99   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
100   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
101   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
102   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
103   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
104   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
105   if (fas->next)      ierr = SNESReset(fas->next);CHKERRQ(ierr);
106 
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "SNESDestroy_FAS"
112 PetscErrorCode SNESDestroy_FAS(SNES snes)
113 {
114   SNES_FAS * fas = (SNES_FAS *)snes->data;
115   PetscErrorCode ierr = 0;
116 
117   PetscFunctionBegin;
118   /* recursively resets and then destroys */
119   ierr = SNESReset(snes);CHKERRQ(ierr);
120   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
121   ierr = PetscFree(fas);CHKERRQ(ierr);
122 
123   PetscFunctionReturn(0);
124 }
125 
126 #undef __FUNCT__
127 #define __FUNCT__ "SNESSetUp_FAS"
128 PetscErrorCode SNESSetUp_FAS(SNES snes)
129 {
130   SNES_FAS                    *fas = (SNES_FAS *) snes->data;
131   PetscErrorCode              ierr;
132   VecScatter                  injscatter;
133   PetscInt                    dm_levels;
134   Vec                         vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
135   SNES                        next;
136   PetscBool                   isFine;
137   SNESLineSearch              linesearch;
138   SNESLineSearch              slinesearch;
139   void                        *lsprectx,*lspostctx;
140   SNESLineSearchPreCheckFunc  lsprefunc;
141   SNESLineSearchPostCheckFunc lspostfunc;
142   PetscFunctionBegin;
143 
144   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
145   if (fas->usedmfornumberoflevels && isFine) {
146     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
147     dm_levels++;
148     if (dm_levels > fas->levels) {
149       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
150       vec_sol = snes->vec_sol;
151       vec_func = snes->vec_func;
152       vec_sol_update = snes->vec_sol_update;
153       vec_rhs = snes->vec_rhs;
154       snes->vec_sol = PETSC_NULL;
155       snes->vec_func = PETSC_NULL;
156       snes->vec_sol_update = PETSC_NULL;
157       snes->vec_rhs = PETSC_NULL;
158 
159       /* reset the number of levels */
160       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
161       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
162 
163       snes->vec_sol = vec_sol;
164       snes->vec_func = vec_func;
165       snes->vec_rhs = vec_rhs;
166       snes->vec_sol_update = vec_sol_update;
167     }
168   }
169   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
170   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
171 
172   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
173     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
174   } else {
175     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
176   }
177 
178   /* set up the smoothers if they haven't already been set up */
179   if (!fas->smoothd) {
180     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
181   }
182   if (!fas->smoothu && fas->level != fas->levels - 1) {
183     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
184   }
185 
186   if (snes->dm) {
187     /* set the smoother DMs properly */
188     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
189     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
190     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
191     if (next) {
192       /* for now -- assume the DM and the evaluation functions have been set externally */
193       if (!next->dm) {
194         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
195         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
196       }
197       /* set the interpolation and restriction from the DM */
198       if (!fas->interpolate) {
199         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
200         if (!fas->restrct) {
201           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
202           fas->restrct = fas->interpolate;
203         }
204       }
205       /* set the injection from the DM */
206       if (!fas->inject) {
207         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
208         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
209         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
210       }
211     }
212   }
213   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
214   if (fas->galerkin) {
215     if (next)
216       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
217     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
218     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
219   }
220 
221   /* sets the down (pre) smoother's default norm and sets it from options */
222   if(fas->smoothd){
223     if (fas->level == 0) {
224       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
225      } else {
226       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
227     }
228     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
229     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
230     ierr = SNESGetSNESLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
231     ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr);
232     ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr);
233     ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr);
234     ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr);
235     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
236   }
237 
238   /* sets the up (post) smoother's default norm and sets it from options */
239   if(fas->smoothu){
240     if (fas->level != fas->levels - 1) {
241       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
242     } else {
243       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
244     }
245     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
246     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
247     ierr = SNESGetSNESLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
248     ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr);
249     ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr);
250     ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr);
251     ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr);
252     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
253   }
254 
255   if (next) {
256     /* gotta set up the solution vector for this to work */
257     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
258     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
259     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
260     ierr = SNESGetSNESLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
261     ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr);
262     ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr);
263     ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr);
264     ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr);
265     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
266     ierr = SNESSetUp(next);CHKERRQ(ierr);
267   }
268   /* setup FAS work vectors */
269   if (fas->galerkin) {
270     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
271     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
272   }
273   PetscFunctionReturn(0);
274 }
275 
276 #undef __FUNCT__
277 #define __FUNCT__ "SNESSetFromOptions_FAS"
278 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
279 {
280   SNES_FAS       *fas = (SNES_FAS *) snes->data;
281   PetscInt       levels = 1;
282   PetscBool      flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE;
283   PetscErrorCode ierr;
284   char           monfilename[PETSC_MAX_PATH_LEN];
285   SNESFASType    fastype;
286   const char     *optionsprefix;
287   SNESLineSearch linesearch;
288   PetscInt       m, n_up, n_down;
289   SNES           next;
290   PetscBool      isFine;
291 
292   PetscFunctionBegin;
293   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
294   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
295 
296   /* number of levels -- only process most options on the finest level */
297   if (isFine) {
298     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
299     if (!flg && snes->dm) {
300       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
301       levels++;
302       fas->usedmfornumberoflevels = PETSC_TRUE;
303     }
304     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
305     fastype = fas->fastype;
306     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
307     if (flg) {
308       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
309     }
310 
311     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
312     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
313     if (flg) {
314       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
315     }
316 
317     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
318     if (flg) {
319       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
320     }
321 
322     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
323 
324     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
325 
326     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
327     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
328   }
329 
330   ierr = PetscOptionsTail();CHKERRQ(ierr);
331   /* setup from the determined types if there is no pointwise procedure or smoother defined */
332   if (upflg) {
333     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
334   }
335   if (downflg) {
336     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
337   }
338 
339   /* set up the default line search for coarse grid corrections */
340   if (fas->fastype == SNES_FAS_ADDITIVE) {
341     if (!snes->linesearch) {
342       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
343       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
344     }
345   }
346 
347   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
348   /* recursive option setting for the smoothers */
349   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNCT__
354 #define __FUNCT__ "SNESView_FAS"
355 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
356 {
357   SNES_FAS       *fas = (SNES_FAS *) snes->data;
358   PetscBool      isFine, iascii;
359   PetscInt       i;
360   PetscErrorCode ierr;
361   SNES           smoothu, smoothd, levelsnes;
362 
363   PetscFunctionBegin;
364   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
365   if (isFine) {
366     ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
367     if (iascii) {
368       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
369       if (fas->galerkin) {
370         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
371       } else {
372         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
373       }
374       for (i=0; i<fas->levels; i++) {
375         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
376         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
377         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
378         if (!i) {
379           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
380         } else {
381           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
382         }
383         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
384         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
385         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
386         if (i && (smoothd == smoothu)) {
387           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
388         } else if (i) {
389           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
390           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
391           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
392           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
393         }
394       }
395     } else {
396       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
397     }
398   }
399   PetscFunctionReturn(0);
400 }
401 
402 #undef __FUNCT__
403 #define __FUNCT__ "SNESFASDownSmooth_Private"
404 /*
405 Defines the action of the downsmoother
406  */
407 PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
408 {
409   PetscErrorCode      ierr = 0;
410   SNESConvergedReason reason;
411   Vec                 FPC;
412   SNES                smoothd;
413   PetscFunctionBegin;
414   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
415   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
416   ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr);
417   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
418   /* check convergence reason for the smoother */
419   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
420   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
421     snes->reason = SNES_DIVERGED_INNER;
422     PetscFunctionReturn(0);
423   }
424   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
425   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
426   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
427   PetscFunctionReturn(0);
428 }
429 
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "SNESFASUpSmooth_Private"
433 /*
434 Defines the action of the upsmoother
435  */
436 PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) {
437   PetscErrorCode      ierr = 0;
438   SNESConvergedReason reason;
439   Vec                 FPC;
440   SNES                smoothu;
441   PetscFunctionBegin;
442 
443   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
444   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
445   /* check convergence reason for the smoother */
446   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
447   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
448     snes->reason = SNES_DIVERGED_INNER;
449     PetscFunctionReturn(0);
450   }
451   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
452   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
453   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
454   PetscFunctionReturn(0);
455 }
456 
457 #undef __FUNCT__
458 #define __FUNCT__ "SNESFASCreateCoarseVec"
459 /*@
460    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
461 
462    Collective
463 
464    Input Arguments:
465 .  snes - SNESFAS
466 
467    Output Arguments:
468 .  Xcoarse - vector on level one coarser than snes
469 
470    Level: developer
471 
472 .seealso: SNESFASSetRestriction(), SNESFASRestrict()
473 @*/
474 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
475 {
476   PetscErrorCode ierr;
477   SNES_FAS       *fas = (SNES_FAS*)snes->data;
478 
479   PetscFunctionBegin;
480   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
481   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
482   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
483   PetscFunctionReturn(0);
484 }
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "SNESFASRestrict"
488 /*@
489    SNESFASRestrict - restrict a Vec to the next coarser level
490 
491    Collective
492 
493    Input Arguments:
494 +  fine - SNES from which to restrict
495 -  Xfine - vector to restrict
496 
497    Output Arguments:
498 .  Xcoarse - result of restriction
499 
500    Level: developer
501 
502 .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
503 @*/
504 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
505 {
506   PetscErrorCode ierr;
507   SNES_FAS       *fas = (SNES_FAS*)fine->data;
508 
509   PetscFunctionBegin;
510   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
511   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
512   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
513   if (fas->inject) {
514     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
515   } else {
516     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
517     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
518   }
519   PetscFunctionReturn(0);
520 }
521 
522 #undef __FUNCT__
523 #define __FUNCT__ "SNESFASCoarseCorrection"
524 /*
525 
526 Performs the FAS coarse correction as:
527 
528 fine problem: F(x) = 0
529 coarse problem: F^c(x) = b^c
530 
531 b^c = F^c(I^c_fx^f - I^c_fF(x))
532 
533  */
534 PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
535   PetscErrorCode      ierr;
536   Vec                 X_c, Xo_c, F_c, B_c;
537   SNESConvergedReason reason;
538   SNES                next;
539   Mat                 restrct, interpolate;
540   PetscFunctionBegin;
541   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
542   if (next) {
543     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
544     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
545 
546     X_c  = next->vec_sol;
547     Xo_c = next->work[0];
548     F_c  = next->vec_func;
549     B_c  = next->vec_rhs;
550 
551     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
552     /* restrict the defect */
553     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
554     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
555     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
556     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
557     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
558     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
559     /* set initial guess of the coarse problem to the projected fine solution */
560     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
561 
562     /* recurse to the next level */
563     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
564     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
565     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
566     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
567       snes->reason = SNES_DIVERGED_INNER;
568       PetscFunctionReturn(0);
569     }
570     /* correct as x <- x + I(x^c - Rx)*/
571     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
572     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
573   }
574   PetscFunctionReturn(0);
575 }
576 
577 #undef __FUNCT__
578 #define __FUNCT__ "SNESFASCycle_Additive"
579 /*
580 
581 The additive cycle looks like:
582 
583 xhat = x
584 xhat = dS(x, b)
585 x = coarsecorrection(xhat, b_d)
586 x = x + nu*(xhat - x);
587 (optional) x = uS(x, b)
588 
589 With the coarse RHS (defect correction) as below.
590 
591  */
592 PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) {
593   Vec                 F, B, Xhat;
594   Vec                 X_c, Xo_c, F_c, B_c;
595   PetscErrorCode      ierr;
596   SNESConvergedReason reason;
597   PetscReal           xnorm, fnorm, ynorm;
598   PetscBool           lssuccess;
599   SNES                next;
600   Mat                 restrct, interpolate;
601   PetscFunctionBegin;
602   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
603   F = snes->vec_func;
604   B = snes->vec_rhs;
605   Xhat = snes->work[1];
606   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
607   /* recurse first */
608   if (next) {
609     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
610     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
611     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
612     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
613     X_c  = next->vec_sol;
614     Xo_c = next->work[0];
615     F_c  = next->vec_func;
616     B_c  = next->vec_rhs;
617 
618     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
619     /* restrict the defect */
620     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
621 
622     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
623     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
624     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
625     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
626     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
627     /* set initial guess of the coarse problem to the projected fine solution */
628     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
629 
630     /* recurse */
631     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
632     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
633 
634     /* smooth on this level */
635     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
636 
637     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
638     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
639       snes->reason = SNES_DIVERGED_INNER;
640       PetscFunctionReturn(0);
641     }
642 
643     /* correct as x <- x + I(x^c - Rx)*/
644     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
645     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
646 
647     /* additive correction of the coarse direction*/
648     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
649     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
650     if (!lssuccess) {
651       if (++snes->numFailures >= snes->maxFailures) {
652         snes->reason = SNES_DIVERGED_LINE_SEARCH;
653         PetscFunctionReturn(0);
654       }
655     }
656     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
657   } else {
658     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
659   }
660   PetscFunctionReturn(0);
661 }
662 
663 #undef __FUNCT__
664 #define __FUNCT__ "SNESFASCycle_Multiplicative"
665 /*
666 
667 Defines the FAS cycle as:
668 
669 fine problem: F(x) = 0
670 coarse problem: F^c(x) = b^c
671 
672 b^c = F^c(I^c_fx^f - I^c_fF(x))
673 
674 correction:
675 
676 x = x + I(x^c - Rx)
677 
678  */
679 PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) {
680 
681   PetscErrorCode      ierr;
682   Vec                 F,B;
683   SNES_FAS            *fas = (SNES_FAS *)snes->data;
684 
685   PetscFunctionBegin;
686   F = snes->vec_func;
687   B = snes->vec_rhs;
688   /* pre-smooth -- just update using the pre-smoother */
689   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
690   if (fas->level != 0) {
691     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
692     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
693   }
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "SNESSolve_FAS"
699 
700 PetscErrorCode SNESSolve_FAS(SNES snes)
701 {
702   PetscErrorCode ierr;
703   PetscInt       i, maxits;
704   Vec            X, F;
705   PetscReal      fnorm;
706   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
707   DM             dm;
708   PetscBool      isFine;
709 
710   PetscFunctionBegin;
711   maxits = snes->max_its;            /* maximum number of iterations */
712   snes->reason = SNES_CONVERGED_ITERATING;
713   X = snes->vec_sol;
714   F = snes->vec_func;
715 
716   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
717   /*norm setup */
718   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
719   snes->iter = 0;
720   snes->norm = 0.;
721   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
722   if (!snes->vec_func_init_set) {
723     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
724     if (snes->domainerror) {
725       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
726       PetscFunctionReturn(0);
727     }
728   } else {
729     snes->vec_func_init_set = PETSC_FALSE;
730   }
731 
732   if (!snes->norm_init_set) {
733     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
734     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
735     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
736   } else {
737     fnorm = snes->norm_init;
738     snes->norm_init_set = PETSC_FALSE;
739   }
740 
741   snes->norm = fnorm;
742   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
743   SNESLogConvHistory(snes,fnorm,0);
744   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
745 
746   /* set parameter for default relative tolerance convergence test */
747   snes->ttol = fnorm*snes->rtol;
748   /* test convergence */
749   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
750   if (snes->reason) PetscFunctionReturn(0);
751 
752 
753   if (isFine) {
754     /* propagate scale-dependent data up the hierarchy */
755     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
756     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
757       DM dmcoarse;
758       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
759       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
760       dm = dmcoarse;
761     }
762   }
763 
764   for (i = 0; i < maxits; i++) {
765     /* Call general purpose update function */
766 
767     if (snes->ops->update) {
768       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
769     }
770     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
771       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
772     } else {
773       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
774     }
775 
776     /* check for FAS cycle divergence */
777     if (snes->reason != SNES_CONVERGED_ITERATING) {
778       PetscFunctionReturn(0);
779     }
780 
781     /* Monitor convergence */
782     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
783     snes->iter = i+1;
784     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
785     SNESLogConvHistory(snes,snes->norm,0);
786     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
787     /* Test for convergence */
788     if (isFine) {
789       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
790       if (snes->reason) break;
791     }
792   }
793   if (i == maxits) {
794     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
795     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
796   }
797   PetscFunctionReturn(0);
798 }
799