xref: /petsc/src/snes/impls/fas/fas.c (revision 8b83055f287c8deb1a16325ae3faee24d0648b7a)
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         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
422         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
423         if (i && (smoothd == smoothu)) {
424           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
425         } else if (i) {
426           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
427           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
428           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
429           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
430         }
431       }
432     } else if (isdraw) {
433       PetscDraw draw;
434       PetscReal x,w,y,bottom,th,wth;
435       SNES_FAS  *curfas = fas;
436       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
437       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
438       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
439       bottom = y - th;
440       while (curfas) {
441         if (!curfas->smoothu) {
442           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
443           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
444           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
445         } else {
446           w    = 0.5*PetscMin(1.0-x,x);
447           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
448           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
449           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
450           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
451           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
452           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
453         }
454         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
455         bottom -= 5*th;
456         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
457         else curfas = NULL;
458       }
459     }
460   }
461   PetscFunctionReturn(0);
462 }
463 
464 #undef __FUNCT__
465 #define __FUNCT__ "SNESFASDownSmooth_Private"
466 /*
467 Defines the action of the downsmoother
468  */
469 PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
470 {
471   PetscErrorCode      ierr = 0;
472   SNESConvergedReason reason;
473   Vec                 FPC;
474   SNES                smoothd;
475   SNES_FAS            *fas = (SNES_FAS*) snes->data;
476 
477   PetscFunctionBegin;
478   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
479   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
480   ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr);
481   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
482   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
483   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
484   /* check convergence reason for the smoother */
485   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
486   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
487     snes->reason = SNES_DIVERGED_INNER;
488     PetscFunctionReturn(0);
489   }
490   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
491   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
492   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
493   PetscFunctionReturn(0);
494 }
495 
496 
497 #undef __FUNCT__
498 #define __FUNCT__ "SNESFASUpSmooth_Private"
499 /*
500 Defines the action of the upsmoother
501  */
502 PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
503 {
504   PetscErrorCode      ierr = 0;
505   SNESConvergedReason reason;
506   Vec                 FPC;
507   SNES                smoothu;
508   SNES_FAS            *fas = (SNES_FAS*) snes->data;
509 
510   PetscFunctionBegin;
511   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
512   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
513   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
514   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
515   /* check convergence reason for the smoother */
516   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
517   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
518     snes->reason = SNES_DIVERGED_INNER;
519     PetscFunctionReturn(0);
520   }
521   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
522   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
523   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "SNESFASCreateCoarseVec"
529 /*@
530    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
531 
532    Collective
533 
534    Input Arguments:
535 .  snes - SNESFAS
536 
537    Output Arguments:
538 .  Xcoarse - vector on level one coarser than snes
539 
540    Level: developer
541 
542 .seealso: SNESFASSetRestriction(), SNESFASRestrict()
543 @*/
544 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
545 {
546   PetscErrorCode ierr;
547   SNES_FAS       *fas = (SNES_FAS*)snes->data;
548 
549   PetscFunctionBegin;
550   if (fas->rscale) {
551     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
552   } else if (fas->inject) {
553     ierr = MatGetVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
554   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
555   PetscFunctionReturn(0);
556 }
557 
558 #undef __FUNCT__
559 #define __FUNCT__ "SNESFASRestrict"
560 /*@
561    SNESFASRestrict - restrict a Vec to the next coarser level
562 
563    Collective
564 
565    Input Arguments:
566 +  fine - SNES from which to restrict
567 -  Xfine - vector to restrict
568 
569    Output Arguments:
570 .  Xcoarse - result of restriction
571 
572    Level: developer
573 
574 .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
575 @*/
576 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
577 {
578   PetscErrorCode ierr;
579   SNES_FAS       *fas = (SNES_FAS*)fine->data;
580 
581   PetscFunctionBegin;
582   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
583   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
584   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
585   if (fas->inject) {
586     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
587   } else {
588     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
589     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
590   }
591   PetscFunctionReturn(0);
592 }
593 
594 #undef __FUNCT__
595 #define __FUNCT__ "SNESFASCoarseCorrection"
596 /*
597 
598 Performs the FAS coarse correction as:
599 
600 fine problem: F(x) = 0
601 coarse problem: F^c(x) = b^c
602 
603 b^c = F^c(I^c_fx^f - I^c_fF(x))
604 
605  */
606 PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
607 {
608   PetscErrorCode      ierr;
609   Vec                 X_c, Xo_c, F_c, B_c;
610   SNESConvergedReason reason;
611   SNES                next;
612   Mat                 restrct, interpolate;
613   SNES_FAS            *fasc;
614 
615   PetscFunctionBegin;
616   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
617   if (next) {
618     fasc = (SNES_FAS*)next->data;
619 
620     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
621     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
622 
623     X_c  = next->vec_sol;
624     Xo_c = next->work[0];
625     F_c  = next->vec_func;
626     B_c  = next->vec_rhs;
627 
628     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
629     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
630     /* restrict the defect */
631     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
632     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
633 
634     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
635     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
636     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
637 
638     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
639     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
640     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
641     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
642     /* set initial guess of the coarse problem to the projected fine solution */
643     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
644 
645     /* recurse to the next level */
646     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
647     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
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     /* correct as x <- x + I(x^c - Rx)*/
654     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
655 
656     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
657     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
658     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
659   }
660   PetscFunctionReturn(0);
661 }
662 
663 #undef __FUNCT__
664 #define __FUNCT__ "SNESFASCycle_Additive"
665 /*
666 
667 The additive cycle looks like:
668 
669 xhat = x
670 xhat = dS(x, b)
671 x = coarsecorrection(xhat, b_d)
672 x = x + nu*(xhat - x);
673 (optional) x = uS(x, b)
674 
675 With the coarse RHS (defect correction) as below.
676 
677  */
678 PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
679 {
680   Vec                 F, B, Xhat;
681   Vec                 X_c, Xo_c, F_c, B_c;
682   PetscErrorCode      ierr;
683   SNESConvergedReason reason;
684   PetscReal           xnorm, fnorm, ynorm;
685   PetscBool           lssuccess;
686   SNES                next;
687   Mat                 restrct, interpolate;
688   SNES_FAS            *fas = (SNES_FAS*)snes->data,*fasc;
689 
690   PetscFunctionBegin;
691   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
692   F    = snes->vec_func;
693   B    = snes->vec_rhs;
694   Xhat = snes->work[1];
695   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
696   /* recurse first */
697   if (next) {
698     fasc = (SNES_FAS*)next->data;
699     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
700     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
701     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
702     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
703     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
704     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
705     X_c  = next->vec_sol;
706     Xo_c = next->work[0];
707     F_c  = next->vec_func;
708     B_c  = next->vec_rhs;
709 
710     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
711     /* restrict the defect */
712     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
713 
714     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
715     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
716     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
717     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
718     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
719     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
720     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
721     /* set initial guess of the coarse problem to the projected fine solution */
722     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
723 
724     /* recurse */
725     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
726     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
727 
728     /* smooth on this level */
729     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
730 
731     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
732     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
733       snes->reason = SNES_DIVERGED_INNER;
734       PetscFunctionReturn(0);
735     }
736 
737     /* correct as x <- x + I(x^c - Rx)*/
738     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
739     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
740 
741     /* additive correction of the coarse direction*/
742     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
743     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
744     if (!lssuccess) {
745       if (++snes->numFailures >= snes->maxFailures) {
746         snes->reason = SNES_DIVERGED_LINE_SEARCH;
747         PetscFunctionReturn(0);
748       }
749     }
750     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
751   } else {
752     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
753   }
754   PetscFunctionReturn(0);
755 }
756 
757 #undef __FUNCT__
758 #define __FUNCT__ "SNESFASCycle_Multiplicative"
759 /*
760 
761 Defines the FAS cycle as:
762 
763 fine problem: F(x) = 0
764 coarse problem: F^c(x) = b^c
765 
766 b^c = F^c(I^c_fx^f - I^c_fF(x))
767 
768 correction:
769 
770 x = x + I(x^c - Rx)
771 
772  */
773 PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
774 {
775 
776   PetscErrorCode ierr;
777   Vec            F,B;
778   SNES_FAS       *fas = (SNES_FAS*)snes->data;
779 
780   PetscFunctionBegin;
781   F = snes->vec_func;
782   B = snes->vec_rhs;
783   /* pre-smooth -- just update using the pre-smoother */
784   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
785   if (fas->level != 0) {
786     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
787     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
788   }
789   PetscFunctionReturn(0);
790 }
791 
792 #undef __FUNCT__
793 #define __FUNCT__ "SNESSolve_FAS"
794 
795 PetscErrorCode SNESSolve_FAS(SNES snes)
796 {
797   PetscErrorCode ierr;
798   PetscInt       i, maxits;
799   Vec            X, F;
800   PetscReal      fnorm;
801   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
802   DM             dm;
803   PetscBool      isFine;
804 
805   PetscFunctionBegin;
806   maxits       = snes->max_its;      /* maximum number of iterations */
807   snes->reason = SNES_CONVERGED_ITERATING;
808   X            = snes->vec_sol;
809   F            = snes->vec_func;
810 
811   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
812   /*norm setup */
813   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
814   snes->iter = 0;
815   snes->norm = 0.;
816   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
817   if (!snes->vec_func_init_set) {
818     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
819     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
820     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
821     if (snes->domainerror) {
822       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
823       PetscFunctionReturn(0);
824     }
825   } else snes->vec_func_init_set = PETSC_FALSE;
826 
827   if (!snes->norm_init_set) {
828     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
829     if (PetscIsInfOrNanReal(fnorm)) {
830       snes->reason = SNES_DIVERGED_FNORM_NAN;
831       PetscFunctionReturn(0);
832     }
833   } else {
834     fnorm               = snes->norm_init;
835     snes->norm_init_set = PETSC_FALSE;
836   }
837 
838   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
839   snes->norm = fnorm;
840   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
841   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
842   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
843 
844   /* set parameter for default relative tolerance convergence test */
845   snes->ttol = fnorm*snes->rtol;
846   /* test convergence */
847   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
848   if (snes->reason) PetscFunctionReturn(0);
849 
850 
851   if (isFine) {
852     /* propagate scale-dependent data up the hierarchy */
853     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
854     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
855       DM dmcoarse;
856       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
857       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
858       dm   = dmcoarse;
859     }
860   }
861 
862   for (i = 0; i < maxits; i++) {
863     /* Call general purpose update function */
864 
865     if (snes->ops->update) {
866       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
867     }
868     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
869       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
870     } else {
871       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
872     }
873 
874     /* check for FAS cycle divergence */
875     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
876 
877     /* Monitor convergence */
878     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
879     snes->iter = i+1;
880     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
881     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
882     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
883     /* Test for convergence */
884     if (isFine) {
885       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
886       if (snes->reason) break;
887     }
888   }
889   if (i == maxits) {
890     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
891     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
892   }
893   PetscFunctionReturn(0);
894 }
895