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