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