xref: /petsc/src/snes/impls/fas/fas.c (revision 6dbf99738a075cf632ae62c52b954b8df5871b59)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnes.h"  I*/
3 
4 const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","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,full,kaskade>  -  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 .   -snes_fas_full_downsweepsmooth<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
32 .   -fas_levels_snes_ -  SNES options for all smoothers
33 .   -fas_levels_cycle_snes_ -  SNES options for all cycles
34 .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
35 .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
36 -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
37 
38 Notes:
39    The organization of the FAS solver is slightly different from the organization of PCMG
40    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
41    The cycle SNES instance may be used for monitoring convergence on a particular level.
42 
43 Level: beginner
44 
45 .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
46 M*/
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "SNESCreate_FAS"
50 PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
51 {
52   SNES_FAS       *fas;
53   PetscErrorCode ierr;
54 
55   PetscFunctionBegin;
56   snes->ops->destroy        = SNESDestroy_FAS;
57   snes->ops->setup          = SNESSetUp_FAS;
58   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
59   snes->ops->view           = SNESView_FAS;
60   snes->ops->solve          = SNESSolve_FAS;
61   snes->ops->reset          = SNESReset_FAS;
62 
63   snes->usesksp = PETSC_FALSE;
64   snes->usespc  = PETSC_FALSE;
65 
66   if (!snes->tolerancesset) {
67     snes->max_funcs = 30000;
68     snes->max_its   = 10000;
69   }
70 
71   ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr);
72 
73   snes->data                  = (void*) fas;
74   fas->level                  = 0;
75   fas->levels                 = 1;
76   fas->n_cycles               = 1;
77   fas->max_up_it              = 1;
78   fas->max_down_it            = 1;
79   fas->smoothu                = NULL;
80   fas->smoothd                = NULL;
81   fas->next                   = NULL;
82   fas->previous               = NULL;
83   fas->fine                   = snes;
84   fas->interpolate            = NULL;
85   fas->restrct                = NULL;
86   fas->inject                 = NULL;
87   fas->monitor                = NULL;
88   fas->usedmfornumberoflevels = PETSC_FALSE;
89   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
90   fas->full_downsweep         = PETSC_FALSE;
91 
92   fas->eventsmoothsetup    = 0;
93   fas->eventsmoothsolve    = 0;
94   fas->eventresidual       = 0;
95   fas->eventinterprestrict = 0;
96   PetscFunctionReturn(0);
97 }
98 
99 #undef __FUNCT__
100 #define __FUNCT__ "SNESReset_FAS"
101 PetscErrorCode SNESReset_FAS(SNES snes)
102 {
103   PetscErrorCode ierr  = 0;
104   SNES_FAS       * fas = (SNES_FAS*)snes->data;
105 
106   PetscFunctionBegin;
107   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
108   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
109   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
110   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
111   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
112   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
113   if (fas->galerkin) {
114     ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr);
115     ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr);
116   }
117   if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);}
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "SNESDestroy_FAS"
123 PetscErrorCode SNESDestroy_FAS(SNES snes)
124 {
125   SNES_FAS       * fas = (SNES_FAS*)snes->data;
126   PetscErrorCode ierr  = 0;
127 
128   PetscFunctionBegin;
129   /* recursively resets and then destroys */
130   ierr = SNESReset(snes);CHKERRQ(ierr);
131   if (fas->next) {
132     ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
133   }
134   ierr = PetscFree(fas);CHKERRQ(ierr);
135   PetscFunctionReturn(0);
136 }
137 
138 #undef __FUNCT__
139 #define __FUNCT__ "SNESSetUp_FAS"
140 PetscErrorCode SNESSetUp_FAS(SNES snes)
141 {
142   SNES_FAS       *fas = (SNES_FAS*) snes->data;
143   PetscErrorCode ierr;
144   VecScatter     injscatter;
145   PetscInt       dm_levels;
146   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
147   SNES           next;
148   PetscBool      isFine;
149   SNESLineSearch linesearch;
150   SNESLineSearch slinesearch;
151   void           *lsprectx,*lspostctx;
152   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
153   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
154 
155   PetscFunctionBegin;
156   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
157   if (fas->usedmfornumberoflevels && isFine) {
158     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
159     dm_levels++;
160     if (dm_levels > fas->levels) {
161       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
162       vec_sol              = snes->vec_sol;
163       vec_func             = snes->vec_func;
164       vec_sol_update       = snes->vec_sol_update;
165       vec_rhs              = snes->vec_rhs;
166       snes->vec_sol        = NULL;
167       snes->vec_func       = NULL;
168       snes->vec_sol_update = NULL;
169       snes->vec_rhs        = NULL;
170 
171       /* reset the number of levels */
172       ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr);
173       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
174 
175       snes->vec_sol        = vec_sol;
176       snes->vec_func       = vec_func;
177       snes->vec_rhs        = vec_rhs;
178       snes->vec_sol_update = vec_sol_update;
179     }
180   }
181   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
182   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
183 
184   ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
185 
186   /* set up the smoothers if they haven't already been set up */
187   if (!fas->smoothd) {
188     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
189   }
190 
191   if (snes->dm) {
192     /* set the smoother DMs properly */
193     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
194     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
195     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
196     if (next) {
197       /* for now -- assume the DM and the evaluation functions have been set externally */
198       if (!next->dm) {
199         ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr);
200         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
201       }
202       /* set the interpolation and restriction from the DM */
203       if (!fas->interpolate) {
204         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
205         if (!fas->restrct) {
206           ierr         = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
207           fas->restrct = fas->interpolate;
208         }
209       }
210       /* set the injection from the DM */
211       if (!fas->inject) {
212         ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr);
213       }
214     }
215   }
216   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
217   if (fas->galerkin) {
218     if (next) {
219       ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
220     }
221     if (fas->smoothd && fas->level != fas->levels - 1) {
222       ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
223     }
224     if (fas->smoothu && fas->level != fas->levels - 1) {
225       ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
226     }
227   }
228 
229   /* sets the down (pre) smoother's default norm and sets it from options */
230   if (fas->smoothd) {
231     if (fas->level == 0 && fas->levels != 1) {
232       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
233     } else {
234       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
235     }
236     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr);
237     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
238     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
239     ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
240     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
241     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
242     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
243     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
244     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
245 
246     fas->smoothd->vec_sol        = snes->vec_sol;
247     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
248     fas->smoothd->vec_sol_update = snes->vec_sol_update;
249     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
250     fas->smoothd->vec_func       = snes->vec_func;
251     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
252 
253     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
254     ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr);
255     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
256   }
257 
258   /* sets the up (post) smoother's default norm and sets it from options */
259   if (fas->smoothu) {
260     if (fas->level != fas->levels - 1) {
261       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
262     } else {
263       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
264     }
265     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr);
266     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
267     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
268     ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
269     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
270     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
271     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
272     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
273     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
274 
275     fas->smoothu->vec_sol        = snes->vec_sol;
276     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
277     fas->smoothu->vec_sol_update = snes->vec_sol_update;
278     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
279     fas->smoothu->vec_func       = snes->vec_func;
280     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
281 
282     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
283     ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr);
284     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
285 
286   }
287 
288   if (next) {
289     /* gotta set up the solution vector for this to work */
290     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
291     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
292     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr);
293     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
294     ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
295     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
296     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
297     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
298     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
299     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
300     ierr = SNESSetUp(next);CHKERRQ(ierr);
301   }
302   /* setup FAS work vectors */
303   if (fas->galerkin) {
304     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
305     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
306   }
307   PetscFunctionReturn(0);
308 }
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "SNESSetFromOptions_FAS"
312 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
313 {
314   SNES_FAS       *fas   = (SNES_FAS*) snes->data;
315   PetscInt       levels = 1;
316   PetscBool      flg    = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE;
317   PetscErrorCode ierr;
318   char           monfilename[PETSC_MAX_PATH_LEN];
319   SNESFASType    fastype;
320   const char     *optionsprefix;
321   SNESLineSearch linesearch;
322   PetscInt       m, n_up, n_down;
323   SNES           next;
324   PetscBool      isFine;
325 
326   PetscFunctionBegin;
327   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
328   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
329 
330   /* number of levels -- only process most options on the finest level */
331   if (isFine) {
332     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
333     if (!flg && snes->dm) {
334       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
335       levels++;
336       fas->usedmfornumberoflevels = PETSC_TRUE;
337     }
338     ierr    = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr);
339     fastype = fas->fastype;
340     ierr    = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
341     if (flg) {
342       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
343     }
344 
345     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
346     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
347     if (flg) {
348       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
349     }
350     ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr);
351     if (flg) {
352       ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr);
353     }
354 
355     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
356     if (flg) {
357       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
358     }
359 
360     if (fas->fastype == SNES_FAS_FULL) {
361       ierr   = PetscOptionsBool("-snes_fas_full_downsweep","Smooth on the initial upsweep for full FAS cycles","SNESFASFullSetDownSweep",fas->full_downsweep,&fas->full_downsweep,&flg);CHKERRQ(ierr);
362       if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);}
363     }
364 
365     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
366 
367     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
368 
369     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
370     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
371 
372     flg    = PETSC_FALSE;
373     monflg = PETSC_TRUE;
374     ierr   = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
375     if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
376   }
377 
378   ierr = PetscOptionsTail();CHKERRQ(ierr);
379   /* setup from the determined types if there is no pointwise procedure or smoother defined */
380   if (upflg) {
381     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
382   }
383   if (downflg) {
384     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
385   }
386 
387   /* set up the default line search for coarse grid corrections */
388   if (fas->fastype == SNES_FAS_ADDITIVE) {
389     if (!snes->linesearch) {
390       ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
391       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
392     }
393   }
394 
395   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
396   /* recursive option setting for the smoothers */
397   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
398   PetscFunctionReturn(0);
399 }
400 
401 #include <petscdraw.h>
402 #undef __FUNCT__
403 #define __FUNCT__ "SNESView_FAS"
404 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
405 {
406   SNES_FAS       *fas = (SNES_FAS*) snes->data;
407   PetscBool      isFine,iascii,isdraw;
408   PetscInt       i;
409   PetscErrorCode ierr;
410   SNES           smoothu, smoothd, levelsnes;
411 
412   PetscFunctionBegin;
413   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
414   if (isFine) {
415     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
416     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
417     if (iascii) {
418       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
419       if (fas->galerkin) {
420         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
421       } else {
422         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
423       }
424       for (i=0; i<fas->levels; i++) {
425         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
426         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
427         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
428         if (!i) {
429           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
430         } else {
431           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
432         }
433         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
434         if (smoothd) {
435           ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
436         } else {
437           ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
438         }
439         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
440         if (i && (smoothd == smoothu)) {
441           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
442         } else if (i) {
443           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
444           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
445           if (smoothu) {
446             ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
447           } else {
448             ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
449           }
450           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
451         }
452       }
453     } else if (isdraw) {
454       PetscDraw draw;
455       PetscReal x,w,y,bottom,th,wth;
456       SNES_FAS  *curfas = fas;
457       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
458       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
459       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
460       bottom = y - th;
461       while (curfas) {
462         if (!curfas->smoothu) {
463           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
464           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
465           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
466         } else {
467           w    = 0.5*PetscMin(1.0-x,x);
468           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
469           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
470           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
471           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
472           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
473           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
474         }
475         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
476         bottom -= 5*th;
477         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
478         else curfas = NULL;
479       }
480     }
481   }
482   PetscFunctionReturn(0);
483 }
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "SNESFASDownSmooth_Private"
487 /*
488 Defines the action of the downsmoother
489  */
490 PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
491 {
492   PetscErrorCode      ierr = 0;
493   SNESConvergedReason reason;
494   Vec                 FPC;
495   SNES                smoothd;
496   SNES_FAS            *fas = (SNES_FAS*) snes->data;
497 
498   PetscFunctionBegin;
499   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
500   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
501   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
502   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
503   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
504   /* check convergence reason for the smoother */
505   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
506   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
507     snes->reason = SNES_DIVERGED_INNER;
508     PetscFunctionReturn(0);
509   }
510   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
511   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
512   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
513   PetscFunctionReturn(0);
514 }
515 
516 
517 #undef __FUNCT__
518 #define __FUNCT__ "SNESFASUpSmooth_Private"
519 /*
520 Defines the action of the upsmoother
521  */
522 PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
523 {
524   PetscErrorCode      ierr = 0;
525   SNESConvergedReason reason;
526   Vec                 FPC;
527   SNES                smoothu;
528   SNES_FAS            *fas = (SNES_FAS*) snes->data;
529 
530   PetscFunctionBegin;
531   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
532   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
533   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
534   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
535   /* check convergence reason for the smoother */
536   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
537   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
538     snes->reason = SNES_DIVERGED_INNER;
539     PetscFunctionReturn(0);
540   }
541   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
542   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
543   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
544   PetscFunctionReturn(0);
545 }
546 
547 #undef __FUNCT__
548 #define __FUNCT__ "SNESFASCreateCoarseVec"
549 /*@
550    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
551 
552    Collective
553 
554    Input Arguments:
555 .  snes - SNESFAS
556 
557    Output Arguments:
558 .  Xcoarse - vector on level one coarser than snes
559 
560    Level: developer
561 
562 .seealso: SNESFASSetRestriction(), SNESFASRestrict()
563 @*/
564 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
565 {
566   PetscErrorCode ierr;
567   SNES_FAS       *fas = (SNES_FAS*)snes->data;
568 
569   PetscFunctionBegin;
570   if (fas->rscale) {
571     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
572   } else if (fas->inject) {
573     ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
574   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
575   PetscFunctionReturn(0);
576 }
577 
578 #undef __FUNCT__
579 #define __FUNCT__ "SNESFASRestrict"
580 /*@
581    SNESFASRestrict - restrict a Vec to the next coarser level
582 
583    Collective
584 
585    Input Arguments:
586 +  fine - SNES from which to restrict
587 -  Xfine - vector to restrict
588 
589    Output Arguments:
590 .  Xcoarse - result of restriction
591 
592    Level: developer
593 
594 .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
595 @*/
596 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
597 {
598   PetscErrorCode ierr;
599   SNES_FAS       *fas = (SNES_FAS*)fine->data;
600 
601   PetscFunctionBegin;
602   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
603   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
604   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
605   if (fas->inject) {
606     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
607   } else {
608     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
609     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
610   }
611   PetscFunctionReturn(0);
612 }
613 
614 #undef __FUNCT__
615 #define __FUNCT__ "SNESFASCoarseCorrection"
616 /*
617 
618 Performs the FAS coarse correction as:
619 
620 fine problem:   F(x) = b
621 coarse problem: F^c(x^c) = b^c
622 
623 b^c = F^c(Rx) - R(F(x) - b)
624 
625  */
626 PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
627 {
628   PetscErrorCode      ierr;
629   Vec                 X_c, Xo_c, F_c, B_c;
630   SNESConvergedReason reason;
631   SNES                next;
632   Mat                 restrct, interpolate;
633   SNES_FAS            *fasc;
634 
635   PetscFunctionBegin;
636   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
637   if (next) {
638     fasc = (SNES_FAS*)next->data;
639 
640     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
641     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
642 
643     X_c  = next->vec_sol;
644     Xo_c = next->work[0];
645     F_c  = next->vec_func;
646     B_c  = next->vec_rhs;
647 
648     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
649     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
650     /* restrict the defect: R(F(x) - b) */
651     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
652     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
653 
654     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
655     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
656     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
657     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
658 
659     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
660     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
661     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
662     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
663     /* set initial guess of the coarse problem to the projected fine solution */
664     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
665 
666     /* recurse to the next level */
667     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
668     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
669     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
670     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
671       snes->reason = SNES_DIVERGED_INNER;
672       PetscFunctionReturn(0);
673     }
674     /* correct as x <- x + I(x^c - Rx)*/
675     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
676 
677     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
678     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
679     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
680   }
681   PetscFunctionReturn(0);
682 }
683 
684 #undef __FUNCT__
685 #define __FUNCT__ "SNESFASCycle_Additive"
686 /*
687 
688 The additive cycle looks like:
689 
690 xhat = x
691 xhat = dS(x, b)
692 x = coarsecorrection(xhat, b_d)
693 x = x + nu*(xhat - x);
694 (optional) x = uS(x, b)
695 
696 With the coarse RHS (defect correction) as below.
697 
698  */
699 PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
700 {
701   Vec                 F, B, Xhat;
702   Vec                 X_c, Xo_c, F_c, B_c;
703   PetscErrorCode      ierr;
704   SNESConvergedReason reason;
705   PetscReal           xnorm, fnorm, ynorm;
706   PetscBool           lssuccess;
707   SNES                next;
708   Mat                 restrct, interpolate;
709   SNES_FAS            *fas = (SNES_FAS*)snes->data,*fasc;
710 
711   PetscFunctionBegin;
712   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
713   F    = snes->vec_func;
714   B    = snes->vec_rhs;
715   Xhat = snes->work[1];
716   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
717   /* recurse first */
718   if (next) {
719     fasc = (SNES_FAS*)next->data;
720     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
721     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
722     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
723     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
724     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
725     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
726     X_c  = next->vec_sol;
727     Xo_c = next->work[0];
728     F_c  = next->vec_func;
729     B_c  = next->vec_rhs;
730 
731     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
732     /* restrict the defect */
733     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
734 
735     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
736     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
737     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
738     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
739     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
740     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
741     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
742     /* set initial guess of the coarse problem to the projected fine solution */
743     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
744 
745     /* recurse */
746     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
747     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
748 
749     /* smooth on this level */
750     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
751 
752     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
753     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
754       snes->reason = SNES_DIVERGED_INNER;
755       PetscFunctionReturn(0);
756     }
757 
758     /* correct as x <- x + I(x^c - Rx)*/
759     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
760     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
761 
762     /* additive correction of the coarse direction*/
763     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
764     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
765     if (!lssuccess) {
766       if (++snes->numFailures >= snes->maxFailures) {
767         snes->reason = SNES_DIVERGED_LINE_SEARCH;
768         PetscFunctionReturn(0);
769       }
770     }
771     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
772   } else {
773     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
774   }
775   PetscFunctionReturn(0);
776 }
777 
778 #undef __FUNCT__
779 #define __FUNCT__ "SNESFASCycle_Multiplicative"
780 /*
781 
782 Defines the FAS cycle as:
783 
784 fine problem: F(x) = b
785 coarse problem: F^c(x) = b^c
786 
787 b^c = F^c(Rx) - R(F(x) - b)
788 
789 correction:
790 
791 x = x + I(x^c - Rx)
792 
793  */
794 PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
795 {
796 
797   PetscErrorCode ierr;
798   Vec            F,B;
799   SNES           next;
800 
801   PetscFunctionBegin;
802   F = snes->vec_func;
803   B = snes->vec_rhs;
804   /* pre-smooth -- just update using the pre-smoother */
805   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
806   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
807   if (next) {
808     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
809     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
810   }
811   PetscFunctionReturn(0);
812 }
813 
814 #undef __FUNCT__
815 #define __FUNCT__ "SNESFASCycleSetupPhase_Full"
816 PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
817 {
818   SNES           next;
819   SNES_FAS       *fas = (SNES_FAS*)snes->data;
820   PetscBool      isFine;
821   PetscErrorCode ierr;
822 
823   PetscFunctionBegin;
824   /* pre-smooth -- just update using the pre-smoother */
825   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
826   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
827   fas->full_stage = 0;
828   if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);}
829   PetscFunctionReturn(0);
830 }
831 
832 #undef __FUNCT__
833 #define __FUNCT__ "SNESFASCycle_Full"
834 PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
835 {
836   PetscErrorCode ierr;
837   Vec            F,B;
838   SNES_FAS       *fas = (SNES_FAS*)snes->data;
839   PetscBool      isFine;
840   SNES           next;
841 
842   PetscFunctionBegin;
843   F = snes->vec_func;
844   B = snes->vec_rhs;
845   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
846   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
847 
848   if (isFine) {
849     ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr);
850   }
851 
852   if (fas->full_stage == 0) {
853     /* downsweep */
854     if (next) {
855       if (fas->level != 1) next->max_its += 1;
856       if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
857       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
858       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
859       if (fas->level != 1) next->max_its -= 1;
860     } else {
861       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
862     }
863     fas->full_stage = 1;
864   } else if (fas->full_stage == 1) {
865     if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
866     if (next) {
867       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
868       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
869     }
870   }
871   /* final v-cycle */
872   if (isFine) {
873     if (next) {
874       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
875       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
876     }
877   }
878   PetscFunctionReturn(0);
879 }
880 
881 #undef __FUNCT__
882 #define __FUNCT__ "SNESFASCycle_Kaskade"
883 PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
884 {
885   PetscErrorCode ierr;
886   Vec            F,B;
887   SNES           next;
888 
889   PetscFunctionBegin;
890   F = snes->vec_func;
891   B = snes->vec_rhs;
892   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
893   if (next) {
894     ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
895     ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
896   } else {
897     ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
898   }
899   PetscFunctionReturn(0);
900 }
901 
902 PetscBool SNEScite = PETSC_FALSE;
903 const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
904                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
905                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
906                             "  year = 2013,\n"
907                             "  type = Preprint,\n"
908                             "  number = {ANL/MCS-P2010-0112},\n"
909                             "  institution = {Argonne National Laboratory}\n}\n";
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "SNESSolve_FAS"
913 PetscErrorCode SNESSolve_FAS(SNES snes)
914 {
915   PetscErrorCode ierr;
916   PetscInt       i, maxits;
917   Vec            X, F;
918   PetscReal      fnorm;
919   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
920   DM             dm;
921   PetscBool      isFine;
922 
923   PetscFunctionBegin;
924   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
925   maxits       = snes->max_its;      /* maximum number of iterations */
926   snes->reason = SNES_CONVERGED_ITERATING;
927   X            = snes->vec_sol;
928   F            = snes->vec_func;
929 
930   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
931   /*norm setup */
932   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
933   snes->iter = 0;
934   snes->norm = 0.;
935   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
936   if (!snes->vec_func_init_set) {
937     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
938     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
939     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
940     if (snes->domainerror) {
941       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
942       PetscFunctionReturn(0);
943     }
944   } else snes->vec_func_init_set = PETSC_FALSE;
945 
946   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
947   if (PetscIsInfOrNanReal(fnorm)) {
948     snes->reason = SNES_DIVERGED_FNORM_NAN;
949     PetscFunctionReturn(0);
950   }
951 
952   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
953   snes->norm = fnorm;
954   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
955   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
956   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
957 
958   /* test convergence */
959   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
960   if (snes->reason) PetscFunctionReturn(0);
961 
962 
963   if (isFine) {
964     /* propagate scale-dependent data up the hierarchy */
965     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
966     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
967       DM dmcoarse;
968       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
969       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
970       dm   = dmcoarse;
971     }
972   }
973 
974   for (i = 0; i < maxits; i++) {
975     /* Call general purpose update function */
976 
977     if (snes->ops->update) {
978       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
979     }
980     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
981       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
982     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
983       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
984     } else if (fas->fastype == SNES_FAS_FULL) {
985       ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr);
986     } else if (fas->fastype ==SNES_FAS_KASKADE) {
987       ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr);
988     } else {
989       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");
990     }
991 
992     /* check for FAS cycle divergence */
993     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
994 
995     /* Monitor convergence */
996     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
997     snes->iter = i+1;
998     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
999     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
1000     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1001     /* Test for convergence */
1002     if (isFine) {
1003       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1004       if (snes->reason) break;
1005     }
1006   }
1007   if (i == maxits) {
1008     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1009     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1010   }
1011   PetscFunctionReturn(0);
1012 }
1013