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