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