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