xref: /petsc/src/snes/impls/fas/fas.c (revision e9923e8da335e25b76d579f89846e60d092e1117)
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 = SNESSetUp(fas->next);CHKERRQ(ierr);
738   }
739 
740   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
741   if (fas->upsmooth) {
742     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
743       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
744     }
745     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
746       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
747     }
748    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
749       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
750     }
751    ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
752   }
753   if (fas->downsmooth) {
754     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
755       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
756     }
757     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
758       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
759     }
760    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
761      ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
762     }
763     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
764   }
765 
766   /* setup FAS work vectors */
767   if (fas->galerkin) {
768     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
769     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
770   }
771 
772 
773   /* got to set them all up at once */
774   PetscFunctionReturn(0);
775 }
776 
777 #undef __FUNCT__
778 #define __FUNCT__ "SNESSetFromOptions_FAS"
779 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
780 {
781   SNES_FAS       *fas = (SNES_FAS *) snes->data;
782   PetscInt       levels = 1;
783   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg;
784   PetscErrorCode ierr;
785   char           monfilename[PETSC_MAX_PATH_LEN];
786   SNESFASType    fastype;
787   const char     *optionsprefix;
788   const char     *prefix;
789 
790   PetscFunctionBegin;
791   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
792 
793   /* number of levels -- only process on the finest level */
794   if (fas->levels - 1 == fas->level) {
795     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
796     if (!flg && snes->dm) {
797       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
798       levels++;
799       fas->usedmfornumberoflevels = PETSC_TRUE;
800     }
801     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
802   }
803   fastype = fas->fastype;
804   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
805   if (flg) {
806     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
807   }
808 
809   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
810 
811   /* smoother setup options */
812   ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr);
813   ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr);
814   ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr);
815   if (fas->level == 0) {
816     ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr);
817   }
818   /* options for the number of preconditioning cycles and cycle type */
819 
820   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
821 
822   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
823 
824   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
825 
826 
827   ierr = PetscOptionsTail();CHKERRQ(ierr);
828   /* setup from the determined types if there is no pointwise procedure or smoother defined */
829 
830   if ((!fas->downsmooth) && smoothcoarseflg) {
831     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
832     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
833     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
834     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
835     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
836   }
837 
838   if ((!fas->downsmooth) && smoothdownflg) {
839     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
840     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
841     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
842     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr);
843     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
844   }
845 
846   if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) {
847     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
848     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
849     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
850     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr);
851     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
852   }
853 
854   if ((!fas->downsmooth) && smoothflg) {
855     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
856     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
857     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
858     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr);
859     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
860   }
861 
862   if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) {
863     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
864     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
865     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
866     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr);
867     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
868   }
869 
870   if (fas->upsmooth) {
871     ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr);
872   }
873 
874   if (fas->downsmooth) {
875     ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr);
876   }
877 
878   if (fas->level != fas->levels - 1) {
879     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr);
880   }
881 
882   /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */
883   if (!fas->downsmooth) {
884     ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
885     ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
886     if (fas->level == 0) {
887       ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
888     }
889   }
890 
891   if (!fas->upsmooth) {
892     ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
893     ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
894   }
895 
896   if (monflg) {
897     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
898     /* set the monitors for the upsmoother and downsmoother */
899     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
900     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
901     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
902     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
903   } else {
904     /* unset the monitors on the coarse levels */
905     if (fas->level != fas->levels - 1) {
906       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
907     }
908   }
909 
910   /* recursive option setting for the smoothers */
911   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
912   PetscFunctionReturn(0);
913 }
914 
915 #undef __FUNCT__
916 #define __FUNCT__ "SNESView_FAS"
917 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
918 {
919   SNES_FAS   *fas = (SNES_FAS *) snes->data;
920   PetscBool      iascii;
921   PetscErrorCode ierr;
922 
923   PetscFunctionBegin;
924   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
925   if (iascii) {
926     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n",  fas->levels);CHKERRQ(ierr);
927     ierr = PetscViewerASCIIPushTab(viewer);
928     ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n",  fas->level);CHKERRQ(ierr);
929     if (fas->upsmooth) {
930       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
931       ierr = PetscViewerASCIIPushTab(viewer);
932       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
933       ierr = PetscViewerASCIIPopTab(viewer);
934     } else {
935       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
936     }
937     if (fas->downsmooth) {
938       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
939       ierr = PetscViewerASCIIPushTab(viewer);
940       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
941       ierr = PetscViewerASCIIPopTab(viewer);
942     } else {
943       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
944     }
945     ierr = PetscViewerASCIIPopTab(viewer);
946   } else {
947     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
948   }
949   PetscFunctionReturn(0);
950 }
951 
952 #undef __FUNCT__
953 #define __FUNCT__ "FASDownSmooth"
954 /*
955 Defines the action of the downsmoother
956  */
957 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
958   PetscErrorCode      ierr = 0;
959   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
960   SNESConvergedReason reason;
961   SNES_FAS            *fas = (SNES_FAS *)snes->data;
962   Vec                 G, W, Y, FPC;
963   PetscBool           lssuccess;
964   PetscInt            k;
965   PetscFunctionBegin;
966   G = snes->work[1];
967   W = snes->work[2];
968   Y = snes->work[3];
969   if (fas->downsmooth) {
970     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
971     /* check convergence reason for the smoother */
972     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
973     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
974       snes->reason = SNES_DIVERGED_INNER;
975       PetscFunctionReturn(0);
976     }
977     ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
978     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
979   } else {
980     /* basic richardson smoother */
981     for (k = 0; k < fas->max_up_it; k++) {
982       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
983       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
984       ierr = VecCopy(F, Y);CHKERRQ(ierr);
985       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
986       if (!lssuccess) {
987         if (++snes->numFailures >= snes->maxFailures) {
988           snes->reason = SNES_DIVERGED_LINE_SEARCH;
989           PetscFunctionReturn(0);
990         }
991       }
992       ierr = VecCopy(W, X);CHKERRQ(ierr);
993       ierr = VecCopy(G, F);CHKERRQ(ierr);
994       fnorm = gnorm;
995     }
996   }
997   PetscFunctionReturn(0);
998 }
999 
1000 
1001 #undef __FUNCT__
1002 #define __FUNCT__ "FASUpSmooth"
1003 /*
1004 Defines the action of the upsmoother
1005  */
1006 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
1007   PetscErrorCode      ierr = 0;
1008   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
1009   SNESConvergedReason reason;
1010   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1011   Vec                 G, W, Y, FPC;
1012   PetscBool           lssuccess;
1013   PetscInt            k;
1014   PetscFunctionBegin;
1015   G = snes->work[1];
1016   W = snes->work[2];
1017   Y = snes->work[3];
1018   if (fas->upsmooth) {
1019     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
1020     /* check convergence reason for the smoother */
1021     ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr);
1022     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1023       snes->reason = SNES_DIVERGED_INNER;
1024       PetscFunctionReturn(0);
1025     }
1026     ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
1027     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
1028   } else {
1029     /* basic richardson smoother */
1030     for (k = 0; k < fas->max_up_it; k++) {
1031       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1032       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1033       ierr = VecCopy(F, Y);CHKERRQ(ierr);
1034       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
1035       if (!lssuccess) {
1036         if (++snes->numFailures >= snes->maxFailures) {
1037           snes->reason = SNES_DIVERGED_LINE_SEARCH;
1038           PetscFunctionReturn(0);
1039         }
1040       }
1041       ierr = VecCopy(W, X);CHKERRQ(ierr);
1042       ierr = VecCopy(G, F);CHKERRQ(ierr);
1043       fnorm = gnorm;
1044     }
1045   }
1046   PetscFunctionReturn(0);
1047 }
1048 
1049 #undef __FUNCT__
1050 #define __FUNCT__ "SNESFASCreateCoarseVec"
1051 /*@
1052    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
1053 
1054    Collective
1055 
1056    Input Arguments:
1057 .  snes - SNESFAS
1058 
1059    Output Arguments:
1060 .  Xcoarse - vector on level one coarser than snes
1061 
1062    Level: developer
1063 
1064 .seealso: SNESFASSetRestriction(), SNESFASRestrict()
1065 @*/
1066 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
1067 {
1068   PetscErrorCode ierr;
1069   SNES_FAS       *fas = (SNES_FAS*)snes->data;
1070 
1071   PetscFunctionBegin;
1072   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
1073   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
1074   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
1075   PetscFunctionReturn(0);
1076 }
1077 
1078 #undef __FUNCT__
1079 #define __FUNCT__ "SNESFASRestrict"
1080 /*@
1081    SNESFASRestrict - restrict a Vec to the next coarser level
1082 
1083    Collective
1084 
1085    Input Arguments:
1086 +  fine - SNES from which to restrict
1087 -  Xfine - vector to restrict
1088 
1089    Output Arguments:
1090 .  Xcoarse - result of restriction
1091 
1092    Level: developer
1093 
1094 .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
1095 @*/
1096 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
1097 {
1098   PetscErrorCode ierr;
1099   SNES_FAS       *fas = (SNES_FAS*)fine->data;
1100 
1101   PetscFunctionBegin;
1102   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
1103   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
1104   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
1105   if (fas->inject) {
1106     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
1107   } else {
1108     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
1109     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
1110   }
1111   PetscFunctionReturn(0);
1112 }
1113 
1114 #undef __FUNCT__
1115 #define __FUNCT__ "FASCoarseCorrection"
1116 /*
1117 
1118 Performs the FAS coarse correction as:
1119 
1120 fine problem: F(x) = 0
1121 coarse problem: F^c(x) = b^c
1122 
1123 b^c = F^c(I^c_fx^f - I^c_fF(x))
1124 
1125 with correction:
1126 
1127 
1128 
1129  */
1130 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
1131   PetscErrorCode      ierr;
1132   Vec                 X_c, Xo_c, F_c, B_c;
1133   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1134   SNESConvergedReason reason;
1135   PetscFunctionBegin;
1136   if (fas->next) {
1137     X_c  = fas->next->vec_sol;
1138     Xo_c = fas->next->work[0];
1139     F_c  = fas->next->vec_func;
1140     B_c  = fas->next->vec_rhs;
1141 
1142     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
1143     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1144 
1145     /* restrict the defect */
1146     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1147 
1148     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1149     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1150     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1151 
1152     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1153 
1154     /* set initial guess of the coarse problem to the projected fine solution */
1155     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1156 
1157     /* recurse to the next level */
1158     fas->next->vec_rhs = B_c;
1159     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1160     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1161     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1162     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1163       snes->reason = SNES_DIVERGED_INNER;
1164       PetscFunctionReturn(0);
1165     }
1166     /* fas->next->vec_rhs = PETSC_NULL; */
1167 
1168     /* correct as x <- x + I(x^c - Rx)*/
1169     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1170     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1171   }
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 #undef __FUNCT__
1176 #define __FUNCT__ "FASCycle_Additive"
1177 /*
1178 
1179 The additive cycle looks like:
1180 
1181 xhat = x
1182 xhat = dS(x, b)
1183 x = coarsecorrection(xhat, b_d)
1184 x = x + nu*(xhat - x);
1185 (optional) x = uS(x, b)
1186 
1187 With the coarse RHS (defect correction) as below.
1188 
1189  */
1190 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
1191   Vec                 F, B, Xhat;
1192   Vec                 X_c, Xo_c, F_c, B_c, G, W;
1193   PetscErrorCode      ierr;
1194   SNES_FAS *          fas = (SNES_FAS *)snes->data;
1195   SNESConvergedReason reason;
1196   PetscReal           xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.;
1197   PetscBool           lssucceed;
1198   PetscFunctionBegin;
1199 
1200   F = snes->vec_func;
1201   B = snes->vec_rhs;
1202   Xhat = snes->work[3];
1203   G    = snes->work[1];
1204   W    = snes->work[2];
1205   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
1206   /* recurse first */
1207   if (fas->next) {
1208     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
1209 
1210     X_c  = fas->next->vec_sol;
1211     Xo_c = fas->next->work[0];
1212     F_c  = fas->next->vec_func;
1213     B_c  = fas->next->vec_rhs;
1214 
1215     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
1216     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1217 
1218     /* restrict the defect */
1219     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1220 
1221     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1222     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1223     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1224 
1225     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1226 
1227     /* set initial guess of the coarse problem to the projected fine solution */
1228     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1229 
1230     /* recurse */
1231     fas->next->vec_rhs = B_c;
1232     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1233 
1234     /* smooth on this level */
1235     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1236 
1237     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1238     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1239       snes->reason = SNES_DIVERGED_INNER;
1240       PetscFunctionReturn(0);
1241     }
1242 
1243     /* correct as x <- x + I(x^c - Rx)*/
1244     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1245     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
1246 
1247     /* additive correction of the coarse direction*/
1248     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1249     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1250     ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr);
1251     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1252     ierr = VecCopy(W, X);CHKERRQ(ierr);
1253     ierr = VecCopy(G, F);CHKERRQ(ierr);
1254     fnorm = gnorm;
1255   } else {
1256     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1257   }
1258   PetscFunctionReturn(0);
1259 }
1260 
1261 #undef __FUNCT__
1262 #define __FUNCT__ "FASCycle_Multiplicative"
1263 /*
1264 
1265 Defines the FAS cycle as:
1266 
1267 fine problem: F(x) = 0
1268 coarse problem: F^c(x) = b^c
1269 
1270 b^c = F^c(I^c_fx^f - I^c_fF(x))
1271 
1272 correction:
1273 
1274 x = x + I(x^c - Rx)
1275 
1276  */
1277 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
1278 
1279   PetscErrorCode      ierr;
1280   Vec                 F,B;
1281   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1282 
1283   PetscFunctionBegin;
1284   F = snes->vec_func;
1285   B = snes->vec_rhs;
1286   /* pre-smooth -- just update using the pre-smoother */
1287   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1288 
1289   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
1290 
1291   if (fas->level != 0) {
1292     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1293   }
1294   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1295 
1296   PetscFunctionReturn(0);
1297 }
1298 
1299 #undef __FUNCT__
1300 #define __FUNCT__ "SNESSolve_FAS"
1301 
1302 PetscErrorCode SNESSolve_FAS(SNES snes)
1303 {
1304   PetscErrorCode ierr;
1305   PetscInt       i, maxits;
1306   Vec            X, F;
1307   PetscReal      fnorm;
1308   SNES_FAS       *fas = (SNES_FAS *)snes->data;
1309   PetscFunctionBegin;
1310   maxits = snes->max_its;            /* maximum number of iterations */
1311   snes->reason = SNES_CONVERGED_ITERATING;
1312   X = snes->vec_sol;
1313   F = snes->vec_func;
1314 
1315   /*norm setup */
1316   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1317   snes->iter = 0;
1318   snes->norm = 0.;
1319   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1320   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1321   if (snes->domainerror) {
1322     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1323     PetscFunctionReturn(0);
1324   }
1325   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1326   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1327   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1328   snes->norm = fnorm;
1329   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1330   SNESLogConvHistory(snes,fnorm,0);
1331   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1332 
1333   /* set parameter for default relative tolerance convergence test */
1334   snes->ttol = fnorm*snes->rtol;
1335   /* test convergence */
1336   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1337   if (snes->reason) PetscFunctionReturn(0);
1338   for (i = 0; i < maxits; i++) {
1339     /* Call general purpose update function */
1340 
1341     if (snes->ops->update) {
1342       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1343     }
1344     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
1345       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
1346     } else {
1347       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
1348     }
1349 
1350     /* check for FAS cycle divergence */
1351     if (snes->reason != SNES_CONVERGED_ITERATING) {
1352       PetscFunctionReturn(0);
1353     }
1354     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1355     /* Monitor convergence */
1356     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1357     snes->iter = i+1;
1358     snes->norm = fnorm;
1359     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1360     SNESLogConvHistory(snes,snes->norm,0);
1361     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1362     /* Test for convergence */
1363     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1364     if (snes->reason) break;
1365   }
1366   if (i == maxits) {
1367     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1368     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1369   }
1370   PetscFunctionReturn(0);
1371 }
1372