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