xref: /petsc/src/snes/impls/fas/fas.c (revision ed020824cf2871ebf0db89bae9162301a51fd7bc)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>
3 
4 /*MC
5 Full Approximation Scheme nonlinear multigrid solver.
6 
7 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
8 
9 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
10 M*/
11 
12 extern PetscErrorCode SNESDestroy_FAS(SNES snes);
13 extern PetscErrorCode SNESSetUp_FAS(SNES snes);
14 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
15 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
16 extern PetscErrorCode SNESSolve_FAS(SNES snes);
17 extern PetscErrorCode SNESReset_FAS(SNES snes);
18 
19 EXTERN_C_BEGIN
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "SNESCreate_FAS"
23 PetscErrorCode SNESCreate_FAS(SNES snes)
24 {
25   SNES_FAS * fas;
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   snes->ops->destroy        = SNESDestroy_FAS;
30   snes->ops->setup          = SNESSetUp_FAS;
31   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
32   snes->ops->view           = SNESView_FAS;
33   snes->ops->solve          = SNESSolve_FAS;
34   snes->ops->reset          = SNESReset_FAS;
35 
36   snes->usesksp             = PETSC_FALSE;
37   snes->usespc              = PETSC_FALSE;
38 
39   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
40   snes->data                = (void*) fas;
41   fas->level                = 0;
42   fas->levels               = 1;
43   fas->n_cycles             = 1;
44   fas->max_up_it            = 1;
45   fas->max_down_it          = 1;
46   fas->upsmooth             = PETSC_NULL;
47   fas->downsmooth           = PETSC_NULL;
48   fas->next                 = PETSC_NULL;
49   fas->interpolate          = PETSC_NULL;
50   fas->restrct              = PETSC_NULL;
51 
52   PetscFunctionReturn(0);
53 }
54 EXTERN_C_END
55 
56 #undef __FUNCT__
57 #define __FUNCT__ "SNESFASGetLevels"
58 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
59   SNES_FAS * fas = (SNES_FAS *)snes->data;
60   PetscFunctionBegin;
61   *levels = fas->levels;
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "SNESFASGetSNES"
67 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
68   SNES_FAS * fas = (SNES_FAS *)snes->data;
69   PetscInt levels = fas->level;
70   PetscInt i;
71   PetscFunctionBegin;
72   *lsnes = snes;
73   if (fas->level < level) {
74     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
75   }
76   if (level > levels - 1) {
77     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
78   }
79   for (i = fas->level; i > level; i--) {
80     *lsnes = fas->next;
81     fas = (SNES_FAS *)(*lsnes)->data;
82   }
83   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
84   PetscFunctionReturn(0);
85 }
86 
87 #undef __FUNCT__
88 #define __FUNCT__ "SNESFASSetLevels"
89 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
90   PetscErrorCode ierr;
91   PetscInt i;
92   SNES_FAS * fas = (SNES_FAS *)snes->data;
93   MPI_Comm comm;
94   PetscFunctionBegin;
95   comm = ((PetscObject)snes)->comm;
96   if (levels == fas->levels) {
97     if (!comms)
98       PetscFunctionReturn(0);
99   }
100   /* user has changed the number of levels; reset */
101   ierr = SNESReset(snes);CHKERRQ(ierr);
102   /* destroy any coarser levels if necessary */
103   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
104   fas->next = PETSC_NULL;
105   /* setup the finest level */
106   for (i = levels - 1; i >= 0; i--) {
107     if (comms) comm = comms[i];
108     fas->level = i;
109     fas->levels = levels;
110     fas->next = PETSC_NULL;
111     if (i > 0) {
112       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
113       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
114       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
115       fas = (SNES_FAS *)fas->next->data;
116     }
117   }
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "SNESFASSetInterpolation"
123 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
124   SNES_FAS * fas =  (SNES_FAS *)snes->data;
125   PetscInt top_level = fas->level,i;
126 
127   PetscFunctionBegin;
128   if (level > top_level)
129     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
130   /* get to the correct level */
131   for (i = fas->level; i > level; i--) {
132     fas = (SNES_FAS *)fas->next->data;
133   }
134   if (fas->level != level)
135     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
136   fas->interpolate = mat;
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "SNESFASSetRestriction"
142 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
143   SNES_FAS * fas =  (SNES_FAS *)snes->data;
144   PetscInt top_level = fas->level,i;
145 
146   PetscFunctionBegin;
147   if (level > top_level)
148     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
149   /* get to the correct level */
150   for (i = fas->level; i > level; i--) {
151     fas = (SNES_FAS *)fas->next->data;
152   }
153   if (fas->level != level)
154     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
155   fas->restrct = mat;
156   PetscFunctionReturn(0);
157 }
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "SNESFASSetRScale"
161 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
162   SNES_FAS * fas =  (SNES_FAS *)snes->data;
163   PetscInt top_level = fas->level,i;
164 
165   PetscFunctionBegin;
166   if (level > top_level)
167     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
168   /* get to the correct level */
169   for (i = fas->level; i > level; i--) {
170     fas = (SNES_FAS *)fas->next->data;
171   }
172   if (fas->level != level)
173     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
174   fas->rscale = rscale;
175   PetscFunctionReturn(0);
176 }
177 
178 #undef __FUNCT__
179 #define __FUNCT__ "SNESReset_FAS"
180 PetscErrorCode SNESReset_FAS(SNES snes)
181 {
182   PetscErrorCode ierr = 0;
183   SNES_FAS * fas = (SNES_FAS *)snes->data;
184 
185   PetscFunctionBegin;
186   /* destroy local data created in SNESSetup_FAS */
187   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
188   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
189   if (fas->interpolate == fas->restrct) {
190     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
191     fas->restrct = PETSC_NULL;
192   } else {
193     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
194     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
195   }
196   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
197   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
198   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
199   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "SNESDestroy_FAS"
205 PetscErrorCode SNESDestroy_FAS(SNES snes)
206 {
207   SNES_FAS * fas = (SNES_FAS *)snes->data;
208   PetscErrorCode ierr;
209 
210   PetscFunctionBegin;
211   /* recursively resets and then destroys */
212   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
213   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
214   ierr = PetscFree(fas);CHKERRQ(ierr);
215 
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "SNESSetUp_FAS"
221 PetscErrorCode SNESSetUp_FAS(SNES snes)
222 {
223   SNES_FAS   *fas = (SNES_FAS *) snes->data;
224   PetscErrorCode ierr;
225 
226   PetscFunctionBegin;
227   ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
228   /* gets the solver ready for solution */
229   if (snes->dm) {
230     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
231     if (fas->next) {
232       /* for now -- assume the DM and the evaluation functions have been set externally */
233         if (!fas->next->dm) {
234           ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
235           ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
236         }
237         /* set the interpolation and restriction from the DM */
238         if (!fas->interpolate) {
239           ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
240           fas->restrct = fas->interpolate;
241         }
242     }
243     /* set the DMs of the pre and post-smoothers here */
244     if (fas->upsmooth)  SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);
245     if (fas->downsmooth)SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);
246   }
247   if (fas->next) {
248     /* gotta set up the solution vector for this to work */
249     ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);
250     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
251   }
252   /* got to set them all up at once */
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "SNESSetFromOptions_FAS"
258 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
259 {
260   SNES_FAS   *fas = (SNES_FAS *) snes->data;
261   PetscInt levels = 1;
262   PetscBool flg, monflg;
263   PetscErrorCode ierr;
264   const char * def_smooth = SNESNRICHARDSON;
265   char pre_type[256];
266   char post_type[256];
267   char                    monfilename[PETSC_MAX_PATH_LEN];
268 
269   PetscFunctionBegin;
270   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
271 
272   /* number of levels -- only process on the finest level */
273   if (fas->levels - 1 == fas->level) {
274     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
275     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
276   }
277 
278   /* type of pre and/or post smoothers -- set both at once */
279   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
280   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
281   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
282   if (flg) {
283     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
284   } else {
285     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
286     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
287   }
288 
289   /* options for the number of preconditioning cycles and cycle type */
290   ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
291   ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
292   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
293 
294   ierr = PetscOptionsString("-snes_fas_monitor","Monitor for smoothers","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
295 
296   /* other options for the coarsest level */
297   if (fas->level == 0) {
298     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
299     if (flg) {
300       ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
301     } else {
302       ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
303       ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
304     }
305   }
306 
307   ierr = PetscOptionsTail();CHKERRQ(ierr);
308   /* setup from the determined types if the smoothers don't exist */
309   if (!fas->upsmooth) {
310     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
311     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
312     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
313   }
314   if (fas->upsmooth) ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
315   if (!fas->downsmooth && fas->level != 0) {
316     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
317     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
318    ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
319   }
320   if (fas->downsmooth) ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
321 
322   if (monflg) {
323     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
324     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
325   }
326 
327   /* recursive option setting for the smoothers */
328   if (fas->next)ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);
329   PetscFunctionReturn(0);
330 }
331 
332 #undef __FUNCT__
333 #define __FUNCT__ "SNESView_FAS"
334 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
335 {
336   SNES_FAS   *fas = (SNES_FAS *) snes->data;
337   PetscBool      iascii;
338   PetscErrorCode ierr;
339   PetscInt levels = fas->levels;
340   PetscInt i;
341 
342   PetscFunctionBegin;
343   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
344   if (iascii) {
345     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
346     ierr = PetscViewerASCIIPushTab(viewer);
347     for (i = levels - 1; i >= 0; i--) {
348       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
349       if (fas->upsmooth) {
350         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
351         ierr = PetscViewerASCIIPushTab(viewer);
352         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
353         ierr = PetscViewerASCIIPopTab(viewer);
354       } else {
355         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
356       }
357       if (fas->downsmooth) {
358         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
359         ierr = PetscViewerASCIIPushTab(viewer);
360         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
361         ierr = PetscViewerASCIIPopTab(viewer);
362       } else {
363         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
364       }
365       if (fas->next) fas = (SNES_FAS *)fas->next->data;
366     }
367     ierr = PetscViewerASCIIPopTab(viewer);
368   } else {
369     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
370   }
371   PetscFunctionReturn(0);
372 }
373 
374 /*
375 
376 Defines the FAS cycle as:
377 
378 fine problem: F(x) = 0
379 coarse problem: F^c(x) = b^c
380 
381 b^c = F^c(I^c_fx^f - I^c_fF(x))
382 
383 correction:
384 
385 x = x + I(x^c - Rx)
386 
387  */
388 
389 #undef __FUNCT__
390 #define __FUNCT__ "FASCycle_Private"
391 PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X, Vec F) {
392 
393   PetscErrorCode ierr;
394   Vec X_c, Xo_c, F_c, B_c;
395   SNES_FAS * fas = (SNES_FAS *)snes->data;
396   PetscInt i;
397 
398   PetscFunctionBegin;
399   /* pre-smooth -- just update using the pre-smoother */
400   if (fas->upsmooth) {
401     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
402   } else if (snes->pc) {
403     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
404   }
405   if (fas->next) {
406     for (i = 0; i < fas->n_cycles; i++) {
407       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
408       X_c  = fas->next->vec_sol;
409       Xo_c = fas->next->work[0];
410       F_c  = fas->next->vec_func;
411       B_c  = fas->next->work[1];
412       /* inject the solution to coarse */
413       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
414       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
415       if (B) {
416         ierr = VecAYPX(F, -1.0, B);CHKERRQ(ierr); /* F = B - F (defect) */
417       } else {
418         ierr = VecScale(F, -1.0);CHKERRQ(ierr);
419       }
420 
421       /* restrict the defect */
422       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
423       ierr = VecPointwiseMult(B_c,  fas->rscale, B_c);CHKERRQ(ierr);
424 
425       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
426       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
427       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
428       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
429 
430       /* set initial guess of the coarse problem to the projected fine solution */
431       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
432 
433       /* recurse to the next level */
434       ierr = FASCycle_Private(fas->next, B_c, X_c, F_c);CHKERRQ(ierr);
435 
436       /* correct as x <- x + I(x^c - Rx)*/
437       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
438       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
439     }
440   }
441     /* down-smooth -- just update using the down-smoother */
442   if (fas->level != 0) {
443     if (fas->downsmooth) {
444       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
445     } else if (snes->pc) {
446       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
447     }
448   }
449   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
450 
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "FASInitialGuess_Private"
456 PetscErrorCode FASInitialGuess_Private(SNES snes, Vec B, Vec X) {
457 
458   PetscErrorCode ierr;
459   Vec X_c, B_c;
460   SNES_FAS * fas;
461 
462   PetscFunctionBegin;
463   fas = (SNES_FAS *)snes->data;
464   /* pre-smooth -- just update using the pre-smoother */
465   if (fas->level == 0) {
466     if (fas->upsmooth) {
467       ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
468     } else if (snes->pc) {
469       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
470     }
471   }
472   if (fas->next) {
473     X_c  = fas->next->vec_sol;
474     B_c  = fas->next->work[0];
475     /* inject the solution to coarse */
476     ierr = MatRestrict(fas->restrct, X, X_c);CHKERRQ(ierr);
477     ierr = VecPointwiseMult(X_c, fas->rscale, X_c);CHKERRQ(ierr);
478     if (B) {
479       ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr);
480       ierr = VecPointwiseMult(B_c, fas->rscale, B_c);CHKERRQ(ierr);
481     } else {
482       B_c = PETSC_NULL;
483     }
484     /* recurse to the next level */
485     ierr = FASInitialGuess_Private(fas->next, B_c, X_c);CHKERRQ(ierr);
486     ierr = MatInterpolate(fas->interpolate, X_c, X);CHKERRQ(ierr);
487   }
488   /* down-smooth -- just update using the down-smoother */
489   if (fas->level != 0) {
490     if (fas->downsmooth) {
491       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
492     } else if (snes->pc) {
493       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
494     }
495   }
496   PetscFunctionReturn(0);
497 }
498 
499 #undef __FUNCT__
500 #define __FUNCT__ "SNESSolve_FAS"
501 
502 PetscErrorCode SNESSolve_FAS(SNES snes)
503 {
504   PetscErrorCode ierr;
505   PetscInt i, maxits;
506   Vec X, F, B;
507   PetscReal fnorm;
508   PetscFunctionBegin;
509   maxits = snes->max_its;            /* maximum number of iterations */
510   snes->reason = SNES_CONVERGED_ITERATING;
511   X = snes->vec_sol;
512   F = snes->vec_func;
513   B = snes->vec_rhs;
514 
515   /*norm setup */
516   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
517   snes->iter = 0;
518   snes->norm = 0.;
519   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
520   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
521   if (snes->domainerror) {
522     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
523     PetscFunctionReturn(0);
524   }
525   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
526   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
527   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
528   snes->norm = fnorm;
529   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
530   SNESLogConvHistory(snes,fnorm,0);
531   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
532 
533   /* set parameter for default relative tolerance convergence test */
534   snes->ttol = fnorm*snes->rtol;
535   /* test convergence */
536   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
537   if (snes->reason) PetscFunctionReturn(0);
538   for (i = 0; i < maxits; i++) {
539     /* Call general purpose update function */
540     if (snes->ops->update) {
541       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
542     }
543     ierr = FASCycle_Private(snes, B, X, F);CHKERRQ(ierr);
544     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
545     /* Monitor convergence */
546     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
547     snes->iter = i+1;
548     snes->norm = fnorm;
549     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
550     SNESLogConvHistory(snes,snes->norm,0);
551     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
552     /* Test for convergence */
553     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
554     if (snes->reason) break;
555   }
556   if (i == maxits) {
557     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
558     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
559   }
560   PetscFunctionReturn(0);
561 }
562