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