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