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