xref: /petsc/src/snes/impls/fas/fas.c (revision eff52c0e83802002f81775ce8e94c1333cc52f8e)
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   fas->inject                 = PETSC_NULL;
52   fas->monitor                = PETSC_NULL;
53   fas->usedmfornumberoflevels = PETSC_FALSE;
54   fas->useGS                  = PETSC_FALSE;
55   PetscFunctionReturn(0);
56 }
57 EXTERN_C_END
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "SNESFASGetLevels"
61 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
62   SNES_FAS * fas = (SNES_FAS *)snes->data;
63   PetscFunctionBegin;
64   *levels = fas->levels;
65   PetscFunctionReturn(0);
66 }
67 
68 #undef __FUNCT__
69 #define __FUNCT__ "SNESFASSetCycles"
70 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
71   SNES_FAS * fas = (SNES_FAS *)snes->data;
72   PetscErrorCode ierr;
73   PetscFunctionBegin;
74   fas->n_cycles = cycles;
75   if (fas->next) {
76     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
77   }
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "SNESFASSetGS"
83 PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
84   PetscErrorCode ierr = 0;
85   SNES_FAS       *fas = (SNES_FAS *)snes->data;
86   PetscFunctionBegin;
87 
88   /* use or don't use it according to user wishes*/
89   fas->useGS = use_gs;
90   if (gsfunc) {
91     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
92     /* push the provided GS up the tree */
93     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
94   } else if (snes->ops->computegs) {
95     /* assume that the user has set the GS solver at this level */
96     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
97   } else if (use_gs) {
98     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level);
99   }
100   PetscFunctionReturn(0);
101 }
102 
103 #undef __FUNCT__
104 #define __FUNCT__ "SNESFASSetGSOnLevel"
105 PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
106   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
107   PetscErrorCode ierr;
108   PetscInt       top_level = fas->level,i;
109   SNES           cur_snes = snes;
110   PetscFunctionBegin;
111   if (level > top_level)
112     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
113   /* get to the correct level */
114   for (i = fas->level; i > level; i--) {
115     fas = (SNES_FAS *)fas->next->data;
116     cur_snes = fas->next;
117   }
118   if (fas->level != level)
119     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
120   fas->useGS = use_gs;
121   if (gsfunc) {
122     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
123   }
124   PetscFunctionReturn(0);
125 }
126 
127 
128 #undef __FUNCT__
129 #define __FUNCT__ "SNESFASSetCyclesOnLevel"
130 PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
131   SNES_FAS * fas =  (SNES_FAS *)snes->data;
132   PetscInt top_level = fas->level,i;
133 
134   PetscFunctionBegin;
135   if (level > top_level)
136     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
137   /* get to the correct level */
138   for (i = fas->level; i > level; i--) {
139     fas = (SNES_FAS *)fas->next->data;
140   }
141   if (fas->level != level)
142     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
143   fas->n_cycles = cycles;
144   PetscFunctionReturn(0);
145 }
146 
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "SNESFASGetSNES"
150 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
151   SNES_FAS * fas = (SNES_FAS *)snes->data;
152   PetscInt levels = fas->level;
153   PetscInt i;
154   PetscFunctionBegin;
155   *lsnes = snes;
156   if (fas->level < level) {
157     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
158   }
159   if (level > levels - 1) {
160     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
161   }
162   for (i = fas->level; i > level; i--) {
163     *lsnes = fas->next;
164     fas = (SNES_FAS *)(*lsnes)->data;
165   }
166   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
167   PetscFunctionReturn(0);
168 }
169 
170 #undef __FUNCT__
171 #define __FUNCT__ "SNESFASSetLevels"
172 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
173   PetscErrorCode ierr;
174   PetscInt i;
175   SNES_FAS * fas = (SNES_FAS *)snes->data;
176   MPI_Comm comm;
177   PetscFunctionBegin;
178    comm = ((PetscObject)snes)->comm;
179   if (levels == fas->levels) {
180     if (!comms)
181       PetscFunctionReturn(0);
182   }
183   /* user has changed the number of levels; reset */
184   ierr = SNESReset(snes);CHKERRQ(ierr);
185   /* destroy any coarser levels if necessary */
186   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
187   fas->next = PETSC_NULL;
188   /* setup the finest level */
189   for (i = levels - 1; i >= 0; i--) {
190     if (comms) comm = comms[i];
191     fas->level = i;
192     fas->levels = levels;
193     fas->next = PETSC_NULL;
194     if (i > 0) {
195       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
196       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
197       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
198       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
199       fas = (SNES_FAS *)fas->next->data;
200     }
201   }
202   /* by default set the GS smoothers up the chain if one exists on the finest level */
203   if (snes->ops->computegs)
204     ierr = SNESFASSetGS(snes, snes->ops->computegs, snes->gsP, fas->useGS);CHKERRQ(ierr);
205   PetscFunctionReturn(0);
206 }
207 
208 #undef __FUNCT__
209 #define __FUNCT__ "SNESFASSetInterpolation"
210 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
211   SNES_FAS * fas =  (SNES_FAS *)snes->data;
212   PetscInt top_level = fas->level,i;
213 
214   PetscFunctionBegin;
215   if (level > top_level)
216     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
217   /* get to the correct level */
218   for (i = fas->level; i > level; i--) {
219     fas = (SNES_FAS *)fas->next->data;
220   }
221   if (fas->level != level)
222     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
223   fas->interpolate = mat;
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "SNESFASSetRestriction"
229 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
230   SNES_FAS * fas =  (SNES_FAS *)snes->data;
231   PetscInt top_level = fas->level,i;
232 
233   PetscFunctionBegin;
234   if (level > top_level)
235     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
236   /* get to the correct level */
237   for (i = fas->level; i > level; i--) {
238     fas = (SNES_FAS *)fas->next->data;
239   }
240   if (fas->level != level)
241     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
242   fas->restrct = mat;
243   PetscFunctionReturn(0);
244 }
245 
246 
247 #undef __FUNCT__
248 #define __FUNCT__ "SNESFASSetRScale"
249 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
250   SNES_FAS * fas =  (SNES_FAS *)snes->data;
251   PetscInt top_level = fas->level,i;
252 
253   PetscFunctionBegin;
254   if (level > top_level)
255     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
256   /* get to the correct level */
257   for (i = fas->level; i > level; i--) {
258     fas = (SNES_FAS *)fas->next->data;
259   }
260   if (fas->level != level)
261     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
262   fas->rscale = rscale;
263   PetscFunctionReturn(0);
264 }
265 
266 
267 #undef __FUNCT__
268 #define __FUNCT__ "SNESFASSetInjection"
269 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
270   SNES_FAS * fas =  (SNES_FAS *)snes->data;
271   PetscInt top_level = fas->level,i;
272 
273   PetscFunctionBegin;
274   if (level > top_level)
275     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
276   /* get to the correct level */
277   for (i = fas->level; i > level; i--) {
278     fas = (SNES_FAS *)fas->next->data;
279   }
280   if (fas->level != level)
281     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
282   fas->inject = mat;
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "SNESReset_FAS"
288 PetscErrorCode SNESReset_FAS(SNES snes)
289 {
290   PetscErrorCode ierr = 0;
291   SNES_FAS * fas = (SNES_FAS *)snes->data;
292 
293   PetscFunctionBegin;
294   /* destroy local data created in SNESSetup_FAS */
295 #if 0
296   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
297 #endif
298   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
299 #if 0
300 #endif
301   PetscFunctionReturn(0);
302 }
303 
304 #undef __FUNCT__
305 #define __FUNCT__ "SNESDestroy_FAS"
306 PetscErrorCode SNESDestroy_FAS(SNES snes)
307 {
308   SNES_FAS * fas = (SNES_FAS *)snes->data;
309   PetscErrorCode ierr;
310 
311   PetscFunctionBegin;
312   /* recursively resets and then destroys */
313   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
314   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
315   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
316   if (fas->inject) {
317     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
318   }
319   if (fas->interpolate == fas->restrct) {
320     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
321     fas->restrct = PETSC_NULL;
322   } else {
323     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
324     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
325   }
326   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
327   if (snes->work)        ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);
328   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
329   ierr = PetscFree(fas);CHKERRQ(ierr);
330 
331   PetscFunctionReturn(0);
332 }
333 
334 #undef __FUNCT__
335 #define __FUNCT__ "SNESSetUp_FAS"
336 PetscErrorCode SNESSetUp_FAS(SNES snes)
337 {
338   SNES_FAS       *fas = (SNES_FAS *) snes->data,*tmp;
339   PetscErrorCode ierr;
340   VecScatter     injscatter;
341   PetscInt       dm_levels;
342 
343 
344   PetscFunctionBegin;
345 
346   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
347     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
348     dm_levels++;
349     if (dm_levels > fas->levels) {
350       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
351       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
352     }
353   }
354 
355   if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */
356 
357   /* should call the SNESSetFromOptions() only when appropriate */
358   tmp = fas;
359   while (tmp) {
360     if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);}
361     if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);}
362     tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0;
363   }
364 
365   /* gets the solver ready for solution */
366   if (fas->next) {
367     if (snes->ops->computegs) {
368       ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
369     }
370   }
371   if (snes->dm) {
372     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
373     if (fas->next) {
374       /* for now -- assume the DM and the evaluation functions have been set externally */
375       if (!fas->next->dm) {
376         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
377         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
378       }
379       /* set the interpolation and restriction from the DM */
380       if (!fas->interpolate) {
381         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
382         fas->restrct = fas->interpolate;
383       }
384       /* set the injection from the DM */
385       if (!fas->inject) {
386         ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
387         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
388         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
389       }
390     }
391     /* set the DMs of the pre and post-smoothers here */
392     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
393     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
394   }
395 
396   if (fas->next) {
397    /* gotta set up the solution vector for this to work */
398     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
399     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
400   }
401   /* got to set them all up at once */
402   PetscFunctionReturn(0);
403 }
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "SNESSetFromOptions_FAS"
407 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
408 {
409   SNES_FAS   *fas = (SNES_FAS *) snes->data;
410   PetscInt levels = 1;
411   PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg;
412   PetscErrorCode ierr;
413   const char * def_smooth = SNESNRICHARDSON;
414   char pre_type[256];
415   char post_type[256];
416   char                    monfilename[PETSC_MAX_PATH_LEN];
417 
418   PetscFunctionBegin;
419   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
420 
421   /* number of levels -- only process on the finest level */
422   if (fas->levels - 1 == fas->level) {
423     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
424     if (!flg && snes->dm) {
425       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
426       levels++;
427       fas->usedmfornumberoflevels = PETSC_TRUE;
428     }
429     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
430   }
431 
432   /* type of pre and/or post smoothers -- set both at once */
433   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
434   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
435   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
436   if (smoothflg) {
437     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
438   } else {
439     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr);
440     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr);
441   }
442 
443   /* options for the number of preconditioning cycles and cycle type */
444   ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
445   ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
446   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
447 
448   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
449 
450   ierr = PetscOptionsBool("-snes_fas_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr);
451 
452   /* other options for the coarsest level */
453   if (fas->level == 0) {
454     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
455     ierr = PetscOptionsBool("-snes_fas_coarse_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr);
456   }
457 
458   ierr = PetscOptionsTail();CHKERRQ(ierr);
459   /* setup from the determined types if there is no pointwise procedure or smoother defined */
460 
461   if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !fas->useGS)) {
462     const char     *prefix;
463     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
464     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
465     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
466     if (fas->level || (fas->levels == 1)) {
467       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr);
468     } else {
469       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
470     }
471     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
472     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
473     if (snes->ops->computefunction) {
474       ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
475     }
476   }
477 
478   if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !fas->useGS)) {
479     const char     *prefix;
480     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
481     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
482     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
483     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr);
484     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
485     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
486     if (snes->ops->computefunction) {
487       ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
488     }
489   }
490   if (fas->upsmooth) {
491     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
492   }
493 
494   if (fas->downsmooth) {
495     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
496   }
497 
498   if (monflg) {
499     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
500     /* set the monitors for the upsmoother and downsmoother */
501     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
502     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
503   }
504 
505   /* recursive option setting for the smoothers */
506   if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);}
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "SNESView_FAS"
512 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
513 {
514   SNES_FAS   *fas = (SNES_FAS *) snes->data;
515   PetscBool      iascii;
516   PetscErrorCode ierr;
517   PetscInt levels = fas->levels;
518   PetscInt i;
519 
520   PetscFunctionBegin;
521   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
522   if (iascii) {
523     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
524     ierr = PetscViewerASCIIPushTab(viewer);
525     for (i = levels - 1; i >= 0; i--) {
526       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
527       if (fas->upsmooth) {
528         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
529         ierr = PetscViewerASCIIPushTab(viewer);
530         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
531         ierr = PetscViewerASCIIPopTab(viewer);
532       } else {
533         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
534       }
535       if (fas->downsmooth) {
536         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
537         ierr = PetscViewerASCIIPushTab(viewer);
538         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
539         ierr = PetscViewerASCIIPopTab(viewer);
540       } else {
541         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
542       }
543       if (fas->useGS) {
544         ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n",  fas->level);CHKERRQ(ierr);
545       }
546       if (fas->next) fas = (SNES_FAS *)fas->next->data;
547     }
548     ierr = PetscViewerASCIIPopTab(viewer);
549   } else {
550     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
551   }
552   PetscFunctionReturn(0);
553 }
554 
555 /*
556 
557 Defines the FAS cycle as:
558 
559 fine problem: F(x) = 0
560 coarse problem: F^c(x) = b^c
561 
562 b^c = F^c(I^c_fx^f - I^c_fF(x))
563 
564 correction:
565 
566 x = x + I(x^c - Rx)
567 
568  */
569 
570 #undef __FUNCT__
571 #define __FUNCT__ "FASCycle_Private"
572 PetscErrorCode FASCycle_Private(SNES snes, Vec X) {
573 
574   PetscErrorCode ierr;
575   Vec X_c, Xo_c, F_c, B_c,F,B;
576   SNES_FAS * fas = (SNES_FAS *)snes->data;
577   PetscInt i, k;
578   PetscReal fnorm;
579 
580   PetscFunctionBegin;
581   F = snes->vec_func;
582   B = snes->vec_rhs;
583   /* pre-smooth -- just update using the pre-smoother */
584   if (fas->downsmooth) {
585     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
586   } else if (snes->ops->computegs) {
587     if (fas->monitor) {
588       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
589       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
590       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
591       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
592       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
593     }
594     for (k = 0; k < fas->max_up_it; k++) {
595       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
596       if (fas->monitor) {
597         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
598         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
599         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
600         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
601         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
602       }
603     }
604   } else if (snes->pc) {
605     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
606   }
607   if (fas->next) {
608     for (i = 0; i < fas->n_cycles; i++) {
609       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
610 
611       X_c  = fas->next->vec_sol;
612       Xo_c = fas->next->work[0];
613       F_c  = fas->next->vec_func;
614       B_c  = fas->next->work[1];
615 
616       /* inject the solution */
617       if (fas->inject) {
618         ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
619       } else {
620         ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
621         ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
622       }
623       ierr = VecScale(F, -1.0);CHKERRQ(ierr);
624 
625       /* restrict the defect */
626       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
627 
628       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
629       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
630       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
631       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
632 
633       /* set initial guess of the coarse problem to the projected fine solution */
634       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
635 
636       /* recurse to the next level */
637       fas->next->vec_rhs = B_c;
638       ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr);
639       fas->next->vec_rhs = PETSC_NULL;
640 
641       /* correct as x <- x + I(x^c - Rx)*/
642       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
643       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
644     }
645   }
646     /* up-smooth -- just update using the down-smoother */
647   if (fas->level != 0) {
648     if (fas->upsmooth) {
649       ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
650     } else if (snes->ops->computegs) {
651       if (fas->monitor) {
652         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
653         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
654         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
655         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
656         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
657       }
658       for (k = 0; k < fas->max_down_it; k++) {
659         ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
660         if (fas->monitor) {
661           ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
662           ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
663           ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
664           ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
665           ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
666         }
667       }
668     } else if (snes->pc) {
669       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
670     }
671   }
672   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
673 
674   PetscFunctionReturn(0);
675 }
676 
677 #undef __FUNCT__
678 #define __FUNCT__ "SNESSolve_FAS"
679 
680 PetscErrorCode SNESSolve_FAS(SNES snes)
681 {
682   PetscErrorCode ierr;
683   PetscInt i, maxits;
684   Vec X, B,F;
685   PetscReal fnorm;
686   PetscFunctionBegin;
687   maxits = snes->max_its;            /* maximum number of iterations */
688   snes->reason = SNES_CONVERGED_ITERATING;
689   X = snes->vec_sol;
690   B = snes->vec_rhs;
691   F = snes->vec_func;
692 
693   /*norm setup */
694   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
695   snes->iter = 0;
696   snes->norm = 0.;
697   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
698   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
699   if (snes->domainerror) {
700     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
701     PetscFunctionReturn(0);
702   }
703   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
704   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
705   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
706   snes->norm = fnorm;
707   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
708   SNESLogConvHistory(snes,fnorm,0);
709   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
710 
711   /* set parameter for default relative tolerance convergence test */
712   snes->ttol = fnorm*snes->rtol;
713   /* test convergence */
714   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
715   if (snes->reason) PetscFunctionReturn(0);
716   for (i = 0; i < maxits; i++) {
717     /* Call general purpose update function */
718 
719     if (snes->ops->update) {
720       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
721     }
722     ierr = FASCycle_Private(snes, X);CHKERRQ(ierr);
723     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
724     /* Monitor convergence */
725     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
726     snes->iter = i+1;
727     snes->norm = fnorm;
728     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
729     SNESLogConvHistory(snes,snes->norm,0);
730     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
731     /* Test for convergence */
732     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
733     if (snes->reason) break;
734   }
735   if (i == maxits) {
736     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
737     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
738   }
739   PetscFunctionReturn(0);
740 }
741