xref: /petsc/src/snes/impls/fas/fas.c (revision b51b94fa95cb114ddc4bf62a84dee696dbc4c345)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3 
4 const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0};
5 
6 /*MC
7 
8 SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
9 
10 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
11 
12 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
13 M*/
14 
15 extern PetscErrorCode SNESDestroy_FAS(SNES snes);
16 extern PetscErrorCode SNESSetUp_FAS(SNES snes);
17 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
18 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
19 extern PetscErrorCode SNESSolve_FAS(SNES snes);
20 extern PetscErrorCode SNESReset_FAS(SNES snes);
21 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
22 
23 EXTERN_C_BEGIN
24 #undef __FUNCT__
25 #define __FUNCT__ "SNESLineSearchSetType_FAS"
26 PetscErrorCode  SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type)
27 {
28   PetscErrorCode ierr;
29   PetscFunctionBegin;
30 
31   switch (type) {
32   case SNES_LS_BASIC:
33     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
34     break;
35   case SNES_LS_BASIC_NONORMS:
36     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
37     break;
38   case SNES_LS_QUADRATIC:
39     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
40     break;
41   default:
42     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type.");
43     break;
44   }
45   snes->ls_type = type;
46   PetscFunctionReturn(0);
47 }
48 EXTERN_C_END
49 
50 EXTERN_C_BEGIN
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "SNESCreate_FAS"
54 PetscErrorCode SNESCreate_FAS(SNES snes)
55 {
56   SNES_FAS * fas;
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   snes->ops->destroy        = SNESDestroy_FAS;
61   snes->ops->setup          = SNESSetUp_FAS;
62   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
63   snes->ops->view           = SNESView_FAS;
64   snes->ops->solve          = SNESSolve_FAS;
65   snes->ops->reset          = SNESReset_FAS;
66 
67   snes->usesksp             = PETSC_FALSE;
68   snes->usespc              = PETSC_FALSE;
69 
70   snes->max_funcs = 30000;
71   snes->max_its   = 10000;
72 
73   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
74   snes->data                  = (void*) fas;
75   fas->level                  = 0;
76   fas->levels                 = 1;
77   fas->n_cycles               = 1;
78   fas->max_up_it              = 1;
79   fas->max_down_it            = 1;
80   fas->upsmooth               = PETSC_NULL;
81   fas->downsmooth             = PETSC_NULL;
82   fas->next                   = PETSC_NULL;
83   fas->previous               = PETSC_NULL;
84   fas->interpolate            = PETSC_NULL;
85   fas->restrct                = PETSC_NULL;
86   fas->inject                 = PETSC_NULL;
87   fas->monitor                = PETSC_NULL;
88   fas->usedmfornumberoflevels = PETSC_FALSE;
89   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
90 
91   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr);
92   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
93 
94   PetscFunctionReturn(0);
95 }
96 EXTERN_C_END
97 
98 #undef __FUNCT__
99 #define __FUNCT__ "SNESFASGetLevels"
100 /*@
101    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
102 
103    Input Parameter:
104 .  snes - the nonlinear solver context
105 
106    Output parameter:
107 .  levels - the number of levels
108 
109    Level: advanced
110 
111 .keywords: MG, get, levels, multigrid
112 
113 .seealso: SNESFASSetLevels(), PCMGGetLevels()
114 @*/
115 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
116   SNES_FAS * fas = (SNES_FAS *)snes->data;
117   PetscFunctionBegin;
118   *levels = fas->levels;
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "SNESFASSetCycles"
124 /*@
125    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
126    complicated cycling.
127 
128    Logically Collective on SNES
129 
130    Input Parameters:
131 +  snes   - the multigrid context
132 -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
133 
134    Options Database Key:
135 $  -snes_fas_cycles 1 or 2
136 
137    Level: advanced
138 
139 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
140 
141 .seealso: SNESFASSetCyclesOnLevel()
142 @*/
143 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
144   SNES_FAS * fas = (SNES_FAS *)snes->data;
145   PetscErrorCode ierr;
146   PetscFunctionBegin;
147   fas->n_cycles = cycles;
148   if (fas->next) {
149     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
150   }
151   PetscFunctionReturn(0);
152 }
153 
154 #undef __FUNCT__
155 #define __FUNCT__ "SNESFASSetCyclesOnLevel"
156 /*@
157    SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level.
158 
159    Logically Collective on SNES
160 
161    Input Parameters:
162 +  snes   - the multigrid context
163 .  level  - the level to set the number of cycles on
164 -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
165 
166    Level: advanced
167 
168 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
169 
170 .seealso: SNESFASSetCycles()
171 @*/
172 PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
173   SNES_FAS * fas =  (SNES_FAS *)snes->data;
174   PetscInt top_level = fas->level,i;
175 
176   PetscFunctionBegin;
177   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D cycle count from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
178   /* get to the correct level */
179   for (i = fas->level; i > level; i--) {
180     fas = (SNES_FAS *)fas->next->data;
181   }
182   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
183   fas->n_cycles = cycles;
184   PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "SNESFASSetGS"
189 /*@C
190    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
191    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
192    and nonlinear preconditioners.
193 
194    Logically Collective on SNES
195 
196    Input Parameters:
197 +  snes    - the multigrid context
198 .  gsfunc  - the nonlinear smoother function
199 .  ctx     - the user context for the nonlinear smoother
200 -  use_gs  - whether to use the nonlinear smoother or not
201 
202    Level: advanced
203 
204 .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
205 
206 .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
207 @*/
208 PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
209   PetscErrorCode ierr = 0;
210   SNES_FAS       *fas = (SNES_FAS *)snes->data;
211   PetscFunctionBegin;
212 
213   if (gsfunc) {
214     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
215     /* push the provided GS up the tree */
216     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
217   } else {
218     ierr = SNESGetGS(snes,&gsfunc,&ctx);CHKERRQ(ierr);
219     if (gsfunc) {
220       /* assume that the user has set the GS solver at this level */
221       if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
222     } else if (use_gs) {
223       SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level);
224     }
225   }
226   PetscFunctionReturn(0);
227 }
228 
229 #undef __FUNCT__
230 #define __FUNCT__ "SNESFASSetGSOnLevel"
231 /*@C
232    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
233 
234    Logically Collective on SNES
235 
236    Input Parameters:
237 +  snes    - the multigrid context
238 .  level   - the level to set the nonlinear smoother on
239 .  gsfunc  - the nonlinear smoother function
240 .  ctx     - the user context for the nonlinear smoother
241 -  use_gs  - whether to use the nonlinear smoother or not
242 
243    Level: advanced
244 
245 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
246 
247 .seealso: SNESSetGS(), SNESFASSetGS()
248 @*/
249 PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
250   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
251   PetscErrorCode ierr;
252   PetscInt       top_level = fas->level,i;
253   SNES           cur_snes = snes;
254   PetscFunctionBegin;
255   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D GS from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
256   /* get to the correct level */
257   for (i = fas->level; i > level; i--) {
258     fas = (SNES_FAS *)fas->next->data;
259     cur_snes = fas->next;
260   }
261   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
262   if (gsfunc) {
263     ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr);
264   }
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "SNESFASGetSNES"
270 /*@
271    SNESFASGetSNES - Gets the SNES corresponding to a particular
272    level of the FAS hierarchy.
273 
274    Input Parameters:
275 +  snes    - the multigrid context
276    level   - the level to get
277 -  lsnes   - whether to use the nonlinear smoother or not
278 
279    Level: advanced
280 
281 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
282 
283 .seealso: SNESFASSetLevels(), SNESFASGetLevels()
284 @*/
285 PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) {
286   SNES_FAS *fas = (SNES_FAS*)snes->data;
287   PetscInt i;
288 
289   PetscFunctionBegin;
290   if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
291   if (fas->level < level) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS on level %D, must call from SNES at least as fine as level",level,fas->level);
292   *lsnes = snes;
293   for (i = fas->level; i > level; i--) {
294     *lsnes = fas->next;
295     fas = (SNES_FAS *)(*lsnes)->data;
296   }
297   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
298   PetscFunctionReturn(0);
299 }
300 
301 #undef __FUNCT__
302 #define __FUNCT__ "SNESFASSetType"
303 /*@
304 SNESFASSetType - Sets the update and correction type used for FAS.
305 e
306 
307 
308 @*/
309 PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
310 {
311   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
312 
313   PetscFunctionBegin;
314   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
315   PetscValidLogicalCollectiveEnum(snes,fastype,2);
316   fas->fastype = fastype;
317   PetscFunctionReturn(0);
318 }
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "SNESFASSetLevels"
322 /*@C
323    SNESFASSetLevels - Sets the number of levels to use with FAS.
324    Must be called before any other FAS routine.
325 
326    Input Parameters:
327 +  snes   - the snes context
328 .  levels - the number of levels
329 -  comms  - optional communicators for each level; this is to allow solving the coarser
330             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
331             Fortran.
332 
333    Level: intermediate
334 
335    Notes:
336      If the number of levels is one then the multigrid uses the -fas_levels prefix
337   for setting the level options rather than the -fas_coarse prefix.
338 
339 .keywords: FAS, MG, set, levels, multigrid
340 
341 .seealso: SNESFASGetLevels()
342 @*/
343 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
344   PetscErrorCode ierr;
345   PetscInt i;
346   SNES_FAS * fas = (SNES_FAS *)snes->data;
347   SNES prevsnes;
348   MPI_Comm comm;
349   PetscFunctionBegin;
350   comm = ((PetscObject)snes)->comm;
351   if (levels == fas->levels) {
352     if (!comms)
353       PetscFunctionReturn(0);
354   }
355   /* user has changed the number of levels; reset */
356   ierr = SNESReset(snes);CHKERRQ(ierr);
357   /* destroy any coarser levels if necessary */
358   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
359   fas->next = PETSC_NULL;
360   fas->previous = PETSC_NULL;
361   prevsnes = snes;
362   /* setup the finest level */
363   for (i = levels - 1; i >= 0; i--) {
364     if (comms) comm = comms[i];
365     fas->level = i;
366     fas->levels = levels;
367     fas->next = PETSC_NULL;
368     if (i > 0) {
369       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
370       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
371       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
372       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
373       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
374       prevsnes = fas->next;
375       fas = (SNES_FAS *)prevsnes->data;
376     }
377   }
378   PetscFunctionReturn(0);
379 }
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "SNESFASSetNumberSmoothUp"
383 /*@
384    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
385    use on all levels.
386 
387    Logically Collective on SNES
388 
389    Input Parameters:
390 +  snes - the multigrid context
391 -  n    - the number of smoothing steps
392 
393    Options Database Key:
394 .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
395 
396    Level: advanced
397 
398 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
399 
400 .seealso: SNESFASSetNumberSmoothDown()
401 @*/
402 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
403   SNES_FAS * fas =  (SNES_FAS *)snes->data;
404   PetscErrorCode ierr = 0;
405   PetscFunctionBegin;
406   fas->max_up_it = n;
407   if (fas->next) {
408     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
409   }
410   PetscFunctionReturn(0);
411 }
412 
413 #undef __FUNCT__
414 #define __FUNCT__ "SNESFASSetNumberSmoothDown"
415 /*@
416    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
417    use on all levels.
418 
419    Logically Collective on SNES
420 
421    Input Parameters:
422 +  snes - the multigrid context
423 -  n    - the number of smoothing steps
424 
425    Options Database Key:
426 .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
427 
428    Level: advanced
429 
430 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
431 
432 .seealso: SNESFASSetNumberSmoothUp()
433 @*/
434 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
435   SNES_FAS * fas =  (SNES_FAS *)snes->data;
436   PetscErrorCode ierr = 0;
437   PetscFunctionBegin;
438   fas->max_down_it = n;
439   if (fas->next) {
440     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
441   }
442   PetscFunctionReturn(0);
443 }
444 
445 #undef __FUNCT__
446 #define __FUNCT__ "SNESFASSetInterpolation"
447 /*@
448    SNESFASSetInterpolation - Sets the function to be used to calculate the
449    interpolation from l-1 to the lth level
450 
451    Input Parameters:
452 +  snes      - the multigrid context
453 .  mat       - the interpolation operator
454 -  level     - the level (0 is coarsest) to supply [do not supply 0]
455 
456    Level: advanced
457 
458    Notes:
459           Usually this is the same matrix used also to set the restriction
460     for the same level.
461 
462           One can pass in the interpolation matrix or its transpose; PETSc figures
463     out from the matrix size which one it is.
464 
465 .keywords:  FAS, multigrid, set, interpolate, level
466 
467 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
468 @*/
469 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
470   SNES_FAS * fas =  (SNES_FAS *)snes->data;
471   PetscInt top_level = fas->level,i;
472   PetscErrorCode ierr;
473 
474   PetscFunctionBegin;
475   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D interpolation from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
476   /* get to the correct level */
477   for (i = fas->level; i > level; i--) {
478     fas = (SNES_FAS *)fas->next->data;
479   }
480   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
481   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
482   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
483   fas->interpolate = mat;
484   PetscFunctionReturn(0);
485 }
486 
487 #undef __FUNCT__
488 #define __FUNCT__ "SNESFASSetRestriction"
489 /*@
490    SNESFASSetRestriction - Sets the function to be used to restrict the defect
491    from level l to l-1.
492 
493    Input Parameters:
494 +  snes  - the multigrid context
495 .  mat   - the restriction matrix
496 -  level - the level (0 is coarsest) to supply [Do not supply 0]
497 
498    Level: advanced
499 
500    Notes:
501           Usually this is the same matrix used also to set the interpolation
502     for the same level.
503 
504           One can pass in the interpolation matrix or its transpose; PETSc figures
505     out from the matrix size which one it is.
506 
507          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
508     is used.
509 
510 .keywords: FAS, MG, set, multigrid, restriction, level
511 
512 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
513 @*/
514 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
515   SNES_FAS * fas =  (SNES_FAS *)snes->data;
516   PetscInt top_level = fas->level,i;
517   PetscErrorCode ierr;
518 
519   PetscFunctionBegin;
520   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D restriction from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
521   /* get to the correct level */
522   for (i = fas->level; i > level; i--) {
523     fas = (SNES_FAS *)fas->next->data;
524   }
525   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
526   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
527   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
528   fas->restrct = mat;
529   PetscFunctionReturn(0);
530 }
531 
532 #undef __FUNCT__
533 #define __FUNCT__ "SNESFASSetRScale"
534 /*@
535    SNESFASSetRScale - Sets the scaling factor of the restriction
536    operator from level l to l-1.
537 
538    Input Parameters:
539 +  snes   - the multigrid context
540 .  rscale - the restriction scaling
541 -  level  - the level (0 is coarsest) to supply [Do not supply 0]
542 
543    Level: advanced
544 
545    Notes:
546          This is only used in the case that the injection is not set.
547 
548 .keywords: FAS, MG, set, multigrid, restriction, level
549 
550 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
551 @*/
552 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
553   SNES_FAS * fas =  (SNES_FAS *)snes->data;
554   PetscInt top_level = fas->level,i;
555   PetscErrorCode ierr;
556 
557   PetscFunctionBegin;
558   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D rscale from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
559   /* get to the correct level */
560   for (i = fas->level; i > level; i--) {
561     fas = (SNES_FAS *)fas->next->data;
562   }
563   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
564   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
565   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
566   fas->rscale = rscale;
567   PetscFunctionReturn(0);
568 }
569 
570 
571 #undef __FUNCT__
572 #define __FUNCT__ "SNESFASSetInjection"
573 /*@
574    SNESFASSetInjection - Sets the function to be used to inject the solution
575    from level l to l-1.
576 
577    Input Parameters:
578 +  snes  - the multigrid context
579 .  mat   - the restriction matrix
580 -  level - the level (0 is coarsest) to supply [Do not supply 0]
581 
582    Level: advanced
583 
584    Notes:
585          If you do not set this, the restriction and rscale is used to
586    project the solution instead.
587 
588 .keywords: FAS, MG, set, multigrid, restriction, level
589 
590 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
591 @*/
592 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
593   SNES_FAS * fas =  (SNES_FAS *)snes->data;
594   PetscInt top_level = fas->level,i;
595   PetscErrorCode ierr;
596 
597   PetscFunctionBegin;
598   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D injection from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
599   /* get to the correct level */
600   for (i = fas->level; i > level; i--) {
601     fas = (SNES_FAS *)fas->next->data;
602   }
603   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
604   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
605   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
606   fas->inject = mat;
607   PetscFunctionReturn(0);
608 }
609 
610 #undef __FUNCT__
611 #define __FUNCT__ "SNESReset_FAS"
612 PetscErrorCode SNESReset_FAS(SNES snes)
613 {
614   PetscErrorCode ierr = 0;
615   SNES_FAS * fas = (SNES_FAS *)snes->data;
616 
617   PetscFunctionBegin;
618   ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
619   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
620   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
621   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
622   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
623   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
624 
625   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
626   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
627   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
628   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","",PETSC_NULL);CHKERRQ(ierr);
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "SNESDestroy_FAS"
634 PetscErrorCode SNESDestroy_FAS(SNES snes)
635 {
636   SNES_FAS * fas = (SNES_FAS *)snes->data;
637   PetscErrorCode ierr = 0;
638 
639   PetscFunctionBegin;
640   /* recursively resets and then destroys */
641   ierr = SNESReset(snes);CHKERRQ(ierr);
642   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
643   ierr = PetscFree(fas);CHKERRQ(ierr);
644 
645   PetscFunctionReturn(0);
646 }
647 
648 #undef __FUNCT__
649 #define __FUNCT__ "SNESSetUp_FAS"
650 PetscErrorCode SNESSetUp_FAS(SNES snes)
651 {
652   SNES_FAS       *fas = (SNES_FAS *) snes->data;
653   PetscErrorCode ierr;
654   VecScatter     injscatter;
655   PetscInt       dm_levels;
656   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
657 
658   PetscFunctionBegin;
659 
660   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
661     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
662     dm_levels++;
663     if (dm_levels > fas->levels) {
664 
665       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
666       vec_sol = snes->vec_sol;
667       vec_func = snes->vec_func;
668       vec_sol_update = snes->vec_sol_update;
669       vec_rhs = snes->vec_rhs;
670       snes->vec_sol = PETSC_NULL;
671       snes->vec_func = PETSC_NULL;
672       snes->vec_sol_update = PETSC_NULL;
673       snes->vec_rhs = PETSC_NULL;
674 
675       /* reset the number of levels */
676       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
677       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
678 
679       snes->vec_sol = vec_sol;
680       snes->vec_func = vec_func;
681       snes->vec_rhs = vec_rhs;
682       snes->vec_sol_update = vec_sol_update;
683     }
684   }
685 
686   if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
687 
688   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
689     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
690   } else {
691     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
692   }
693 
694   if (snes->dm) {
695     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
696     if (fas->next) {
697       /* for now -- assume the DM and the evaluation functions have been set externally */
698       if (!fas->next->dm) {
699         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
700         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
701       }
702       /* set the interpolation and restriction from the DM */
703       if (!fas->interpolate) {
704         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
705         if (!fas->restrct) {
706           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
707           fas->restrct = fas->interpolate;
708         }
709       }
710       /* set the injection from the DM */
711       if (!fas->inject) {
712         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
713         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
714         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
715       }
716     }
717     /* set the DMs of the pre and post-smoothers here */
718     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
719     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
720   }
721   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
722 
723  if (fas->next) {
724     if (fas->galerkin) {
725       ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr);
726     }
727   }
728 
729   if (fas->next) {
730     /* gotta set up the solution vector for this to work */
731     if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);}
732     if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);}
733     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
734   }
735 
736   /* setup the pre and post smoothers */
737   if (fas->upsmooth) {ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);}
738   if (fas->downsmooth) {ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);}
739 
740   /* setup FAS work vectors */
741   if (fas->galerkin) {
742     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
743     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
744   }
745   PetscFunctionReturn(0);
746 }
747 
748 #undef __FUNCT__
749 #define __FUNCT__ "SNESSetFromOptions_FAS"
750 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
751 {
752   SNES_FAS       *fas = (SNES_FAS *) snes->data;
753   PetscInt       levels = 1;
754   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg;
755   PetscErrorCode ierr;
756   char           monfilename[PETSC_MAX_PATH_LEN];
757   SNESFASType    fastype;
758   const char     *optionsprefix;
759   const char     *prefix;
760 
761   PetscFunctionBegin;
762   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
763 
764   /* number of levels -- only process on the finest level */
765   if (fas->levels - 1 == fas->level) {
766     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
767     if (!flg && snes->dm) {
768       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
769       levels++;
770       fas->usedmfornumberoflevels = PETSC_TRUE;
771     }
772     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
773   }
774   fastype = fas->fastype;
775   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
776   if (flg) {
777     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
778   }
779 
780   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
781 
782   /* smoother setup options */
783   ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr);
784   ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr);
785   ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr);
786   if (fas->level == 0) {
787     ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr);
788   }
789   /* options for the number of preconditioning cycles and cycle type */
790 
791   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
792 
793   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
794 
795   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
796 
797 
798   ierr = PetscOptionsTail();CHKERRQ(ierr);
799   /* setup from the determined types if there is no pointwise procedure or smoother defined */
800 
801   if ((!fas->downsmooth) && smoothcoarseflg) {
802     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
803     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
804     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
805     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
806     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
807   }
808 
809   if ((!fas->downsmooth) && smoothdownflg) {
810     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
811     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
812     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
813     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr);
814     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
815   }
816 
817   if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) {
818     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
819     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
820     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
821     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr);
822     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
823   }
824 
825   if ((!fas->downsmooth) && smoothflg) {
826     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
827     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
828     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
829     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr);
830     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
831   }
832 
833   if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) {
834     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
835     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
836     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
837     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr);
838     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
839   }
840 
841   if (fas->upsmooth) {
842     ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr);
843   }
844 
845   if (fas->downsmooth) {
846     ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr);
847   }
848 
849   if (fas->level != fas->levels - 1) {
850     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr);
851   }
852 
853   /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */
854   if (!fas->downsmooth) {
855     ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
856     ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
857     if (fas->level == 0) {
858       ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
859     }
860   }
861 
862   if (!fas->upsmooth) {
863     ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
864     ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
865   }
866 
867   if (monflg) {
868     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
869     /* set the monitors for the upsmoother and downsmoother */
870     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
871     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
872     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
873     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
874   } else {
875     /* unset the monitors on the coarse levels */
876     if (fas->level != fas->levels - 1) {
877       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
878     }
879   }
880 
881   /* recursive option setting for the smoothers */
882   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
883   PetscFunctionReturn(0);
884 }
885 
886 #undef __FUNCT__
887 #define __FUNCT__ "SNESView_FAS"
888 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
889 {
890   SNES_FAS   *fas = (SNES_FAS *) snes->data;
891   PetscBool      iascii;
892   PetscErrorCode ierr;
893 
894   PetscFunctionBegin;
895   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
896   if (iascii) {
897     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n",  fas->levels);CHKERRQ(ierr);
898     ierr = PetscViewerASCIIPushTab(viewer);
899     ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n",  fas->level);CHKERRQ(ierr);
900     if (fas->upsmooth) {
901       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
902       ierr = PetscViewerASCIIPushTab(viewer);
903       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
904       ierr = PetscViewerASCIIPopTab(viewer);
905     } else {
906       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
907     }
908     if (fas->downsmooth) {
909       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
910       ierr = PetscViewerASCIIPushTab(viewer);
911       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
912       ierr = PetscViewerASCIIPopTab(viewer);
913     } else {
914       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
915     }
916     ierr = PetscViewerASCIIPopTab(viewer);
917   } else {
918     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
919   }
920   PetscFunctionReturn(0);
921 }
922 
923 #undef __FUNCT__
924 #define __FUNCT__ "FASDownSmooth"
925 /*
926 Defines the action of the downsmoother
927  */
928 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
929   PetscErrorCode      ierr = 0;
930   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
931   SNESConvergedReason reason;
932   SNES_FAS            *fas = (SNES_FAS *)snes->data;
933   Vec                 G, W, Y, FPC;
934   PetscBool           lssuccess;
935   PetscInt            k;
936   PetscFunctionBegin;
937   G = snes->work[1];
938   W = snes->work[2];
939   Y = snes->work[3];
940   if (fas->downsmooth) {
941     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
942     /* check convergence reason for the smoother */
943     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
944     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
945       snes->reason = SNES_DIVERGED_INNER;
946       PetscFunctionReturn(0);
947     }
948     ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
949     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
950   } else {
951     /* basic richardson smoother */
952     for (k = 0; k < fas->max_up_it; k++) {
953       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
954       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
955       ierr = VecCopy(F, Y);CHKERRQ(ierr);
956       ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);CHKERRQ(ierr);
957       ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
958       if (!lssuccess) {
959         if (++snes->numFailures >= snes->maxFailures) {
960           snes->reason = SNES_DIVERGED_LINE_SEARCH;
961           PetscFunctionReturn(0);
962         }
963       }
964       ierr = VecCopy(W, X);CHKERRQ(ierr);
965       ierr = VecCopy(G, F);CHKERRQ(ierr);
966       fnorm = gnorm;
967     }
968   }
969   PetscFunctionReturn(0);
970 }
971 
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "FASUpSmooth"
975 /*
976 Defines the action of the upsmoother
977  */
978 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
979   PetscErrorCode      ierr = 0;
980   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
981   SNESConvergedReason reason;
982   SNES_FAS            *fas = (SNES_FAS *)snes->data;
983   Vec                 G, W, Y, FPC;
984   PetscBool           lssuccess;
985   PetscInt            k;
986   PetscFunctionBegin;
987   G = snes->work[1];
988   W = snes->work[2];
989   Y = snes->work[3];
990   if (fas->upsmooth) {
991     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
992     /* check convergence reason for the smoother */
993     ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr);
994     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
995       snes->reason = SNES_DIVERGED_INNER;
996       PetscFunctionReturn(0);
997     }
998     ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
999     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
1000   } else {
1001     /* basic richardson smoother */
1002     for (k = 0; k < fas->max_up_it; k++) {
1003       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1004       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1005       ierr = VecCopy(F, Y);CHKERRQ(ierr);
1006       ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);CHKERRQ(ierr);
1007       ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
1008       if (!lssuccess) {
1009         if (++snes->numFailures >= snes->maxFailures) {
1010           snes->reason = SNES_DIVERGED_LINE_SEARCH;
1011           PetscFunctionReturn(0);
1012         }
1013       }
1014       ierr = VecCopy(W, X);CHKERRQ(ierr);
1015       ierr = VecCopy(G, F);CHKERRQ(ierr);
1016       fnorm = gnorm;
1017     }
1018   }
1019   PetscFunctionReturn(0);
1020 }
1021 
1022 #undef __FUNCT__
1023 #define __FUNCT__ "SNESFASCreateCoarseVec"
1024 /*@
1025    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
1026 
1027    Collective
1028 
1029    Input Arguments:
1030 .  snes - SNESFAS
1031 
1032    Output Arguments:
1033 .  Xcoarse - vector on level one coarser than snes
1034 
1035    Level: developer
1036 
1037 .seealso: SNESFASSetRestriction(), SNESFASRestrict()
1038 @*/
1039 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
1040 {
1041   PetscErrorCode ierr;
1042   SNES_FAS       *fas = (SNES_FAS*)snes->data;
1043 
1044   PetscFunctionBegin;
1045   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
1046   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
1047   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 #undef __FUNCT__
1052 #define __FUNCT__ "SNESFASRestrict"
1053 /*@
1054    SNESFASRestrict - restrict a Vec to the next coarser level
1055 
1056    Collective
1057 
1058    Input Arguments:
1059 +  fine - SNES from which to restrict
1060 -  Xfine - vector to restrict
1061 
1062    Output Arguments:
1063 .  Xcoarse - result of restriction
1064 
1065    Level: developer
1066 
1067 .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
1068 @*/
1069 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
1070 {
1071   PetscErrorCode ierr;
1072   SNES_FAS       *fas = (SNES_FAS*)fine->data;
1073 
1074   PetscFunctionBegin;
1075   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
1076   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
1077   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
1078   if (fas->inject) {
1079     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
1080   } else {
1081     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
1082     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
1083   }
1084   PetscFunctionReturn(0);
1085 }
1086 
1087 #undef __FUNCT__
1088 #define __FUNCT__ "FASCoarseCorrection"
1089 /*
1090 
1091 Performs the FAS coarse correction as:
1092 
1093 fine problem: F(x) = 0
1094 coarse problem: F^c(x) = b^c
1095 
1096 b^c = F^c(I^c_fx^f - I^c_fF(x))
1097 
1098 with correction:
1099 
1100 
1101 
1102  */
1103 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
1104   PetscErrorCode      ierr;
1105   Vec                 X_c, Xo_c, F_c, B_c;
1106   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1107   SNESConvergedReason reason;
1108   PetscFunctionBegin;
1109   if (fas->next) {
1110     X_c  = fas->next->vec_sol;
1111     Xo_c = fas->next->work[0];
1112     F_c  = fas->next->vec_func;
1113     B_c  = fas->next->vec_rhs;
1114 
1115     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
1116     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1117 
1118     /* restrict the defect */
1119     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1120 
1121     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1122     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1123     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1124 
1125     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1126 
1127     /* set initial guess of the coarse problem to the projected fine solution */
1128     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1129 
1130     /* recurse to the next level */
1131     fas->next->vec_rhs = B_c;
1132     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1133     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1134     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1135     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1136       snes->reason = SNES_DIVERGED_INNER;
1137       PetscFunctionReturn(0);
1138     }
1139     /* fas->next->vec_rhs = PETSC_NULL; */
1140 
1141     /* correct as x <- x + I(x^c - Rx)*/
1142     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1143     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1144   }
1145   PetscFunctionReturn(0);
1146 }
1147 
1148 #undef __FUNCT__
1149 #define __FUNCT__ "FASCycle_Additive"
1150 /*
1151 
1152 The additive cycle looks like:
1153 
1154 xhat = x
1155 xhat = dS(x, b)
1156 x = coarsecorrection(xhat, b_d)
1157 x = x + nu*(xhat - x);
1158 (optional) x = uS(x, b)
1159 
1160 With the coarse RHS (defect correction) as below.
1161 
1162  */
1163 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
1164   Vec                 F, B, Xhat;
1165   Vec                 X_c, Xo_c, F_c, B_c, G, W;
1166   PetscErrorCode      ierr;
1167   SNES_FAS *          fas = (SNES_FAS *)snes->data;
1168   SNESConvergedReason reason;
1169   PetscReal           xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.;
1170   PetscBool           lssucceed;
1171   PetscFunctionBegin;
1172 
1173   F = snes->vec_func;
1174   B = snes->vec_rhs;
1175   Xhat = snes->work[3];
1176   G    = snes->work[1];
1177   W    = snes->work[2];
1178   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
1179   /* recurse first */
1180   if (fas->next) {
1181     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
1182 
1183     X_c  = fas->next->vec_sol;
1184     Xo_c = fas->next->work[0];
1185     F_c  = fas->next->vec_func;
1186     B_c  = fas->next->vec_rhs;
1187 
1188     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
1189     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1190 
1191     /* restrict the defect */
1192     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1193 
1194     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1195     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1196     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1197 
1198     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1199 
1200     /* set initial guess of the coarse problem to the projected fine solution */
1201     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1202 
1203     /* recurse */
1204     fas->next->vec_rhs = B_c;
1205     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1206 
1207     /* smooth on this level */
1208     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1209 
1210     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1211     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1212       snes->reason = SNES_DIVERGED_INNER;
1213       PetscFunctionReturn(0);
1214     }
1215 
1216     /* correct as x <- x + I(x^c - Rx)*/
1217     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1218     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
1219 
1220     /* additive correction of the coarse direction*/
1221     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1222     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1223     ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr);
1224     ierr = SNESLineSearchPreCheckApply(snes, X,Xhat,PETSC_NULL);CHKERRQ(ierr);
1225     ierr = SNESLineSearchApply(snes,X,F,Xhat,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1226     ierr = VecCopy(W, X);CHKERRQ(ierr);
1227     ierr = VecCopy(G, F);CHKERRQ(ierr);
1228     fnorm = gnorm;
1229   } else {
1230     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1231   }
1232   PetscFunctionReturn(0);
1233 }
1234 
1235 #undef __FUNCT__
1236 #define __FUNCT__ "FASCycle_Multiplicative"
1237 /*
1238 
1239 Defines the FAS cycle as:
1240 
1241 fine problem: F(x) = 0
1242 coarse problem: F^c(x) = b^c
1243 
1244 b^c = F^c(I^c_fx^f - I^c_fF(x))
1245 
1246 correction:
1247 
1248 x = x + I(x^c - Rx)
1249 
1250  */
1251 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
1252 
1253   PetscErrorCode      ierr;
1254   Vec                 F,B;
1255   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1256 
1257   PetscFunctionBegin;
1258   F = snes->vec_func;
1259   B = snes->vec_rhs;
1260   /* pre-smooth -- just update using the pre-smoother */
1261   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1262 
1263   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
1264 
1265   if (fas->level != 0) {
1266     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1267   }
1268   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1269 
1270   PetscFunctionReturn(0);
1271 }
1272 
1273 #undef __FUNCT__
1274 #define __FUNCT__ "SNESSolve_FAS"
1275 
1276 PetscErrorCode SNESSolve_FAS(SNES snes)
1277 {
1278   PetscErrorCode ierr;
1279   PetscInt       i, maxits;
1280   Vec            X, F;
1281   PetscReal      fnorm;
1282   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
1283   DM             dm;
1284 
1285   PetscFunctionBegin;
1286   maxits = snes->max_its;            /* maximum number of iterations */
1287   snes->reason = SNES_CONVERGED_ITERATING;
1288   X = snes->vec_sol;
1289   F = snes->vec_func;
1290 
1291   /*norm setup */
1292   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1293   snes->iter = 0;
1294   snes->norm = 0.;
1295   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1296   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1297   if (snes->domainerror) {
1298     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1299     PetscFunctionReturn(0);
1300   }
1301   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1302   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1303   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1304   snes->norm = fnorm;
1305   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1306   SNESLogConvHistory(snes,fnorm,0);
1307   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1308 
1309   /* set parameter for default relative tolerance convergence test */
1310   snes->ttol = fnorm*snes->rtol;
1311   /* test convergence */
1312   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1313   if (snes->reason) PetscFunctionReturn(0);
1314 
1315   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1316   for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
1317     DM dmcoarse;
1318     ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
1319     ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
1320     dm = dmcoarse;
1321   }
1322 
1323   for (i = 0; i < maxits; i++) {
1324     /* Call general purpose update function */
1325 
1326     if (snes->ops->update) {
1327       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1328     }
1329     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
1330       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
1331     } else {
1332       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
1333     }
1334 
1335     /* check for FAS cycle divergence */
1336     if (snes->reason != SNES_CONVERGED_ITERATING) {
1337       PetscFunctionReturn(0);
1338     }
1339     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1340     /* Monitor convergence */
1341     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1342     snes->iter = i+1;
1343     snes->norm = fnorm;
1344     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1345     SNESLogConvHistory(snes,snes->norm,0);
1346     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1347     /* Test for convergence */
1348     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1349     if (snes->reason) break;
1350   }
1351   if (i == maxits) {
1352     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1353     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1354   }
1355   PetscFunctionReturn(0);
1356 }
1357